source: trunk/src/readavailable.f90 @ 28

Last change on this file since 28 was 27, checked in by hasod, 10 years ago
  • Implemented optional namelist input for COMMAND, RELEASES, SPECIES, AGECLASSES,OUTGRID,OUTGRID_NEST,RECEPTORS
  • Implemented com_mod switch nmlout to write input files as namelist to the output directory (.true. by default)
  • Proposed updated startup and runtime output (may change back to previous info if desired)
  • Property svn:executable set to *
File size: 11.3 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine readavailable
23
24  !*****************************************************************************
25  !                                                                            *
26  !   This routine reads the dates and times for which windfields are          *
27  !   available.                                                               *
28  !                                                                            *
29  !     Authors: A. Stohl                                                      *
30  !                                                                            *
31  !     6 February 1994                                                        *
32  !     8 February 1999, Use of nested fields, A. Stohl                        *
33  !                                                                            *
34  !*****************************************************************************
35  !                                                                            *
36  ! Variables:                                                                 *
37  ! bdate                beginning date as Julian date                         *
38  ! beg                  beginning date for windfields                         *
39  ! end                  ending date for windfields                            *
40  ! fname                filename of wind field, help variable                 *
41  ! ideltas [s]          duration of modelling period                          *
42  ! idiff                time difference between 2 wind fields                 *
43  ! idiffnorm            normal time difference between 2 wind fields          *
44  ! idiffmax [s]         maximum allowable time between 2 wind fields          *
45  ! jul                  julian date, help variable                            *
46  ! numbwf               actual number of wind fields                          *
47  ! wfname(maxwf)        file names of needed wind fields                      *
48  ! wfspec(maxwf)        file specifications of wind fields (e.g., if on disc) *
49  ! wftime(maxwf) [s]times of wind fields relative to beginning time           *
50  ! wfname1,wfspec1,wftime1 = same as above, but only local (help variables)   *
51  !                                                                            *
52  ! Constants:                                                                 *
53  ! maxwf                maximum number of wind fields                         *
54  ! unitavailab          unit connected to file AVAILABLE                      *
55  !                                                                            *
56  !*****************************************************************************
57
58  use par_mod
59  use com_mod
60
61  implicit none
62
63  integer :: i,idiff,ldat,ltim,wftime1(maxwf),numbwfn(maxnests),k
64  integer :: wftime1n(maxnests,maxwf),wftimen(maxnests,maxwf)
65  real(kind=dp) :: juldate,jul,beg,end
66  character(len=255) :: fname,spec,wfname1(maxwf),wfspec1(maxwf)
67  character(len=255) :: wfname1n(maxnests,maxwf)
68  character(len=40) :: wfspec1n(maxnests,maxwf)
69
70
71  ! Windfields are only used, if they are within the modelling period.
72  ! However, 1 additional day at the beginning and at the end is used for
73  ! interpolation. -> Compute beginning and ending date for the windfields.
74  !************************************************************************
75
76  if (ideltas.gt.0) then         ! forward trajectories
77    beg=bdate-1._dp
78    end=bdate+real(ideltas,kind=dp)/86400._dp+real(idiffmax,kind=dp)/ &
79         86400._dp
80  else                           ! backward trajectories
81    beg=bdate+real(ideltas,kind=dp)/86400._dp-real(idiffmax,kind=dp)/ &
82         86400._dp
83    end=bdate+1._dp
84  endif
85
86  ! Open the wind field availability file and read available wind fields
87  ! within the modelling period.
88  !*********************************************************************
89
90  open(unitavailab,file=path(4)(1:length(4)),status='old', &
91       err=999)
92
93  do i=1,3
94    read(unitavailab,*)
95  end do
96
97  numbwf=0
98100   read(unitavailab,'(i8,1x,i6,2(6x,a255))',end=99) &
99           ldat,ltim,fname,spec
100    jul=juldate(ldat,ltim)
101    if ((jul.ge.beg).and.(jul.le.end)) then
102      numbwf=numbwf+1
103      if (numbwf.gt.maxwf) then      ! check exceedance of dimension
104       write(*,*) 'Number of wind fields needed is too great.'
105       write(*,*) 'Reduce modelling period (file "COMMAND") or'
106       write(*,*) 'reduce number of wind fields (file "AVAILABLE").'
107       stop
108      endif
109
110      wfname1(numbwf)=fname(1:index(fname,' '))
111      wfspec1(numbwf)=spec
112      wftime1(numbwf)=nint((jul-bdate)*86400._dp)
113    endif
114    goto 100       ! next wind field
115
11699   continue
117
118  close(unitavailab)
119
120  ! Open the wind field availability file and read available wind fields
121  ! within the modelling period (nested grids)
122  !*********************************************************************
123
124  do k=1,numbnests
125  !print*,length(numpath+2*(k-1)+1),length(numpath+2*(k-1)+2),length(4),length(3)
126  !print*,path(numpath+2*(k-1)+2)(1:length(numpath+2*(k-1)+2))
127    open(unitavailab,file=path(numpath+2*(k-1)+2) &
128         (1:length(numpath+2*(k-1)+2)),status='old',err=998)
129
130    do i=1,3
131      read(unitavailab,*)
132    end do
133
134    numbwfn(k)=0
135700   read(unitavailab,'(i8,1x,i6,2(6x,a255))',end=699) ldat, &
136           ltim,fname,spec
137      jul=juldate(ldat,ltim)
138      if ((jul.ge.beg).and.(jul.le.end)) then
139        numbwfn(k)=numbwfn(k)+1
140        if (numbwfn(k).gt.maxwf) then      ! check exceedance of dimension
141       write(*,*) 'Number of nested wind fields is too great.'
142       write(*,*) 'Reduce modelling period (file "COMMAND") or'
143       write(*,*) 'reduce number of wind fields (file "AVAILABLE").'
144          stop
145        endif
146
147        wfname1n(k,numbwfn(k))=fname
148        wfspec1n(k,numbwfn(k))=spec
149        wftime1n(k,numbwfn(k))=nint((jul-bdate)*86400._dp)
150      endif
151      goto 700       ! next wind field
152
153699   continue
154
155    close(unitavailab)
156  end do
157
158
159  ! Check wind field times of file AVAILABLE (expected to be in temporal order)
160  !****************************************************************************
161
162  if (numbwf.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,numbwf
169    if (wftime1(i).le.wftime1(i-1)) then
170      write(*,*) 'FLEXPART ERROR: FILE AVAILABLE IS CORRUPT.'
171      write(*,*) 'THE WIND FIELDS ARE NOT IN TEMPORAL ORDER.'
172      write(*,*) 'PLEASE CHECK FIELD ',wfname1(i)
173      stop
174    endif
175  end do
176
177  ! Check wind field times of file AVAILABLE for the nested fields
178  ! (expected to be in temporal order)
179  !***************************************************************
180
181  do k=1,numbnests
182    if (numbwfn(k).eq.0) then
183      write(*,*) '#### FLEXPART MODEL ERROR! NO WIND FIELDS  ####'
184      write(*,*) '#### AVAILABLE FOR SELECTED TIME PERIOD.   ####'
185      stop
186    endif
187
188    do i=2,numbwfn(k)
189      if (wftime1n(k,i).le.wftime1n(k,i-1)) then
190      write(*,*) 'FLEXPART ERROR: FILE AVAILABLE IS CORRUPT. '
191      write(*,*) 'THE NESTED WIND FIELDS ARE NOT IN TEMPORAL ORDER.'
192      write(*,*) 'PLEASE CHECK FIELD ',wfname1n(k,i)
193      write(*,*) 'AT NESTING LEVEL ',k
194      stop
195      endif
196    end do
197
198  end do
199
200
201  ! For backward trajectories, reverse the order of the windfields
202  !***************************************************************
203
204  if (ideltas.ge.0) then
205    do i=1,numbwf
206      wfname(i)=wfname1(i)
207      wfspec(i)=wfspec1(i)
208      wftime(i)=wftime1(i)
209    end do
210    do k=1,numbnests
211      do i=1,numbwfn(k)
212        wfnamen(k,i)=wfname1n(k,i)
213        wfspecn(k,i)=wfspec1n(k,i)
214        wftimen(k,i)=wftime1n(k,i)
215      end do
216    end do
217  else
218    do i=1,numbwf
219      wfname(numbwf-i+1)=wfname1(i)
220      wfspec(numbwf-i+1)=wfspec1(i)
221      wftime(numbwf-i+1)=wftime1(i)
222    end do
223    do k=1,numbnests
224      do i=1,numbwfn(k)
225        wfnamen(k,numbwfn(k)-i+1)=wfname1n(k,i)
226        wfspecn(k,numbwfn(k)-i+1)=wfspec1n(k,i)
227        wftimen(k,numbwfn(k)-i+1)=wftime1n(k,i)
228      end do
229    end do
230  endif
231
232  ! Check the time difference between the wind fields. If it is big,
233  ! write a warning message. If it is too big, terminate the trajectory.
234  !*********************************************************************
235
236  do i=2,numbwf
237    idiff=abs(wftime(i)-wftime(i-1))
238    if (idiff.gt.idiffmax) then
239      write(*,*) 'FLEXPART WARNING: TIME DIFFERENCE BETWEEN TWO'
240      write(*,*) 'WIND FIELDS IS TOO BIG FOR TRANSPORT CALCULATION.&
241           &'
242      write(*,*) 'THEREFORE, TRAJECTORIES HAVE TO BE SKIPPED.'
243    else if (idiff.gt.idiffnorm) then
244      write(*,*) 'FLEXPART WARNING: TIME DIFFERENCE BETWEEN TWO'
245      write(*,*) 'WIND FIELDS IS BIG. THIS MAY CAUSE A DEGRADATION'
246      write(*,*) 'OF SIMULATION QUALITY.'
247    endif
248  end do
249
250  do k=1,numbnests
251    if (numbwfn(k).ne.numbwf) then
252      write(*,*) 'FLEXPART ERROR: THE AVAILABLE FILES FOR THE'
253      write(*,*) 'NESTED WIND FIELDS ARE NOT CONSISTENT WITH'
254      write(*,*) 'THE AVAILABLE FILE OF THE MOTHER DOMAIN.  '
255      write(*,*) 'ERROR AT NEST LEVEL: ',k
256      stop
257    endif
258    do i=1,numbwf
259      if (wftimen(k,i).ne.wftime(i)) then
260        write(*,*) 'FLEXPART ERROR: THE AVAILABLE FILES FOR THE'
261        write(*,*) 'NESTED WIND FIELDS ARE NOT CONSISTENT WITH'
262        write(*,*) 'THE AVAILABLE FILE OF THE MOTHER DOMAIN.  '
263        write(*,*) 'ERROR AT NEST LEVEL: ',k
264        stop
265      endif
266    end do
267  end do
268
269  ! Reset the times of the wind fields that are kept in memory to no time
270  !**********************************************************************
271
272  do i=1,2
273    memind(i)=i
274    memtime(i)=999999999
275  end do
276
277  return
278
279998   write(*,*) ' #### FLEXPART MODEL ERROR! AVAILABLE FILE   #### '
280  write(*,'(a)') '     '//path(numpath+2*(k-1)+2) &
281       (1:length(numpath+2*(k-1)+2))
282  write(*,*) ' #### CANNOT BE OPENED             #### '
283  stop
284
285999   write(*,*) ' #### FLEXPART MODEL ERROR! AVAILABLE IILE #### '
286  write(*,'(a)') '     '//path(4)(1:length(4))
287  write(*,*) ' #### CANNOT BE OPENED           #### '
288  stop
289
290end subroutine readavailable
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG