source: flexpart.git/src/readavailable.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

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