source: branches/flexpart91_hasod/src_parallel/readavailable.f90 @ 10

Last change on this file since 10 was 10, checked in by hasod, 11 years ago

ADD: namelist input implemented for all common input files

  • Property svn:executable set to *
File size: 11.0 KB
RevLine 
[8]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  ! wftime(maxwf) [s]times of wind fields relative to beginning time           *
[10]49  ! wfname1,wftime1 = same as above, but only local (help variables)   *
[8]50  !                                                                            *
51  ! Constants:                                                                 *
52  ! maxwf                maximum number of wind fields                         *
53  ! unitavailab          unit connected to file AVAILABLE                      *
54  !                                                                            *
55  !*****************************************************************************
56
57  use par_mod
58  use com_mod
59
60  implicit none
61
62  integer :: i,idiff,ldat,ltim,wftime1(maxwf),numbwfn(maxnests),k
63  integer :: wftime1n(maxnests,maxwf),wftimen(maxnests,maxwf)
64  real(kind=dp) :: juldate,jul,beg,end
[10]65  character(len=255) :: fname,spec,wfname1(maxwf)
[8]66  character(len=255) :: wfname1n(maxnests,maxwf)
67
68
69  ! Windfields are only used, if they are within the modelling period.
70  ! However, 1 additional day at the beginning and at the end is used for
71  ! interpolation. -> Compute beginning and ending date for the windfields.
72  !************************************************************************
73
74  if (ideltas.gt.0) then         ! forward trajectories
75    beg=bdate-1._dp
76    end=bdate+real(ideltas,kind=dp)/86400._dp+real(idiffmax,kind=dp)/ &
77         86400._dp
78  else                           ! backward trajectories
79    beg=bdate+real(ideltas,kind=dp)/86400._dp-real(idiffmax,kind=dp)/ &
80         86400._dp
81    end=bdate+1._dp
82  endif
83
84  ! Open the wind field availability file and read available wind fields
85  ! within the modelling period.
86  !*********************************************************************
87
88  open(unitavailab,file=path(4)(1:length(4)),status='old', &
89       err=999)
90
91  do i=1,3
92    read(unitavailab,*)
93  end do
94
95  numbwf=0
96100   read(unitavailab,'(i8,1x,i6,2(6x,a255))',end=99) &
97           ldat,ltim,fname,spec
98    jul=juldate(ldat,ltim)
99    if ((jul.ge.beg).and.(jul.le.end)) then
100      numbwf=numbwf+1
101      if (numbwf.gt.maxwf) then      ! check exceedance of dimension
102       write(*,*) 'Number of wind fields needed is too great.'
103       write(*,*) 'Reduce modelling period (file "COMMAND") or'
104       write(*,*) 'reduce number of wind fields (file "AVAILABLE").'
105       stop
106      endif
107
108      wfname1(numbwf)=fname(1:index(fname,' '))
109      wftime1(numbwf)=nint((jul-bdate)*86400._dp)
110    endif
111    goto 100       ! next wind field
112
11399   continue
114
115  close(unitavailab)
116
117  ! Open the wind field availability file and read available wind fields
118  ! within the modelling period (nested grids)
119  !*********************************************************************
120
121  do k=1,numbnests
122  print*,length(numpath+2*(k-1)+1),length(numpath+2*(k-1)+2),length(4),length(3)
123  print*,path(numpath+2*(k-1)+2)(1:length(numpath+2*(k-1)+2))
124    open(unitavailab,file=path(numpath+2*(k-1)+2) &
125         (1:length(numpath+2*(k-1)+2)),status='old',err=998)
126
127    do i=1,3
128      read(unitavailab,*)
129    end do
130
131    numbwfn(k)=0
132700   read(unitavailab,'(i8,1x,i6,2(6x,a255))',end=699) ldat, &
133           ltim,fname,spec
134      jul=juldate(ldat,ltim)
135      if ((jul.ge.beg).and.(jul.le.end)) then
136        numbwfn(k)=numbwfn(k)+1
137        if (numbwfn(k).gt.maxwf) then      ! check exceedance of dimension
138       write(*,*) 'Number of nested wind fields is too great.'
139       write(*,*) 'Reduce modelling period (file "COMMAND") or'
140       write(*,*) 'reduce number of wind fields (file "AVAILABLE").'
141          stop
142        endif
143
144        wfname1n(k,numbwfn(k))=fname
145        wftime1n(k,numbwfn(k))=nint((jul-bdate)*86400._dp)
146      endif
147      goto 700       ! next wind field
148
149699   continue
150
151    close(unitavailab)
152  end do
153
154
155  ! Check wind field times of file AVAILABLE (expected to be in temporal order)
156  !****************************************************************************
157
158  if (numbwf.eq.0) then
159    write(*,*) ' #### FLEXPART MODEL ERROR! NO WIND FIELDS    #### '
160    write(*,*) ' #### AVAILABLE FOR SELECTED TIME PERIOD.     #### '
161    stop
162  endif
163
164  do i=2,numbwf
165    if (wftime1(i).le.wftime1(i-1)) then
166      write(*,*) 'FLEXPART ERROR: FILE AVAILABLE IS CORRUPT.'
167      write(*,*) 'THE WIND FIELDS ARE NOT IN TEMPORAL ORDER.'
168      write(*,*) 'PLEASE CHECK FIELD ',wfname1(i)
169      stop
170    endif
171  end do
172
173  ! Check wind field times of file AVAILABLE for the nested fields
174  ! (expected to be in temporal order)
175  !***************************************************************
176
177  do k=1,numbnests
178    if (numbwfn(k).eq.0) then
179      write(*,*) '#### FLEXPART MODEL ERROR! NO WIND FIELDS  ####'
180      write(*,*) '#### AVAILABLE FOR SELECTED TIME PERIOD.   ####'
181      stop
182    endif
183
184    do i=2,numbwfn(k)
185      if (wftime1n(k,i).le.wftime1n(k,i-1)) then
186      write(*,*) 'FLEXPART ERROR: FILE AVAILABLE IS CORRUPT. '
187      write(*,*) 'THE NESTED WIND FIELDS ARE NOT IN TEMPORAL ORDER.'
188      write(*,*) 'PLEASE CHECK FIELD ',wfname1n(k,i)
189      write(*,*) 'AT NESTING LEVEL ',k
190      stop
191      endif
192    end do
193
194  end do
195
196
197  ! For backward trajectories, reverse the order of the windfields
198  !***************************************************************
199
200  if (ideltas.ge.0) then
201    do i=1,numbwf
202      wfname(i)=wfname1(i)
203      wftime(i)=wftime1(i)
204    end do
205    do k=1,numbnests
206      do i=1,numbwfn(k)
207        wfnamen(k,i)=wfname1n(k,i)
208        wftimen(k,i)=wftime1n(k,i)
209      end do
210    end do
211  else
212    do i=1,numbwf
213      wfname(numbwf-i+1)=wfname1(i)
214      wftime(numbwf-i+1)=wftime1(i)
215    end do
216    do k=1,numbnests
217      do i=1,numbwfn(k)
218        wfnamen(k,numbwfn(k)-i+1)=wfname1n(k,i)
219        wftimen(k,numbwfn(k)-i+1)=wftime1n(k,i)
220      end do
221    end do
222  endif
223
224  ! Check the time difference between the wind fields. If it is big,
225  ! write a warning message. If it is too big, terminate the trajectory.
226  !*********************************************************************
227
228  do i=2,numbwf
229    idiff=abs(wftime(i)-wftime(i-1))
230    if (idiff.gt.idiffmax) then
231      write(*,*) 'FLEXPART WARNING: TIME DIFFERENCE BETWEEN TWO'
232      write(*,*) 'WIND FIELDS IS TOO BIG FOR TRANSPORT CALCULATION.&
233           &'
234      write(*,*) 'THEREFORE, TRAJECTORIES HAVE TO BE SKIPPED.'
235    else if (idiff.gt.idiffnorm) then
236      write(*,*) 'FLEXPART WARNING: TIME DIFFERENCE BETWEEN TWO'
237      write(*,*) 'WIND FIELDS IS BIG. THIS MAY CAUSE A DEGRADATION'
238      write(*,*) 'OF SIMULATION QUALITY.'
239    endif
240  end do
241
242  do k=1,numbnests
243    if (numbwfn(k).ne.numbwf) 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    do i=1,numbwf
251      if (wftimen(k,i).ne.wftime(i)) 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    end do
259  end do
260
261  ! Reset the times of the wind fields that are kept in memory to no time
262  !**********************************************************************
263
264  do i=1,2
265    memind(i)=i
266    memtime(i)=999999999
267  end do
268
269  return
270
271998   write(*,*) ' #### FLEXPART MODEL ERROR! AVAILABLE FILE   #### '
272  write(*,'(a)') '     '//path(numpath+2*(k-1)+2) &
273       (1:length(numpath+2*(k-1)+2))
274  write(*,*) ' #### CANNOT BE OPENED             #### '
275  stop
276
277999   write(*,*) ' #### FLEXPART MODEL ERROR! AVAILABLE IILE #### '
278  write(*,'(a)') '     '//path(4)(1:length(4))
279  write(*,*) ' #### CANNOT BE OPENED           #### '
280  stop
281
282end subroutine readavailable
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG