[e200b7a] | 1 | subroutine readavailable |
---|
| 2 | |
---|
| 3 | !***************************************************************************** |
---|
| 4 | ! * |
---|
| 5 | ! This routine reads the dates and times for which windfields are * |
---|
| 6 | ! available. * |
---|
| 7 | ! * |
---|
| 8 | ! Authors: A. Stohl * |
---|
| 9 | ! * |
---|
| 10 | ! 6 February 1994 * |
---|
| 11 | ! 8 February 1999, Use of nested fields, A. Stohl * |
---|
| 12 | ! * |
---|
| 13 | !***************************************************************************** |
---|
| 14 | ! * |
---|
| 15 | ! Variables: * |
---|
| 16 | ! bdate beginning date as Julian date * |
---|
| 17 | ! beg beginning date for windfields * |
---|
| 18 | ! end ending date for windfields * |
---|
| 19 | ! fname filename of wind field, help variable * |
---|
| 20 | ! ideltas [s] duration of modelling period * |
---|
| 21 | ! idiff time difference between 2 wind fields * |
---|
| 22 | ! idiffnorm normal time difference between 2 wind fields * |
---|
| 23 | ! idiffmax [s] maximum allowable time between 2 wind fields * |
---|
| 24 | ! jul julian date, help variable * |
---|
| 25 | ! numbwf actual number of wind fields * |
---|
| 26 | ! wfname(maxwf) file names of needed wind fields * |
---|
| 27 | ! wfspec(maxwf) file specifications of wind fields (e.g., if on disc) * |
---|
| 28 | ! wftime(maxwf) [s]times of wind fields relative to beginning time * |
---|
| 29 | ! wfname1,wfspec1,wftime1 = same as above, but only local (help variables) * |
---|
| 30 | ! * |
---|
| 31 | ! Constants: * |
---|
| 32 | ! maxwf maximum number of wind fields * |
---|
| 33 | ! unitavailab unit connected to file AVAILABLE * |
---|
| 34 | ! * |
---|
| 35 | !***************************************************************************** |
---|
| 36 | |
---|
| 37 | use par_mod |
---|
| 38 | use com_mod |
---|
| 39 | |
---|
| 40 | implicit none |
---|
| 41 | |
---|
| 42 | integer :: i,idiff,ldat,ltim,wftime1(maxwf),numbwfn(maxnests),k |
---|
| 43 | integer :: wftime1n(maxnests,maxwf),wftimen(maxnests,maxwf) |
---|
[b5127f9] | 44 | logical :: lwarntd=.true. |
---|
[e200b7a] | 45 | real(kind=dp) :: juldate,jul,beg,end |
---|
| 46 | character(len=255) :: fname,spec,wfname1(maxwf),wfspec1(maxwf) |
---|
| 47 | character(len=255) :: wfname1n(maxnests,maxwf) |
---|
| 48 | character(len=40) :: wfspec1n(maxnests,maxwf) |
---|
| 49 | |
---|
| 50 | |
---|
| 51 | ! Windfields are only used, if they are within the modelling period. |
---|
| 52 | ! However, 1 additional day at the beginning and at the end is used for |
---|
| 53 | ! interpolation. -> Compute beginning and ending date for the windfields. |
---|
| 54 | !************************************************************************ |
---|
| 55 | |
---|
| 56 | if (ideltas.gt.0) then ! forward trajectories |
---|
| 57 | beg=bdate-1._dp |
---|
| 58 | end=bdate+real(ideltas,kind=dp)/86400._dp+real(idiffmax,kind=dp)/ & |
---|
| 59 | 86400._dp |
---|
| 60 | else ! backward trajectories |
---|
| 61 | beg=bdate+real(ideltas,kind=dp)/86400._dp-real(idiffmax,kind=dp)/ & |
---|
| 62 | 86400._dp |
---|
| 63 | end=bdate+1._dp |
---|
| 64 | endif |
---|
| 65 | |
---|
| 66 | ! Open the wind field availability file and read available wind fields |
---|
| 67 | ! within the modelling period. |
---|
| 68 | !********************************************************************* |
---|
| 69 | |
---|
| 70 | open(unitavailab,file=path(4)(1:length(4)),status='old', & |
---|
| 71 | err=999) |
---|
| 72 | |
---|
| 73 | do i=1,3 |
---|
| 74 | read(unitavailab,*) |
---|
| 75 | end do |
---|
| 76 | |
---|
| 77 | numbwf=0 |
---|
| 78 | 100 read(unitavailab,'(i8,1x,i6,2(6x,a255))',end=99) & |
---|
| 79 | ldat,ltim,fname,spec |
---|
| 80 | jul=juldate(ldat,ltim) |
---|
| 81 | if ((jul.ge.beg).and.(jul.le.end)) then |
---|
| 82 | numbwf=numbwf+1 |
---|
| 83 | if (numbwf.gt.maxwf) then ! check exceedance of dimension |
---|
| 84 | write(*,*) 'Number of wind fields needed is too great.' |
---|
| 85 | write(*,*) 'Reduce modelling period (file "COMMAND") or' |
---|
| 86 | write(*,*) 'reduce number of wind fields (file "AVAILABLE").' |
---|
| 87 | stop |
---|
| 88 | endif |
---|
| 89 | |
---|
| 90 | wfname1(numbwf)=fname(1:index(fname,' ')) |
---|
| 91 | wfspec1(numbwf)=spec |
---|
| 92 | wftime1(numbwf)=nint((jul-bdate)*86400._dp) |
---|
| 93 | endif |
---|
| 94 | goto 100 ! next wind field |
---|
| 95 | |
---|
| 96 | 99 continue |
---|
| 97 | |
---|
| 98 | close(unitavailab) |
---|
| 99 | |
---|
| 100 | ! Open the wind field availability file and read available wind fields |
---|
| 101 | ! within the modelling period (nested grids) |
---|
| 102 | !********************************************************************* |
---|
| 103 | |
---|
| 104 | do k=1,numbnests |
---|
[b4d29ce] | 105 | !print*,length(numpath+2*(k-1)+1),length(numpath+2*(k-1)+2),length(4),length(3) |
---|
| 106 | !print*,path(numpath+2*(k-1)+2)(1:length(numpath+2*(k-1)+2)) |
---|
[e200b7a] | 107 | open(unitavailab,file=path(numpath+2*(k-1)+2) & |
---|
| 108 | (1:length(numpath+2*(k-1)+2)),status='old',err=998) |
---|
| 109 | |
---|
| 110 | do i=1,3 |
---|
| 111 | read(unitavailab,*) |
---|
| 112 | end do |
---|
| 113 | |
---|
| 114 | numbwfn(k)=0 |
---|
| 115 | 700 read(unitavailab,'(i8,1x,i6,2(6x,a255))',end=699) ldat, & |
---|
| 116 | ltim,fname,spec |
---|
| 117 | jul=juldate(ldat,ltim) |
---|
| 118 | if ((jul.ge.beg).and.(jul.le.end)) then |
---|
| 119 | numbwfn(k)=numbwfn(k)+1 |
---|
| 120 | if (numbwfn(k).gt.maxwf) then ! check exceedance of dimension |
---|
| 121 | write(*,*) 'Number of nested wind fields is too great.' |
---|
| 122 | write(*,*) 'Reduce modelling period (file "COMMAND") or' |
---|
| 123 | write(*,*) 'reduce number of wind fields (file "AVAILABLE").' |
---|
| 124 | stop |
---|
| 125 | endif |
---|
| 126 | |
---|
| 127 | wfname1n(k,numbwfn(k))=fname |
---|
| 128 | wfspec1n(k,numbwfn(k))=spec |
---|
| 129 | wftime1n(k,numbwfn(k))=nint((jul-bdate)*86400._dp) |
---|
| 130 | endif |
---|
| 131 | goto 700 ! next wind field |
---|
| 132 | |
---|
| 133 | 699 continue |
---|
| 134 | |
---|
| 135 | close(unitavailab) |
---|
| 136 | end do |
---|
| 137 | |
---|
| 138 | |
---|
| 139 | ! Check wind field times of file AVAILABLE (expected to be in temporal order) |
---|
| 140 | !**************************************************************************** |
---|
| 141 | |
---|
| 142 | if (numbwf.eq.0) then |
---|
| 143 | write(*,*) ' #### FLEXPART MODEL ERROR! NO WIND FIELDS #### ' |
---|
| 144 | write(*,*) ' #### AVAILABLE FOR SELECTED TIME PERIOD. #### ' |
---|
| 145 | stop |
---|
| 146 | endif |
---|
| 147 | |
---|
| 148 | do i=2,numbwf |
---|
| 149 | if (wftime1(i).le.wftime1(i-1)) then |
---|
| 150 | write(*,*) 'FLEXPART ERROR: FILE AVAILABLE IS CORRUPT.' |
---|
| 151 | write(*,*) 'THE WIND FIELDS ARE NOT IN TEMPORAL ORDER.' |
---|
| 152 | write(*,*) 'PLEASE CHECK FIELD ',wfname1(i) |
---|
| 153 | stop |
---|
| 154 | endif |
---|
| 155 | end do |
---|
| 156 | |
---|
| 157 | ! Check wind field times of file AVAILABLE for the nested fields |
---|
| 158 | ! (expected to be in temporal order) |
---|
| 159 | !*************************************************************** |
---|
| 160 | |
---|
| 161 | do k=1,numbnests |
---|
| 162 | if (numbwfn(k).eq.0) then |
---|
| 163 | write(*,*) '#### FLEXPART MODEL ERROR! NO WIND FIELDS ####' |
---|
| 164 | write(*,*) '#### AVAILABLE FOR SELECTED TIME PERIOD. ####' |
---|
| 165 | stop |
---|
| 166 | endif |
---|
| 167 | |
---|
| 168 | do i=2,numbwfn(k) |
---|
| 169 | if (wftime1n(k,i).le.wftime1n(k,i-1)) then |
---|
| 170 | write(*,*) 'FLEXPART ERROR: FILE AVAILABLE IS CORRUPT. ' |
---|
| 171 | write(*,*) 'THE NESTED WIND FIELDS ARE NOT IN TEMPORAL ORDER.' |
---|
| 172 | write(*,*) 'PLEASE CHECK FIELD ',wfname1n(k,i) |
---|
| 173 | write(*,*) 'AT NESTING LEVEL ',k |
---|
| 174 | stop |
---|
| 175 | endif |
---|
| 176 | end do |
---|
| 177 | |
---|
| 178 | end do |
---|
| 179 | |
---|
| 180 | |
---|
| 181 | ! For backward trajectories, reverse the order of the windfields |
---|
| 182 | !*************************************************************** |
---|
| 183 | |
---|
| 184 | if (ideltas.ge.0) then |
---|
| 185 | do i=1,numbwf |
---|
| 186 | wfname(i)=wfname1(i) |
---|
| 187 | wfspec(i)=wfspec1(i) |
---|
| 188 | wftime(i)=wftime1(i) |
---|
| 189 | end do |
---|
| 190 | do k=1,numbnests |
---|
| 191 | do i=1,numbwfn(k) |
---|
| 192 | wfnamen(k,i)=wfname1n(k,i) |
---|
| 193 | wfspecn(k,i)=wfspec1n(k,i) |
---|
| 194 | wftimen(k,i)=wftime1n(k,i) |
---|
| 195 | end do |
---|
| 196 | end do |
---|
| 197 | else |
---|
| 198 | do i=1,numbwf |
---|
| 199 | wfname(numbwf-i+1)=wfname1(i) |
---|
| 200 | wfspec(numbwf-i+1)=wfspec1(i) |
---|
| 201 | wftime(numbwf-i+1)=wftime1(i) |
---|
| 202 | end do |
---|
| 203 | do k=1,numbnests |
---|
| 204 | do i=1,numbwfn(k) |
---|
| 205 | wfnamen(k,numbwfn(k)-i+1)=wfname1n(k,i) |
---|
| 206 | wfspecn(k,numbwfn(k)-i+1)=wfspec1n(k,i) |
---|
| 207 | wftimen(k,numbwfn(k)-i+1)=wftime1n(k,i) |
---|
| 208 | end do |
---|
| 209 | end do |
---|
| 210 | endif |
---|
| 211 | |
---|
| 212 | ! Check the time difference between the wind fields. If it is big, |
---|
| 213 | ! write a warning message. If it is too big, terminate the trajectory. |
---|
| 214 | !********************************************************************* |
---|
| 215 | |
---|
| 216 | do i=2,numbwf |
---|
| 217 | idiff=abs(wftime(i)-wftime(i-1)) |
---|
[4c64400] | 218 | if (idiff.gt.idiffmax.and.lroot) then |
---|
[e200b7a] | 219 | write(*,*) 'FLEXPART WARNING: TIME DIFFERENCE BETWEEN TWO' |
---|
| 220 | write(*,*) 'WIND FIELDS IS TOO BIG FOR TRANSPORT CALCULATION.& |
---|
| 221 | &' |
---|
| 222 | write(*,*) 'THEREFORE, TRAJECTORIES HAVE TO BE SKIPPED.' |
---|
[b5127f9] | 223 | else if (idiff.gt.idiffnorm.and.lroot.and.lwarntd) then |
---|
[e200b7a] | 224 | write(*,*) 'FLEXPART WARNING: TIME DIFFERENCE BETWEEN TWO' |
---|
| 225 | write(*,*) 'WIND FIELDS IS BIG. THIS MAY CAUSE A DEGRADATION' |
---|
| 226 | write(*,*) 'OF SIMULATION QUALITY.' |
---|
[b5127f9] | 227 | lwarntd=.false. ! only issue this warning once |
---|
[e200b7a] | 228 | endif |
---|
| 229 | end do |
---|
| 230 | |
---|
| 231 | do k=1,numbnests |
---|
| 232 | if (numbwfn(k).ne.numbwf) then |
---|
| 233 | write(*,*) 'FLEXPART ERROR: THE AVAILABLE FILES FOR THE' |
---|
| 234 | write(*,*) 'NESTED WIND FIELDS ARE NOT CONSISTENT WITH' |
---|
| 235 | write(*,*) 'THE AVAILABLE FILE OF THE MOTHER DOMAIN. ' |
---|
| 236 | write(*,*) 'ERROR AT NEST LEVEL: ',k |
---|
| 237 | stop |
---|
| 238 | endif |
---|
| 239 | do i=1,numbwf |
---|
| 240 | if (wftimen(k,i).ne.wftime(i)) then |
---|
| 241 | write(*,*) 'FLEXPART ERROR: THE AVAILABLE FILES FOR THE' |
---|
| 242 | write(*,*) 'NESTED WIND FIELDS ARE NOT CONSISTENT WITH' |
---|
| 243 | write(*,*) 'THE AVAILABLE FILE OF THE MOTHER DOMAIN. ' |
---|
| 244 | write(*,*) 'ERROR AT NEST LEVEL: ',k |
---|
| 245 | stop |
---|
| 246 | endif |
---|
| 247 | end do |
---|
| 248 | end do |
---|
| 249 | |
---|
| 250 | ! Reset the times of the wind fields that are kept in memory to no time |
---|
| 251 | !********************************************************************** |
---|
| 252 | |
---|
| 253 | do i=1,2 |
---|
| 254 | memind(i)=i |
---|
| 255 | memtime(i)=999999999 |
---|
| 256 | end do |
---|
| 257 | |
---|
| 258 | return |
---|
| 259 | |
---|
[f13406c] | 260 | 998 write(*,*) ' #### FLEXPART MODEL ERROR! AVAILABLE FILE #### ' |
---|
[e200b7a] | 261 | write(*,'(a)') ' '//path(numpath+2*(k-1)+2) & |
---|
| 262 | (1:length(numpath+2*(k-1)+2)) |
---|
| 263 | write(*,*) ' #### CANNOT BE OPENED #### ' |
---|
| 264 | stop |
---|
| 265 | |
---|
[db712a8] | 266 | 999 write(*,*) ' #### FLEXPART MODEL ERROR! AVAILABLE FILE #### ' |
---|
[e200b7a] | 267 | write(*,'(a)') ' '//path(4)(1:length(4)) |
---|
| 268 | write(*,*) ' #### CANNOT BE OPENED #### ' |
---|
| 269 | stop |
---|
| 270 | |
---|
| 271 | end subroutine readavailable |
---|