1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359
|
!-------------------------------------------------------------------------------
! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2016 EDF S.A.
!
! This program is free software; you can redistribute it and/or modify it under
! the terms of the GNU General Public License as published by the Free Software
! Foundation; either version 2 of the License, or (at your option) any later
! version.
!
! This program is distributed in the hope that it will be useful, but WITHOUT
! ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
! FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
! details.
!
! You should have received a copy of the GNU General Public License along with
! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
! Street, Fifth Floor6b1e433ea7e612c18f94b936f94db2eba9e887f2, Boston, MA 02110-1301, USA.
!-------------------------------------------------------------------------------
!> \file ecrhis.f90
!> \brief Write plot data
!>
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
! mode name role
!------------------------------------------------------------------------------
!> \param[in] modhis 0 or 1: initialize/output; 2: finalize
!______________________________________________________________________________
subroutine ecrhis &
( modhis )
!===============================================================================
! Module files
!===============================================================================
use paramx
use numvar
use entsor
use cstnum
use optcal
use parall
use mesh
use field
!===============================================================================
implicit none
! Arguments
integer modhis
! Local variables
character nompre*300, nomhis*300
integer tplnum, ii, lpre, lnom, lng
integer icap, ncap, ipp, ippf
integer c_id, f_id, f_dim, f_type, n_fields
double precision varcap(ncaptm)
integer, dimension(:), allocatable :: lsttmp
double precision, dimension(:,:), allocatable :: xyztmp
double precision, dimension(:,:), pointer :: val_v
double precision, dimension(:), pointer :: val_s
character(len=3), dimension(3), save :: nomext3 = (/'[X]', '[Y]', '[Z]'/)
character(len=4), dimension(6), save :: nomext6 &
= (/'[XX]', '[YY]', '[ZZ]', '[XY]', '[YZ]', '[XZ]'/)
character(len=4), dimension(9), save :: nomext9 &
= (/'[XX]', '[XY]', '[XZ]', '[YX]', '[YY]', '[YZ]', '[ZX]', '[ZY]', '[ZZ]'/)
integer, save :: keypp = -1
! Time plot number shift (in case multiple routines define plots)
integer nptpl
data nptpl /0/
save nptpl
! Number of passes in this routine
integer ipass
data ipass /0/
save ipass
!===============================================================================
! 0. Local initializations
!===============================================================================
ipass = ipass + 1
!--> If no probe data has been output or there are no probes, return
if ((ipass.eq.1 .and. modhis.eq.2) .or. ncapt.eq.0) return
if (ipass.eq.1) then
call tplnbr(nptpl)
endif
if (keypp.lt.0) then
call field_get_key_id("post_id", keypp)
endif
!===============================================================================
! 1. Search for neighboring nodes -> nodcap
!===============================================================================
if (ipass.eq.1) then
do ii = 1, ncapt
call findpt &
(ncelet, ncel, xyzcen, &
xyzcap(1,ii), xyzcap(2,ii), xyzcap(3,ii), nodcap(ii), ndrcap(ii))
enddo
endif
!===============================================================================
! 2. Initialize output
!===============================================================================
! Create directory if required
if (ipass.eq.1 .and. irangp.le.0) then
call csmkdr(emphis, len(emphis))
endif
if (ipass.eq.1) then
! plot prefix
nompre = trim(adjustl(emphis)) // trim(adjustl(prehis))
lpre = len_trim(nompre)
allocate(lsttmp(ncapt))
allocate(xyztmp(3, ncapt))
! Initialize one output per variable
call field_get_n_fields(n_fields)
do f_id = 0, n_fields - 1
call field_get_key_int(f_id, keypp, ippf)
if (ippf.le.1) cycle
call field_get_dim (f_id, f_dim)
do c_id = 1, min(f_dim, 9)
ipp = ippf + c_id - 1
if (ihisvr(ipp,1) .gt. 0) then
ncap = ihisvr(ipp,1)
do ii=1, ncap
lsttmp(ii) = ihisvr(ipp,ii+1)
if (irangp.lt.0 .or. irangp.eq.ndrcap(ihisvr(ipp, ii+1))) then
xyztmp(1, ii) = xyzcen(1, nodcap(ihisvr(ipp, ii+1)))
xyztmp(2, ii) = xyzcen(2, nodcap(ihisvr(ipp, ii+1)))
xyztmp(3, ii) = xyzcen(3, nodcap(ihisvr(ipp, ii+1)))
endif
if (irangp.ge.0) then
lng = 3
call parbcr(ndrcap(ihisvr(ipp,ii+1)), lng , xyztmp(1, ii))
endif
enddo
else if (ihisvr(ipp,1) .lt. 0) then
ncap = ncapt
do ii = 1, ncap
lsttmp(ii) = ii
if (irangp.lt.0 .or. irangp.eq.ndrcap(ii)) then
xyztmp(1, ii) = xyzcen(1, nodcap(ii))
xyztmp(2, ii) = xyzcen(2, nodcap(ii))
xyztmp(3, ii) = xyzcen(3, nodcap(ii))
endif
if (irangp.ge.0) then
lng = 3
call parbcr(ndrcap(ii), lng , xyztmp(1, ii))
endif
enddo
else
ncap = 0
endif
if (ncap .gt. 0) then
if (irangp.le.0) then
call field_get_label(f_id, nomhis)
nomhis = adjustl(nomhis)
if (f_dim .eq. 3) then
nomhis = trim(nomhis) // nomext3(c_id)
else if (f_dim .eq. 6) then
nomhis = trim(nomhis) // nomext6(c_id)
else if (f_dim .eq. 9) then
nomhis = trim(nomhis) // nomext9(c_id)
endif
lnom = len_trim(nomhis)
tplnum = nptpl + ipp
call tppini(tplnum, nomhis, nompre, tplfmt, idtvar, &
ncap, lsttmp(1), xyzcap(1,1), lnom, lpre)
endif ! (irangp.le.0)
endif
enddo
enddo
deallocate(lsttmp)
deallocate(xyztmp)
endif
!===============================================================================
! 3. Output results
!===============================================================================
if (modhis.eq.0 .or. modhis.eq.1) then
call field_get_n_fields(n_fields)
! Loop on fields
do f_id = 0, n_fields - 1
call field_get_key_int(f_id, keypp, ippf)
if (ippf.lt.2) cycle
call field_get_dim (f_id, f_dim)
call field_get_type(f_id, f_type)
do c_id = 1, min(f_dim, 9)
ipp = ippf + c_id - 1
if (ihisvr(ipp,1).ne.0) then
! Case of 1D fields, including moments
if (f_dim .eq. 1) then
call field_get_val_s(f_id, val_s)
if (ihisvr(ipp,1).lt.0) then
do icap = 1, ncapt
if (irangp.lt.0) then
varcap(icap) = val_s(nodcap(icap))
else
call parhis(nodcap(icap), ndrcap(icap), val_s, varcap(icap))
endif
enddo
ncap = ncapt
else
do icap = 1, ihisvr(ipp,1)
if (irangp.lt.0) then
varcap(icap) = val_s(nodcap(ihisvr(ipp,icap+1)))
else
call parhis(nodcap(ihisvr(ipp,icap+1)), &
ndrcap(ihisvr(ipp,icap+1)), &
val_s, varcap(icap))
endif
enddo
ncap = ihisvr(ipp,1)
endif
if (irangp.le.0 .and. ncap.gt.0) then
tplnum = nptpl + ipp
call tplwri(tplnum, tplfmt, ncap, ntcabs, ttcabs, varcap)
endif
else ! For vector field
call field_get_val_v(f_id, val_v)
if (ihisvr(ipp,1).lt.0) then
do icap = 1, ncapt
if (irangp.lt.0 .or. ndrcap(icap).eq.irangp) then
varcap(icap) = val_v(c_id, nodcap(icap))
endif
if (irangp.ge.0) then
lng = 1
call parbcr(ndrcap(icap), lng, varcap(icap))
endif
enddo
ncap = ncapt
else if (ihisvr(ipp,1).gt.0) then
do icap = 1, ihisvr(ipp,1)
if (irangp.lt.0 .or. ndrcap(icap).eq.irangp) then
varcap(icap) = val_v(c_id, nodcap(ihisvr(ipp,icap+1)))
endif
if (irangp.ge.0) then
lng = 1
call parbcr(ndrcap(icap), lng, varcap(icap))
endif
enddo
ncap = ihisvr(ipp,1)
else
ncap = 0
endif
if (irangp.le.0 .and. ncap.gt.0) then
tplnum = nptpl + ipp
call tplwri(tplnum, tplfmt, ncap, ntcabs, ttcabs, varcap)
endif
endif ! Scalar or vector field
endif ! Output for this component
enddo ! loop on components
enddo ! loop on fields
endif
!===============================================================================
! 4. Close output
!===============================================================================
if (modhis.eq.2) then
call field_get_n_fields(n_fields)
do f_id = 0, n_fields - 1
call field_get_key_int(f_id, keypp, ippf)
if (ippf.le.1) cycle
call field_get_dim (f_id, f_dim)
do c_id = 1, min(f_dim, 3)
ipp = ippf + c_id - 1
ncap = ihisvr(ipp,1)
if (ncap.lt.0) ncap = ncapt
if (ncap.gt.0 .and. irangp.le.0) then
tplnum = nptpl + ipp
call tplend(tplnum, tplfmt)
endif
enddo
enddo
endif
!===============================================================================
! 5. End
!===============================================================================
return
end subroutine
|