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 | |
---|
22 | subroutine getfields(itime,nstop,metdata_format) |
---|
23 | ! i o |
---|
24 | !***************************************************************************** |
---|
25 | ! * |
---|
26 | ! This subroutine manages the 3 data fields to be kept in memory. * |
---|
27 | ! During the first time step of petterssen it has to be fulfilled that the * |
---|
28 | ! first data field must have |wftime|<itime, i.e. the absolute value of * |
---|
29 | ! wftime must be smaller than the absolute value of the current time in [s].* |
---|
30 | ! The other 2 fields are the next in time after the first one. * |
---|
31 | ! Pointers (memind) are used, because otherwise one would have to resort the* |
---|
32 | ! wind fields, which costs a lot of computing time. Here only the pointers * |
---|
33 | ! are resorted. * |
---|
34 | ! * |
---|
35 | ! Author: A. Stohl * |
---|
36 | ! * |
---|
37 | ! 29 April 1994 * |
---|
38 | ! * |
---|
39 | !***************************************************************************** |
---|
40 | ! Changes, Bernd C. Krueger, Feb. 2001: * |
---|
41 | ! Variables tth,qvh,tthn,qvhn (on eta coordinates) in common block. * |
---|
42 | ! Function of nstop extended. * |
---|
43 | ! * |
---|
44 | ! Unified ECMWF and GFS builds * |
---|
45 | ! Marian Harustak, 12.5.2017 * |
---|
46 | ! - Added passing of metdata_format as it was needed by called routines * |
---|
47 | !***************************************************************************** |
---|
48 | ! * |
---|
49 | ! Variables: * |
---|
50 | ! lwindinterval [s] time difference between the two wind fields read in * |
---|
51 | ! indj indicates the number of the wind field to be read in * |
---|
52 | ! indmin remembers the number of wind fields already treated * |
---|
53 | ! memind(2) pointer, on which place the wind fields are stored * |
---|
54 | ! memtime(2) [s] times of the wind fields, which are kept in memory * |
---|
55 | ! itime [s] current time since start date of trajectory calcu- * |
---|
56 | ! lation * |
---|
57 | ! nstop > 0, if trajectory has to be terminated * |
---|
58 | ! nx,ny,nuvz,nwz field dimensions in x,y and z direction * |
---|
59 | ! uu(0:nxmax,0:nymax,nuvzmax,2) wind components in x-direction [m/s] * |
---|
60 | ! vv(0:nxmax,0:nymax,nuvzmax,2) wind components in y-direction [m/s] * |
---|
61 | ! ww(0:nxmax,0:nymax,nwzmax,2) wind components in z-direction [deltaeta/s]* |
---|
62 | ! tt(0:nxmax,0:nymax,nuvzmax,2) temperature [K] * |
---|
63 | ! ps(0:nxmax,0:nymax,2) surface pressure [Pa] * |
---|
64 | ! metdata_format format of metdata (ecmwf/gfs) * |
---|
65 | ! * |
---|
66 | ! Constants: * |
---|
67 | ! idiffmax maximum allowable time difference between 2 wind * |
---|
68 | ! fields * |
---|
69 | ! * |
---|
70 | !***************************************************************************** |
---|
71 | |
---|
72 | use par_mod |
---|
73 | use com_mod |
---|
74 | use class_gribfile |
---|
75 | |
---|
76 | implicit none |
---|
77 | |
---|
78 | integer :: indj,itime,nstop,memaux |
---|
79 | integer :: metdata_format |
---|
80 | |
---|
81 | real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
82 | real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
83 | real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
84 | real :: wwh(0:nxmax-1,0:nymax-1,nwzmax) |
---|
85 | real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
86 | real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
87 | real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests) |
---|
88 | real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests) |
---|
89 | |
---|
90 | integer :: indmin = 1 |
---|
91 | |
---|
92 | |
---|
93 | ! Check, if wind fields are available for the current time step |
---|
94 | !************************************************************** |
---|
95 | |
---|
96 | nstop=0 |
---|
97 | |
---|
98 | if ((ldirect*wftime(1).gt.ldirect*itime).or. & |
---|
99 | (ldirect*wftime(numbwf).lt.ldirect*itime)) then |
---|
100 | write(*,*) 'FLEXPART WARNING: NO WIND FIELDS ARE AVAILABLE.' |
---|
101 | write(*,*) 'A TRAJECTORY HAS TO BE TERMINATED.' |
---|
102 | nstop=4 |
---|
103 | return |
---|
104 | endif |
---|
105 | |
---|
106 | |
---|
107 | if ((ldirect*memtime(1).le.ldirect*itime).and. & |
---|
108 | (ldirect*memtime(2).gt.ldirect*itime)) then |
---|
109 | |
---|
110 | ! The right wind fields are already in memory -> don't do anything |
---|
111 | !***************************************************************** |
---|
112 | |
---|
113 | continue |
---|
114 | |
---|
115 | else if ((ldirect*memtime(2).le.ldirect*itime).and. & |
---|
116 | (memtime(2).ne.999999999)) then |
---|
117 | |
---|
118 | |
---|
119 | ! Current time is after 2nd wind field |
---|
120 | ! -> Resort wind field pointers, so that current time is between 1st and 2nd |
---|
121 | !*************************************************************************** |
---|
122 | |
---|
123 | memaux=memind(1) |
---|
124 | memind(1)=memind(2) |
---|
125 | memind(2)=memaux |
---|
126 | memtime(1)=memtime(2) |
---|
127 | |
---|
128 | |
---|
129 | ! Read a new wind field and store it on place memind(2) |
---|
130 | !****************************************************** |
---|
131 | |
---|
132 | do indj=indmin,numbwf-1 |
---|
133 | if (ldirect*wftime(indj+1).gt.ldirect*itime) then |
---|
134 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
135 | call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh) |
---|
136 | else |
---|
137 | call readwind_gfs(indj+1,memind(2),uuh,vvh,wwh) |
---|
138 | end if |
---|
139 | call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn) |
---|
140 | call calcpar(memind(2),uuh,vvh,pvh,metdata_format) |
---|
141 | call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format) |
---|
142 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
143 | call verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh) |
---|
144 | else |
---|
145 | call verttransform_gfs(memind(2),uuh,vvh,wwh,pvh) |
---|
146 | end if |
---|
147 | call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn) |
---|
148 | memtime(2)=wftime(indj+1) |
---|
149 | nstop = 1 |
---|
150 | goto 40 |
---|
151 | endif |
---|
152 | end do |
---|
153 | 40 indmin=indj |
---|
154 | |
---|
155 | if (WETBKDEP) then |
---|
156 | call writeprecip(itime,memind(1)) |
---|
157 | endif |
---|
158 | |
---|
159 | else |
---|
160 | |
---|
161 | ! No wind fields, which can be used, are currently in memory |
---|
162 | ! -> read both wind fields |
---|
163 | !*********************************************************** |
---|
164 | |
---|
165 | do indj=indmin,numbwf-1 |
---|
166 | if ((ldirect*wftime(indj).le.ldirect*itime).and. & |
---|
167 | (ldirect*wftime(indj+1).gt.ldirect*itime)) then |
---|
168 | memind(1)=1 |
---|
169 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
170 | call readwind_ecmwf(indj,memind(1),uuh,vvh,wwh) |
---|
171 | else |
---|
172 | call readwind_gfs(indj,memind(1),uuh,vvh,wwh) |
---|
173 | end if |
---|
174 | call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn) |
---|
175 | call calcpar(memind(1),uuh,vvh,pvh,metdata_format) |
---|
176 | call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format) |
---|
177 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
178 | call verttransform_ecmwf(memind(1),uuh,vvh,wwh,pvh) |
---|
179 | else |
---|
180 | call verttransform_gfs(memind(1),uuh,vvh,wwh,pvh) |
---|
181 | end if |
---|
182 | call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn) |
---|
183 | memtime(1)=wftime(indj) |
---|
184 | memind(2)=2 |
---|
185 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
186 | call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh) |
---|
187 | else |
---|
188 | call readwind_gfs(indj+1,memind(2),uuh,vvh,wwh) |
---|
189 | end if |
---|
190 | call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn) |
---|
191 | call calcpar(memind(2),uuh,vvh,pvh,metdata_format) |
---|
192 | call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format) |
---|
193 | if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then |
---|
194 | call verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh) |
---|
195 | else |
---|
196 | call verttransform_gfs(memind(2),uuh,vvh,wwh,pvh) |
---|
197 | end if |
---|
198 | call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn) |
---|
199 | memtime(2)=wftime(indj+1) |
---|
200 | nstop = 1 |
---|
201 | goto 60 |
---|
202 | endif |
---|
203 | end do |
---|
204 | 60 indmin=indj |
---|
205 | |
---|
206 | if (WETBKDEP) then |
---|
207 | call writeprecip(itime,memind(1)) |
---|
208 | endif |
---|
209 | |
---|
210 | endif |
---|
211 | |
---|
212 | lwindinterv=abs(memtime(2)-memtime(1)) |
---|
213 | |
---|
214 | if (lwindinterv.gt.idiffmax) nstop=3 |
---|
215 | |
---|
216 | end subroutine getfields |
---|