[16] | 1 | !*********************************************************************** |
---|
| 2 | !* Copyright 2012,2013 * |
---|
| 3 | !* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine, * |
---|
| 4 | !* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,* |
---|
| 5 | !* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso, * |
---|
| 6 | !* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann, * |
---|
| 7 | !* * |
---|
| 8 | !* This file is part of FLEXPART WRF * |
---|
| 9 | !* * |
---|
| 10 | !* FLEXPART is free software: you can redistribute it and/or modify * |
---|
| 11 | !* it under the terms of the GNU General Public License as published by* |
---|
| 12 | !* the Free Software Foundation, either version 3 of the License, or * |
---|
| 13 | !* (at your option) any later version. * |
---|
| 14 | !* * |
---|
| 15 | !* FLEXPART is distributed in the hope that it will be useful, * |
---|
| 16 | !* but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
| 17 | !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
| 18 | !* GNU General Public License for more details. * |
---|
| 19 | !* * |
---|
| 20 | !* You should have received a copy of the GNU General Public License * |
---|
| 21 | !* along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
| 22 | !*********************************************************************** |
---|
| 23 | subroutine getfields(itime,nstop) |
---|
| 24 | ! i o |
---|
| 25 | !******************************************************************************* |
---|
| 26 | ! * |
---|
| 27 | ! This subroutine manages the 3 data fields to be kept in memory. * |
---|
| 28 | ! During the first time step of petterssen it has to be fulfilled that the * |
---|
| 29 | ! first data field must have |wftime|<itime, i.e. the absolute value of wftime* |
---|
| 30 | ! must be smaller than the absolute value of the current time in [s]. * |
---|
| 31 | ! The other 2 fields are the next in time after the first one. * |
---|
| 32 | ! Pointers (memind) are used, because otherwise one would have to resort the * |
---|
| 33 | ! wind fields, which costs a lot of computing time. Here only the pointers are* |
---|
| 34 | ! resorted. * |
---|
| 35 | ! * |
---|
| 36 | ! Author: A. Stohl * |
---|
| 37 | ! * |
---|
| 38 | ! 29 April 1994 * |
---|
| 39 | ! * |
---|
| 40 | ! Changes, Bernd C. Krueger, Feb. 2001: * |
---|
| 41 | ! Variables tth,qvh,tthn,qvhn (on eta coordinates) in common block. * |
---|
| 42 | ! Function of nstop extended. * |
---|
| 43 | ! * |
---|
| 44 | ! Dec 2005, R. Easter - * |
---|
| 45 | ! When "memtime(2) = itime = wftime(numbwf)", do not read a new file. * |
---|
| 46 | ! This allows the ending date/time of the flexpart run to match * |
---|
| 47 | ! the date/time of the last met. file. * |
---|
| 48 | ! * |
---|
| 49 | !******************************************************************************* |
---|
| 50 | ! * |
---|
| 51 | ! Variables: * |
---|
| 52 | ! lwindinterval [s] time difference between the two wind fields read in * |
---|
| 53 | ! indj indicates the number of the wind field to be read in * |
---|
| 54 | ! indmin remembers the number of wind fields already treated * |
---|
| 55 | ! memind(2) pointer, on which place the wind fields are stored * |
---|
| 56 | ! memtime(2) [s] times of the wind fields, which are kept in memory * |
---|
| 57 | ! itime [s] current time since start date of trajectory calculation * |
---|
| 58 | ! nstop > 0, if trajectory has to be terminated * |
---|
| 59 | ! nx,ny,nuvz,nwz field dimensions in x,y and z direction * |
---|
| 60 | ! uu(0:nxmax,0:nymax,nuvzmax,2) wind components in x-direction [m/s] * |
---|
| 61 | ! vv(0:nxmax,0:nymax,nuvzmax,2) wind components in y-direction [m/s] * |
---|
| 62 | ! ww(0:nxmax,0:nymax,nwzmax,2) wind components in z-direction [deltaeta/s] * |
---|
| 63 | ! tt(0:nxmax,0:nymax,nuvzmax,2) temperature [K] * |
---|
| 64 | ! ps(0:nxmax,0:nymax,2) surface pressure [Pa] * |
---|
| 65 | ! * |
---|
| 66 | ! Constants: * |
---|
| 67 | ! idiffmax maximum allowable time difference between 2 wind fields * |
---|
| 68 | ! * |
---|
| 69 | !******************************************************************************* |
---|
| 70 | |
---|
| 71 | use par_mod |
---|
| 72 | use com_mod |
---|
| 73 | |
---|
| 74 | implicit none |
---|
| 75 | |
---|
| 76 | integer :: indj,itime,nstop,memaux |
---|
| 77 | |
---|
| 78 | real(kind=4) :: uuh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 79 | real(kind=4) :: vvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 80 | real(kind=4) :: wwh(0:nxmax-1,0:nymax-1,nwzmax) |
---|
| 81 | real(kind=4) :: divh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 82 | |
---|
| 83 | ! real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 84 | ! real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 85 | real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 86 | ! real :: wwh(0:nxmax-1,0:nymax-1,nwzmax) |
---|
| 87 | real(kind=4) :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
| 88 | real(kind=4) :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
| 89 | real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
| 90 | real(kind=4) :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests) |
---|
| 91 | real :: start, finish |
---|
| 92 | |
---|
| 93 | integer :: indmin = 1 |
---|
| 94 | |
---|
| 95 | |
---|
| 96 | real(kind=4) :: divhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
| 97 | character(len=28) :: name2 |
---|
| 98 | character(len=8) :: chartime |
---|
| 99 | integer :: ix,jy,kz |
---|
| 100 | |
---|
| 101 | |
---|
| 102 | |
---|
| 103 | ! Check, if wind fields are available for the current time step |
---|
| 104 | !************************************************************** |
---|
| 105 | |
---|
| 106 | nstop=0 |
---|
| 107 | |
---|
| 108 | if ((ldirect*wftime(1).gt.ldirect*itime).or. & |
---|
| 109 | (ldirect*wftime(numbwf).lt.ldirect*itime)) then |
---|
| 110 | write(*,*) 'FLEXPART WARNING: NO WIND FIELDS ARE AVAILABLE.' |
---|
| 111 | write(*,*) 'A TRAJECTORY HAS TO BE TERMINATED.' |
---|
| 112 | nstop=4 |
---|
| 113 | return |
---|
| 114 | endif |
---|
| 115 | |
---|
| 116 | |
---|
| 117 | if ((ldirect*memtime(1).le.ldirect*itime).and. & |
---|
| 118 | (ldirect*memtime(2).gt.ldirect*itime)) then |
---|
| 119 | |
---|
| 120 | ! The right wind fields are already in memory -> don't do anything |
---|
| 121 | !***************************************************************** |
---|
| 122 | |
---|
| 123 | continue |
---|
| 124 | |
---|
| 125 | ! FLEXPART_WRF - following change allows the ending date/time |
---|
| 126 | ! of the flexpart run to match that of the last met. file |
---|
| 127 | else if ( (ldirect*memtime(1).lt.ldirect*itime).and. & |
---|
| 128 | (memtime(2).eq.itime) .and. (wftime(numbwf).eq.itime) ) then |
---|
| 129 | |
---|
| 130 | continue |
---|
| 131 | |
---|
| 132 | else if ((ldirect*memtime(2).le.ldirect*itime).and. & |
---|
| 133 | (memtime(2).ne.999999999)) then |
---|
| 134 | |
---|
| 135 | |
---|
| 136 | ! Current time is after 2nd wind field |
---|
| 137 | ! -> Resort wind field pointers, so that current time is between 1st and 2nd |
---|
| 138 | !*************************************************************************** |
---|
| 139 | |
---|
| 140 | memaux=memind(1) |
---|
| 141 | memind(1)=memind(2) |
---|
| 142 | memind(2)=memaux |
---|
| 143 | memtime(1)=memtime(2) |
---|
| 144 | |
---|
| 145 | |
---|
| 146 | ! Read a new wind field and store it on place memind(2) |
---|
| 147 | !****************************************************** |
---|
| 148 | |
---|
| 149 | do indj=indmin,numbwf-1 |
---|
| 150 | if (ldirect*wftime(indj+1).gt.ldirect*itime) then |
---|
| 151 | call cpu_time(start) |
---|
| 152 | if ((time_option.eq.1).and.(wind_option.eq.1)) & |
---|
| 153 | call readwind_timeav(indj+1,memind(2),uuh,vvh,wwh) |
---|
| 154 | call readwind(indj+1,memind(2),uuh,vvh,wwh,divh) |
---|
| 155 | if ((time_option.eq.1).and.(wind_option.eq.1)) & |
---|
| 156 | call readwind_nests_timeav(indj+1,memind(2),uuhn,vvhn,wwhn) |
---|
| 157 | call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn,divhn) |
---|
| 158 | if (option_verbose.gt.1) then |
---|
| 159 | call cpu_time(finish) |
---|
| 160 | print*,'readwind',finish-start |
---|
| 161 | call cpu_time(start) |
---|
| 162 | endif |
---|
| 163 | call calcpar(memind(2),uuh,vvh,pvh) |
---|
| 164 | call calcpar_nests(memind(2),uuhn,vvhn,pvhn) |
---|
| 165 | if (option_verbose.gt.1) then |
---|
| 166 | call cpu_time(finish) |
---|
| 167 | print*,'calcpar',finish-start |
---|
| 168 | call cpu_time(start) |
---|
| 169 | endif |
---|
| 170 | call verttransform(memind(2),uuh,vvh,wwh,pvh,divh) |
---|
| 171 | call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn,divhn) |
---|
| 172 | if (option_verbose.gt.1) then |
---|
| 173 | call cpu_time(finish) |
---|
| 174 | print*,'verttran',finish-start |
---|
| 175 | endif |
---|
| 176 | memtime(2)=wftime(indj+1) |
---|
| 177 | nstop = 1 |
---|
| 178 | goto 40 |
---|
| 179 | endif |
---|
| 180 | enddo |
---|
| 181 | 40 indmin=indj |
---|
| 182 | |
---|
| 183 | else |
---|
| 184 | |
---|
| 185 | ! No wind fields, which can be used, are currently in memory |
---|
| 186 | ! -> read both wind fields |
---|
| 187 | !*********************************************************** |
---|
| 188 | |
---|
| 189 | do indj=indmin,numbwf-1 |
---|
| 190 | if ((ldirect*wftime(indj).le.ldirect*itime).and. & |
---|
| 191 | (ldirect*wftime(indj+1).gt.ldirect*itime)) then |
---|
| 192 | memind(1)=1 |
---|
| 193 | if ((time_option.eq.1).and.(wind_option.eq.1)) & |
---|
| 194 | call readwind_timeav(indj,memind(1),uuh,vvh,wwh) |
---|
| 195 | call readwind(indj,memind(1),uuh,vvh,wwh,divh) |
---|
| 196 | if ((time_option.eq.1).and.(wind_option.eq.1)) & |
---|
| 197 | call readwind_nests_timeav(indj,memind(1),uuhn,vvhn,wwhn) |
---|
| 198 | call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn,divhn) |
---|
| 199 | call calcpar(memind(1),uuh,vvh,pvh) |
---|
| 200 | call calcpar_nests(memind(1),uuhn,vvhn,pvhn) |
---|
| 201 | call verttransform(memind(1),uuh,vvh,wwh,pvh,divh) |
---|
| 202 | call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn,divhn) |
---|
| 203 | memtime(1)=wftime(indj) |
---|
| 204 | memind(2)=2 |
---|
| 205 | if ((time_option.eq.1).and.(wind_option.eq.1)) & |
---|
| 206 | call readwind_timeav(indj+1,memind(2),uuh,vvh,wwh) |
---|
| 207 | call readwind(indj+1,memind(2),uuh,vvh,wwh,divh) |
---|
| 208 | if ((time_option.eq.1).and.(wind_option.eq.1)) & |
---|
| 209 | call readwind_nests_timeav(indj+1,memind(2),uuhn,vvhn,wwhn) |
---|
| 210 | call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn,divhn) |
---|
| 211 | call calcpar(memind(2),uuh,vvh,pvh) |
---|
| 212 | call calcpar_nests(memind(2),uuhn,vvhn,pvhn) |
---|
| 213 | call verttransform(memind(2),uuh,vvh,wwh,pvh,divh) |
---|
| 214 | call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn,divhn) |
---|
| 215 | memtime(2)=wftime(indj+1) |
---|
| 216 | nstop = 1 |
---|
| 217 | goto 60 |
---|
| 218 | endif |
---|
| 219 | end do |
---|
| 220 | 60 indmin=indj |
---|
| 221 | |
---|
| 222 | endif |
---|
| 223 | |
---|
| 224 | lwindinterv=abs(memtime(2)-memtime(1)) |
---|
| 225 | |
---|
| 226 | if (lwindinterv.gt.idiffmax) nstop=3 |
---|
| 227 | |
---|
| 228 | end subroutine getfields |
---|