CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
getfpfields.F90
Go to the documentation of this file.
1 
2 subroutine getfpfields(itime,nstop,metdata_format)
3  ! i o
4  !*****************************************************************************
5  ! *
6  ! This subroutine manages the 3 data fields to be kept in memory. *
7  ! During the first time step of petterssen it has to be fulfilled that the *
8  ! first data field must have |wftime|<itime, i.e. the absolute value of *
9  ! wftime must be smaller than the absolute value of the current time in [s].*
10  ! The other 2 fields are the next in time after the first one. *
11  ! Pointers (memind) are used, because otherwise one would have to resort the*
12  ! wind fields, which costs a lot of computing time. Here only the pointers *
13  ! are resorted. *
14  ! *
15  ! Author: A. Stohl *
16  ! *
17  ! 29 April 1994 *
18  ! *
19  !*****************************************************************************
20  ! Changes, Bernd C. Krueger, Feb. 2001:
21  ! Variables tth,qvh,tthn,qvhn (on eta coordinates) in common block.
22  ! Function of nstop extended.
23  !*****************************************************************************
24  ! *
25  ! Variables: *
26  ! lwindinterval [s] time difference between the two wind fields read in *
27  ! indj indicates the number of the wind field to be read in *
28  ! indmin remembers the number of wind fields already treated *
29  ! memind(2) pointer, on which place the wind fields are stored *
30  ! memtime(2) [s] times of the wind fields, which are kept in memory *
31  ! itime [s] current time since start date of trajectory calcu- *
32  ! lation *
33  ! nstop > 0, if trajectory has to be terminated *
34  ! nx,ny,nuvz,nwz field dimensions in x,y and z direction *
35  ! uu(0:nxmax,0:nymax,nuvzmax,2) wind components in x-direction [m/s] *
36  ! vv(0:nxmax,0:nymax,nuvzmax,2) wind components in y-direction [m/s] *
37  ! ww(0:nxmax,0:nymax,nwzmax,2) wind components in z-direction [deltaeta/s]*
38  ! tt(0:nxmax,0:nymax,nuvzmax,2) temperature [K] *
39  ! ps(0:nxmax,0:nymax,2) surface pressure [Pa] *
40  ! *
41  ! Constants: *
42  ! idiffmax maximum allowable time difference between 2 wind *
43  ! fields *
44  ! *
45  !*****************************************************************************
46 
47  use par_mod
48  use com_mod
49 
50  use fpmetbinary_mod
51 
52  implicit none
53 
54  integer :: indj,itime,nstop,memaux,metdata_format
55 
56  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
57  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
58  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
59  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
60  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
61  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
62  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
63  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
64 
65  integer :: indmin = 1
66 
67  !character(len=120) :: fpdir ! Directory the .fp files are in
68  character(len=512) fpfname ! .fp filename
69 
70 #ifdef PERFTIMER
71  INTEGER millisecs_start, millisecs_stop, count_rate, count_max
72 #endif
73 
74 
75  ! Check, if wind fields are available for the current time step
76  !**************************************************************
77 
78  nstop=0
79 
80  if ((ldirect*wftime(1).gt.ldirect*itime).or. &
81  (ldirect*wftime(numbwf).lt.ldirect*itime)) then
82  write(*,*) 'FLEXPART WARNING: NO WIND FIELDS ARE AVAILABLE.'
83  write(*,*) 'A TRAJECTORY HAS TO BE TERMINATED.'
84  nstop=4
85  return
86  endif
87 
88 
89  if ((ldirect*memtime(1).le.ldirect*itime).and. &
90  (ldirect*memtime(2).gt.ldirect*itime)) then
91 
92  ! The right wind fields are already in memory -> don't do anything
93  !*****************************************************************
94 
95  continue
96 
97  else if ((ldirect*memtime(2).le.ldirect*itime).and. &
98  (memtime(2).ne.999999999)) then
99 
100 
101  ! Current time is after 2nd wind field
102  ! -> Resort wind field pointers, so that current time is between 1st and 2nd
103  !***************************************************************************
104 
105  memaux=memind(1)
106  memind(1)=memind(2)
107  memind(2)=memaux
108  memtime(1)=memtime(2)
109 
110 
111  ! Read a new wind field and store it on place memind(2)
112  !******************************************************
113 
114  do indj=indmin,numbwf-1
115  if (ldirect*wftime(indj+1).gt.ldirect*itime) then
116 
117 #ifdef PERFTIMER
118  CALL system_clock(millisecs_start, count_rate, count_max)
119 #endif
120  ! Read in a single .fp file, placing contents at index 2
121  if ( ldirect.eq.1 ) then
122  fpfname = trim(path(3)) // trim(wfname(indj+1)) // '_fwd.fp'
123  else
124  fpfname = trim(path(3)) // trim(wfname(indj+1)) // '_bwd.fp'
125  endif
126  print *, 'loading... ', trim(fpfname)
127  CALL fpmetbinary_load(trim(fpfname), memind(2))
128  memtime(2)=wftime(indj+1)
129  nstop = 1
130 
131 #ifdef PERFTIMER
132  CALL system_clock(millisecs_stop, count_rate, count_max)
133  print *, 'Wall time to process: ', trim(fpfname), &
134  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
135 #endif
136 
137 
138  goto 40
139 
140  endif
141  end do
142  40 indmin=indj
143 
144  else
145 
146  ! No wind fields, which can be used, are currently in memory
147  ! -> read both wind fields
148  !***********************************************************
149 
150  do indj=indmin,numbwf-1
151  if ((ldirect*wftime(indj).le.ldirect*itime).and. &
152  (ldirect*wftime(indj+1).gt.ldirect*itime)) then
153  memind(1)=1
154 
155 #ifdef PERFTIMER
156  CALL system_clock(millisecs_start, count_rate, count_max)
157 #endif
158  ! Read in first .fp file
159  if ( ldirect.eq.1 ) then
160  fpfname = trim(path(3)) // trim(wfname(indj)) // '_fwd.fp'
161  else
162  fpfname = trim(path(3)) // trim(wfname(indj)) // '_bwd.fp'
163  endif
164  print *, 'loading... ', trim(fpfname)
165  CALL fpmetbinary_load(trim(fpfname), memind(1))
166  memtime(1)=wftime(indj)
167  memind(2)=2
168 #ifdef PERFTIMER
169  CALL system_clock(millisecs_stop, count_rate, count_max)
170  print *, 'Wall time to process: ', trim(fpfname), &
171  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
172 #endif
173 
174 #ifdef PERFTIMER
175  CALL system_clock(millisecs_start, count_rate, count_max)
176 #endif
177  ! Read in second .fp file
178  if ( ldirect.eq.1 ) then
179  fpfname = trim(path(3)) // trim(wfname(indj+1)) // '_fwd.fp'
180  else
181  fpfname = trim(path(3)) // trim(wfname(indj+1)) // '_bwd.fp'
182  endif
183  print *, 'loading... ', trim(fpfname)
184  CALL fpmetbinary_load(trim(fpfname), memind(2))
185  memtime(2)=wftime(indj+1)
186  nstop = 1
187 #ifdef PERFTIMER
188  CALL system_clock(millisecs_stop, count_rate, count_max)
189  print *, 'Wall time to process: ', trim(fpfname), &
190  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
191 #endif
192 
193 
194  goto 60
195  endif
196  end do
197  60 indmin=indj
198 
199  endif
200 
201  lwindinterv=abs(memtime(2)-memtime(1))
202 
203  if (lwindinterv.gt.idiffmax) nstop=3
204 
205 end subroutine getfpfields
206 
207 
subroutine fpmetbinary_load(filename, cm_index)
subroutine getfpfields(itime, nstop, metdata_format)
Definition: getfpfields.F90:2