[92fab65] | 1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
| 2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
[332fbbd] | 3 | |
---|
[6ecb30a] | 4 | subroutine getfields(itime,nstop,memstat,metdata_format) |
---|
[8a65cb0] | 5 | ! i o o |
---|
| 6 | !***************************************************************************** |
---|
| 7 | ! * |
---|
| 8 | ! This subroutine manages the 3 data fields to be kept in memory. * |
---|
| 9 | ! During the first time step of petterssen it has to be fulfilled that the * |
---|
| 10 | ! first data field must have |wftime|<itime, i.e. the absolute value of * |
---|
| 11 | ! wftime must be smaller than the absolute value of the current time in [s].* |
---|
| 12 | ! The other 2 fields are the next in time after the first one. * |
---|
| 13 | ! Pointers (memind) are used, because otherwise one would have to resort the* |
---|
| 14 | ! wind fields, which costs a lot of computing time. Here only the pointers * |
---|
| 15 | ! are resorted. * |
---|
| 16 | ! * |
---|
| 17 | ! Author: A. Stohl * |
---|
| 18 | ! * |
---|
| 19 | ! 29 April 1994 * |
---|
| 20 | ! * |
---|
| 21 | !***************************************************************************** |
---|
| 22 | ! Changes, Bernd C. Krueger, Feb. 2001: |
---|
| 23 | ! Variables tth,qvh,tthn,qvhn (on eta coordinates) in common block. |
---|
| 24 | ! Function of nstop extended. |
---|
| 25 | ! |
---|
| 26 | ! eso 2014: |
---|
[78e62dc] | 27 | ! MPI version. |
---|
| 28 | ! If running with number of processes >= mpi_mod::read_grp_min, |
---|
| 29 | ! only one process (mpi_mod::lmpreader=.true.) does the actual reading, while |
---|
| 30 | ! the rest call this routine only to update memind, memstat etc. |
---|
| 31 | ! |
---|
| 32 | ! If mpi_mod::lmp_sync=.true., uses 3 fields instead of 2, to allow reading |
---|
| 33 | ! the newest in the background. |
---|
[8a65cb0] | 34 | ! |
---|
| 35 | ! Return memstat, which is the sum of |
---|
| 36 | ! |
---|
| 37 | ! memstat=16: one new field to be read |
---|
| 38 | ! memstat=32: two new fields must be read |
---|
| 39 | ! memstat=1,2,3: position(s) in array to read next field |
---|
| 40 | ! memstat=0: no new fields to be read |
---|
| 41 | ! |
---|
[6ecb30a] | 42 | ! Unified ECMWF and GFS builds |
---|
| 43 | ! Marian Harustak, 12.5.2017 |
---|
| 44 | ! - Added passing of metdata_format as it was needed by called routines |
---|
| 45 | ! |
---|
[8a65cb0] | 46 | !***************************************************************************** |
---|
| 47 | ! * |
---|
| 48 | ! Variables: * |
---|
| 49 | ! lwindinterval [s] time difference between the two wind fields read in * |
---|
| 50 | ! indj indicates the number of the wind field to be read in * |
---|
| 51 | ! indmin remembers the number of wind fields already treated * |
---|
[78e62dc] | 52 | ! memind(2[3]) pointer, on which place the wind fields are stored * |
---|
| 53 | ! memtime(2[3]) [s] times of the wind fields, which are kept in memory * |
---|
[8a65cb0] | 54 | ! itime [s] current time since start date of trajectory calcu- * |
---|
| 55 | ! lation * |
---|
| 56 | ! nstop > 0, if trajectory has to be terminated * |
---|
| 57 | ! memstat keep track of change of status for field pointers * |
---|
| 58 | ! nx,ny,nuvz,nwz field dimensions in x,y and z direction * |
---|
| 59 | ! uu(0:nxmax,0:nymax,nuvzmax,numwfmem) wind components in x-direction [m/s]* |
---|
| 60 | ! vv(0:nxmax,0:nymax,nuvzmax,numwfmem) wind components in y-direction [m/s]* |
---|
| 61 | ! ww(0:nxmax,0:nymax,nwzmax,numwfmem) wind components in z-direction * |
---|
| 62 | ! [deltaeta/s] * |
---|
| 63 | ! tt(0:nxmax,0:nymax,nuvzmax,numwfmem) temperature [K] * |
---|
| 64 | ! ps(0:nxmax,0:nymax,numwfmem) surface pressure [Pa] * |
---|
[6ecb30a] | 65 | ! metdata_format format of metdata (ecmwf/gfs) * |
---|
[8a65cb0] | 66 | ! * |
---|
| 67 | ! Constants: * |
---|
| 68 | ! idiffmax maximum allowable time difference between 2 wind * |
---|
| 69 | ! fields * |
---|
| 70 | ! * |
---|
| 71 | !***************************************************************************** |
---|
| 72 | |
---|
| 73 | use par_mod |
---|
| 74 | use com_mod |
---|
| 75 | use mpi_mod, only: lmpreader,lmp_use_reader,lmp_sync |
---|
[6ecb30a] | 76 | use class_gribfile |
---|
[8a65cb0] | 77 | |
---|
| 78 | implicit none |
---|
| 79 | |
---|
[6ecb30a] | 80 | integer :: metdata_format |
---|
[8a65cb0] | 81 | integer :: indj,itime,nstop,memaux,mindread |
---|
| 82 | ! mindread: index of where to read 3rd field |
---|
| 83 | integer, intent(out) :: memstat |
---|
| 84 | ! keep track of 3rd field index. updated when new fields are read |
---|
| 85 | integer, save :: mind3=0 |
---|
| 86 | |
---|
| 87 | real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 88 | real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 89 | real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 90 | real :: wwh(0:nxmax-1,0:nymax-1,nwzmax) |
---|
| 91 | real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
| 92 | real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
| 93 | real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
| 94 | real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests) |
---|
| 95 | |
---|
| 96 | integer :: indmin = 1 |
---|
| 97 | |
---|
| 98 | |
---|
| 99 | ! Check, if wind fields are available for the current time step |
---|
| 100 | !************************************************************** |
---|
| 101 | |
---|
| 102 | nstop=0 |
---|
| 103 | memstat=0 |
---|
| 104 | |
---|
| 105 | if ((ldirect*wftime(1).gt.ldirect*itime).or. & |
---|
| 106 | (ldirect*wftime(numbwf).lt.ldirect*itime)) then |
---|
| 107 | write(*,*) 'FLEXPART WARNING: NO WIND FIELDS ARE AVAILABLE.' |
---|
| 108 | write(*,*) 'A TRAJECTORY HAS TO BE TERMINATED.' |
---|
| 109 | nstop=4 |
---|
| 110 | return |
---|
| 111 | endif |
---|
| 112 | |
---|
| 113 | |
---|
| 114 | if ((ldirect*memtime(1).le.ldirect*itime).and. & |
---|
| 115 | (ldirect*memtime(2).gt.ldirect*itime)) then |
---|
| 116 | |
---|
| 117 | ! The right wind fields are already in memory -> don't do anything |
---|
| 118 | !***************************************************************** |
---|
| 119 | memstat=0 |
---|
| 120 | continue |
---|
| 121 | |
---|
| 122 | else if ((ldirect*memtime(2).le.ldirect*itime).and. & |
---|
| 123 | (memtime(2).ne.999999999)) then |
---|
| 124 | |
---|
| 125 | |
---|
| 126 | ! Current time is after 2nd wind field |
---|
| 127 | ! -> Resort wind field pointers, so that current time is between 1st and 2nd |
---|
| 128 | !*************************************************************************** |
---|
| 129 | |
---|
| 130 | ! 2 fields in memory |
---|
| 131 | !******************* |
---|
| 132 | if (lmp_sync) then |
---|
| 133 | memaux=memind(1) |
---|
| 134 | memind(1)=memind(2) |
---|
| 135 | memind(2)=memaux |
---|
| 136 | memtime(1)=memtime(2) |
---|
| 137 | memstat=16+memind(2) |
---|
| 138 | mindread=memind(2) |
---|
| 139 | |
---|
| 140 | else |
---|
| 141 | ! 3 fields in memory |
---|
| 142 | ! Note that the reader process never needs to keep 3 fields in memory, |
---|
[5f9d14a] | 143 | ! (2 is enough) so it is possible to save some memory |
---|
[8a65cb0] | 144 | !********************************************************************* |
---|
| 145 | if (mind3.eq.32.or.mind3.eq.19) then |
---|
| 146 | if (lmpreader) then |
---|
| 147 | memind(1)=2 |
---|
| 148 | memind(2)=3 |
---|
| 149 | memind(3)=3 |
---|
| 150 | else |
---|
| 151 | memind(1)=2 |
---|
| 152 | memind(2)=3 |
---|
| 153 | memind(3)=1 |
---|
| 154 | end if |
---|
| 155 | memstat=17 |
---|
| 156 | else if (mind3.eq.17) then |
---|
| 157 | if (lmpreader) then |
---|
| 158 | memind(1)=3 |
---|
| 159 | memind(2)=1 |
---|
| 160 | memind(3)=1 |
---|
| 161 | else |
---|
| 162 | memind(1)=3 |
---|
| 163 | memind(2)=1 |
---|
| 164 | memind(3)=2 |
---|
| 165 | end if |
---|
| 166 | memstat=18 |
---|
| 167 | else if (mind3.eq.18) then |
---|
| 168 | if (lmpreader) then |
---|
| 169 | memind(1)=1 |
---|
| 170 | memind(2)=2 |
---|
| 171 | memind(3)=2 |
---|
| 172 | else |
---|
| 173 | memind(1)=1 |
---|
| 174 | memind(2)=2 |
---|
| 175 | memind(3)=3 |
---|
| 176 | end if |
---|
| 177 | memstat=19 |
---|
| 178 | else |
---|
| 179 | write(*,*) '#### getfields_mpi.f90> ERROR: ', & |
---|
| 180 | & 'unknown mind3, exiting.', mind3,' ####' |
---|
| 181 | stop |
---|
| 182 | end if |
---|
| 183 | mind3=memstat |
---|
| 184 | memtime(1)=memtime(2) |
---|
| 185 | mindread=memind(3) |
---|
| 186 | end if |
---|
| 187 | |
---|
| 188 | |
---|
| 189 | ! Read a new wind field and store it on place memind(2) |
---|
| 190 | ! or memind(3) for the 3-field read-ahead version |
---|
| 191 | !****************************************************** |
---|
| 192 | do indj=indmin,numbwf-1 |
---|
| 193 | if (ldirect*wftime(indj+1).gt.ldirect*itime) then |
---|
| 194 | if (lmpreader.or..not. lmp_use_reader) then |
---|
[6ecb30a] | 195 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
| 196 | call readwind_ecmwf(indj+1,mindread,uuh,vvh,wwh) |
---|
| 197 | else |
---|
| 198 | call readwind_gfs(indj+1,mindread,uuh,vvh,wwh) |
---|
| 199 | end if |
---|
[8a65cb0] | 200 | call readwind_nests(indj+1,mindread,uuhn,vvhn,wwhn) |
---|
[6ecb30a] | 201 | call calcpar(mindread,uuh,vvh,pvh,metdata_format) |
---|
| 202 | call calcpar_nests(mindread,uuhn,vvhn,pvhn,metdata_format) |
---|
| 203 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
| 204 | call verttransform_ecmwf(mindread,uuh,vvh,wwh,pvh) |
---|
| 205 | else |
---|
| 206 | call verttransform_gfs(mindread,uuh,vvh,wwh,pvh) |
---|
| 207 | end if |
---|
[8a65cb0] | 208 | call verttransform_nests(mindread,uuhn,vvhn,wwhn,pvhn) |
---|
| 209 | end if |
---|
| 210 | memtime(2)=wftime(indj+1) |
---|
| 211 | nstop = 1 |
---|
| 212 | goto 40 |
---|
| 213 | endif |
---|
| 214 | end do |
---|
| 215 | 40 indmin=indj |
---|
| 216 | |
---|
[20963b1] | 217 | if (WETBKDEP.and.(lmpreader.or.(.not.lmp_use_reader.and.lroot))) then |
---|
| 218 | call writeprecip(itime,memind(1)) |
---|
| 219 | endif |
---|
| 220 | |
---|
[8a65cb0] | 221 | else |
---|
| 222 | |
---|
| 223 | ! No wind fields, which can be used, are currently in memory |
---|
| 224 | ! -> read both wind fields |
---|
| 225 | !*********************************************************** |
---|
| 226 | memstat=32 |
---|
| 227 | |
---|
| 228 | do indj=indmin,numbwf-1 |
---|
| 229 | if ((ldirect*wftime(indj).le.ldirect*itime).and. & |
---|
| 230 | (ldirect*wftime(indj+1).gt.ldirect*itime)) then |
---|
| 231 | memind(1)=1 |
---|
| 232 | if (lmpreader.or..not.lmp_use_reader) then |
---|
[6ecb30a] | 233 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
| 234 | call readwind_ecmwf(indj,memind(1),uuh,vvh,wwh) |
---|
| 235 | else |
---|
| 236 | call readwind_gfs(indj,memind(1),uuh,vvh,wwh) |
---|
| 237 | end if |
---|
[8a65cb0] | 238 | call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn) |
---|
[6ecb30a] | 239 | call calcpar(memind(1),uuh,vvh,pvh,metdata_format) |
---|
| 240 | call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format) |
---|
| 241 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
| 242 | call verttransform_ecmwf(memind(1),uuh,vvh,wwh,pvh) |
---|
| 243 | else |
---|
| 244 | call verttransform_gfs(memind(1),uuh,vvh,wwh,pvh) |
---|
| 245 | end if |
---|
[8a65cb0] | 246 | call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn) |
---|
| 247 | end if |
---|
| 248 | memtime(1)=wftime(indj) |
---|
| 249 | memind(2)=2 |
---|
| 250 | if (lmpreader.or..not.lmp_use_reader) then |
---|
[6ecb30a] | 251 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
| 252 | call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh) |
---|
| 253 | else |
---|
| 254 | call readwind_gfs(indj+1,memind(2),uuh,vvh,wwh) |
---|
| 255 | end if |
---|
[8a65cb0] | 256 | call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn) |
---|
[6ecb30a] | 257 | call calcpar(memind(2),uuh,vvh,pvh,metdata_format) |
---|
| 258 | call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format) |
---|
| 259 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
| 260 | call verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh) |
---|
| 261 | else |
---|
| 262 | call verttransform_gfs(memind(2),uuh,vvh,wwh,pvh) |
---|
| 263 | end if |
---|
[8a65cb0] | 264 | call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn) |
---|
| 265 | end if |
---|
| 266 | memtime(2)=wftime(indj+1) |
---|
[5f9d14a] | 267 | ! DEV: not used |
---|
[8a65cb0] | 268 | if (.not.lmp_sync) memind(3)=3 ! indicate position of next read |
---|
| 269 | nstop = 1 |
---|
| 270 | goto 60 |
---|
| 271 | endif |
---|
| 272 | end do |
---|
| 273 | 60 indmin=indj |
---|
| 274 | |
---|
[20963b1] | 275 | ! if (WETBKDEP.and.lroot) then |
---|
| 276 | if (WETBKDEP.and.(lmpreader.or.(.not.lmp_use_reader.and.lroot))) then |
---|
| 277 | call writeprecip(itime,memind(1)) |
---|
| 278 | endif |
---|
[8a65cb0] | 279 | |
---|
| 280 | endif |
---|
| 281 | |
---|
| 282 | lwindinterv=abs(memtime(2)-memtime(1)) |
---|
| 283 | |
---|
| 284 | if (lwindinterv.gt.idiffmax) nstop=3 |
---|
| 285 | |
---|
| 286 | end subroutine getfields |
---|