1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
3 | |
---|
4 | subroutine readwind_gfs(indj,n,uuh,vvh,wwh) |
---|
5 | |
---|
6 | !*********************************************************************** |
---|
7 | !* * |
---|
8 | !* TRAJECTORY MODEL SUBROUTINE READWIND * |
---|
9 | !* * |
---|
10 | !*********************************************************************** |
---|
11 | !* * |
---|
12 | !* AUTHOR: G. WOTAWA * |
---|
13 | !* DATE: 1997-08-05 * |
---|
14 | !* LAST UPDATE: 2000-10-17, Andreas Stohl * |
---|
15 | !* CHANGE: 01/02/2001, Bernd C. Krueger, Variables tth and * |
---|
16 | !* qvh (on eta coordinates) in common block * |
---|
17 | !* CHANGE: 16/11/2005, Caroline Forster, GFS data * |
---|
18 | !* CHANGE: 11/01/2008, Harald Sodemann, Input of GRIB1/2 * |
---|
19 | !* data with the ECMWF grib_api library * |
---|
20 | !* CHANGE: 03/12/2008, Harald Sodemann, update to f90 with * |
---|
21 | !* ECMWF grib_api * |
---|
22 | ! * |
---|
23 | ! Unified ECMWF and GFS builds * |
---|
24 | ! Marian Harustak, 12.5.2017 * |
---|
25 | ! - Renamed routine from readwind to readwind_gfs * |
---|
26 | !* * |
---|
27 | !*********************************************************************** |
---|
28 | !* * |
---|
29 | !* DESCRIPTION: * |
---|
30 | !* * |
---|
31 | !* READING OF ECMWF METEOROLOGICAL FIELDS FROM INPUT DATA FILES. THE * |
---|
32 | !* INPUT DATA FILES ARE EXPECTED TO BE AVAILABLE IN GRIB CODE * |
---|
33 | !* * |
---|
34 | !* INPUT: * |
---|
35 | !* indj indicates number of the wind field to be read in * |
---|
36 | !* n temporal index for meteorological fields (1 to 3)* |
---|
37 | !* * |
---|
38 | !* IMPORTANT VARIABLES FROM COMMON BLOCK: * |
---|
39 | !* * |
---|
40 | !* wfname File name of data to be read in * |
---|
41 | !* nx,ny,nuvz,nwz expected field dimensions * |
---|
42 | !* nlev_ec number of vertical levels ecmwf model * |
---|
43 | !* uu,vv,ww wind fields * |
---|
44 | !* tt,qv temperature and specific humidity * |
---|
45 | !* ps surface pressure * |
---|
46 | !* * |
---|
47 | !*********************************************************************** |
---|
48 | |
---|
49 | use eccodes |
---|
50 | use par_mod |
---|
51 | use com_mod |
---|
52 | |
---|
53 | implicit none |
---|
54 | |
---|
55 | !HSO new parameters for grib_api |
---|
56 | integer :: ifile |
---|
57 | integer :: iret |
---|
58 | integer :: igrib |
---|
59 | integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl |
---|
60 | !HSO end edits |
---|
61 | real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
62 | real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
63 | real :: wwh(0:nxmax-1,0:nymax-1,nwzmax) |
---|
64 | integer :: ii,indj,i,j,k,n,levdiff2,ifield,iumax,iwmax |
---|
65 | |
---|
66 | ! NCEP |
---|
67 | integer :: numpt,numpu,numpv,numpw,numprh,numpclwch |
---|
68 | real :: help, temp, ew |
---|
69 | real :: elev |
---|
70 | real :: ulev1(0:nxmax-1,0:nymax-1),vlev1(0:nxmax-1,0:nymax-1) |
---|
71 | real :: tlev1(0:nxmax-1,0:nymax-1) |
---|
72 | real :: qvh2(0:nxmax-1,0:nymax-1) |
---|
73 | |
---|
74 | integer :: i179,i180,i181 |
---|
75 | |
---|
76 | ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING |
---|
77 | !HSO kept isec1, isec2 and zsec4 for consistency with gribex GRIB input |
---|
78 | |
---|
79 | integer :: isec1(8),isec2(3) |
---|
80 | real(kind=4) :: zsec4(jpunp) |
---|
81 | real(kind=4) :: xaux,yaux,xaux0,yaux0 |
---|
82 | real(kind=8) :: xauxin,yauxin |
---|
83 | real,parameter :: eps=1.e-4 |
---|
84 | real(kind=4) :: ewss(0:nxmax-1,0:nymax-1),nsss(0:nxmax-1,0:nymax-1) |
---|
85 | real :: plev1,hlev1,ff10m,fflev1 |
---|
86 | |
---|
87 | logical :: hflswitch,strswitch |
---|
88 | |
---|
89 | !HSO for grib api error messages |
---|
90 | character(len=24) :: gribErrorMsg = 'Error reading grib file' |
---|
91 | character(len=20) :: gribFunction = 'readwind_gfs' |
---|
92 | character(len=20) :: shortname |
---|
93 | |
---|
94 | |
---|
95 | hflswitch=.false. |
---|
96 | strswitch=.false. |
---|
97 | levdiff2=nlev_ec-nwz+1 |
---|
98 | iumax=0 |
---|
99 | iwmax=0 |
---|
100 | |
---|
101 | |
---|
102 | ! OPENING OF DATA FILE (GRIB CODE) |
---|
103 | |
---|
104 | !HSO |
---|
105 | call grib_open_file(ifile,path(3)(1:length(3)) & |
---|
106 | //trim(wfname(indj)),'r',iret) |
---|
107 | if (iret.ne.GRIB_SUCCESS) then |
---|
108 | goto 888 ! ERROR DETECTED |
---|
109 | endif |
---|
110 | !turn on support for multi fields messages |
---|
111 | call grib_multi_support_on |
---|
112 | |
---|
113 | numpt=0 |
---|
114 | numpu=0 |
---|
115 | numpv=0 |
---|
116 | numpw=0 |
---|
117 | numprh=0 |
---|
118 | numpclwch=0 |
---|
119 | ifield=0 |
---|
120 | 10 ifield=ifield+1 |
---|
121 | ! |
---|
122 | ! GET NEXT FIELDS |
---|
123 | ! |
---|
124 | call grib_new_from_file(ifile,igrib,iret) |
---|
125 | if (iret.eq.GRIB_END_OF_FILE) then |
---|
126 | goto 50 ! EOF DETECTED |
---|
127 | elseif (iret.ne.GRIB_SUCCESS) then |
---|
128 | goto 888 ! ERROR DETECTED |
---|
129 | endif |
---|
130 | |
---|
131 | !first see if we read GRIB1 or GRIB2 |
---|
132 | call grib_get_int(igrib,'editionNumber',gribVer,iret) |
---|
133 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
134 | |
---|
135 | if (gribVer.eq.1) then ! GRIB Edition 1 |
---|
136 | |
---|
137 | !read the grib1 identifiers |
---|
138 | call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret) |
---|
139 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
140 | call grib_get_int(igrib,'indicatorOfTypeOfLevel',isec1(7),iret) |
---|
141 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
142 | call grib_get_int(igrib,'level',isec1(8),iret) |
---|
143 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
144 | |
---|
145 | else ! GRIB Edition 2 |
---|
146 | |
---|
147 | !read the grib2 identifiers |
---|
148 | call grib_get_string(igrib,'shortName',shortname,iret) |
---|
149 | |
---|
150 | call grib_get_int(igrib,'discipline',discipl,iret) |
---|
151 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
152 | call grib_get_int(igrib,'parameterCategory',parCat,iret) |
---|
153 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
154 | call grib_get_int(igrib,'parameterNumber',parNum,iret) |
---|
155 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
156 | call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret) |
---|
157 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
158 | call grib_get_int(igrib,'scaledValueOfFirstFixedSurface', & |
---|
159 | valSurf,iret) |
---|
160 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
161 | |
---|
162 | ! write(*,*) 'Field: ',ifield,parCat,parNum,typSurf,shortname |
---|
163 | !convert to grib1 identifiers |
---|
164 | isec1(6)=-1 |
---|
165 | isec1(7)=-1 |
---|
166 | isec1(8)=-1 |
---|
167 | if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.100)) then ! T |
---|
168 | isec1(6)=11 ! indicatorOfParameter |
---|
169 | isec1(7)=100 ! indicatorOfTypeOfLevel |
---|
170 | isec1(8)=valSurf/100 ! level, convert to hPa |
---|
171 | elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.100)) then ! U |
---|
172 | isec1(6)=33 ! indicatorOfParameter |
---|
173 | isec1(7)=100 ! indicatorOfTypeOfLevel |
---|
174 | isec1(8)=valSurf/100 ! level, convert to hPa |
---|
175 | elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.100)) then ! V |
---|
176 | isec1(6)=34 ! indicatorOfParameter |
---|
177 | isec1(7)=100 ! indicatorOfTypeOfLevel |
---|
178 | isec1(8)=valSurf/100 ! level, convert to hPa |
---|
179 | elseif ((parCat.eq.2).and.(parNum.eq.8).and.(typSurf.eq.100)) then ! W |
---|
180 | isec1(6)=39 ! indicatorOfParameter |
---|
181 | isec1(7)=100 ! indicatorOfTypeOfLevel |
---|
182 | isec1(8)=valSurf/100 ! level, convert to hPa |
---|
183 | elseif ((parCat.eq.1).and.(parNum.eq.1).and.(typSurf.eq.100)) then ! RH |
---|
184 | isec1(6)=52 ! indicatorOfParameter |
---|
185 | isec1(7)=100 ! indicatorOfTypeOfLevel |
---|
186 | isec1(8)=valSurf/100 ! level, convert to hPa |
---|
187 | elseif ((parCat.eq.1).and.(parNum.eq.1).and.(typSurf.eq.103)) then ! RH2 |
---|
188 | isec1(6)=52 ! indicatorOfParameter |
---|
189 | isec1(7)=105 ! indicatorOfTypeOfLevel |
---|
190 | isec1(8)=2 |
---|
191 | elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! T2 |
---|
192 | isec1(6)=11 ! indicatorOfParameter |
---|
193 | isec1(7)=105 ! indicatorOfTypeOfLevel |
---|
194 | isec1(8)=2 |
---|
195 | elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! U10 |
---|
196 | isec1(6)=33 ! indicatorOfParameter |
---|
197 | isec1(7)=105 ! indicatorOfTypeOfLevel |
---|
198 | isec1(8)=10 |
---|
199 | elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! V10 |
---|
200 | isec1(6)=34 ! indicatorOfParameter |
---|
201 | isec1(7)=105 ! indicatorOfTypeOfLevel |
---|
202 | isec1(8)=10 |
---|
203 | elseif ((parCat.eq.1).and.(parNum.eq.22).and.(typSurf.eq.100)) then ! CLWMR Cloud Mixing Ratio [kg/kg]: |
---|
204 | isec1(6)=153 ! indicatorOfParameter |
---|
205 | isec1(7)=100 ! indicatorOfTypeOfLevel |
---|
206 | isec1(8)=valSurf/100 ! level, convert to hPa |
---|
207 | elseif ((parCat.eq.3).and.(parNum.eq.1).and.(typSurf.eq.101)) then ! SLP |
---|
208 | isec1(6)=2 ! indicatorOfParameter |
---|
209 | isec1(7)=102 ! indicatorOfTypeOfLevel |
---|
210 | isec1(8)=0 |
---|
211 | elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then ! SP |
---|
212 | isec1(6)=1 ! indicatorOfParameter |
---|
213 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
214 | isec1(8)=0 |
---|
215 | elseif ((parCat.eq.1).and.(parNum.eq.13).and.(typSurf.eq.1)) then ! SNOW |
---|
216 | isec1(6)=66 ! indicatorOfParameter |
---|
217 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
218 | isec1(8)=0 |
---|
219 | elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.104)) then ! T sigma 0 |
---|
220 | isec1(6)=11 ! indicatorOfParameter |
---|
221 | isec1(7)=107 ! indicatorOfTypeOfLevel |
---|
222 | isec1(8)=0.995 ! lowest sigma level |
---|
223 | elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.104)) then ! U sigma 0 |
---|
224 | isec1(6)=33 ! indicatorOfParameter |
---|
225 | isec1(7)=107 ! indicatorOfTypeOfLevel |
---|
226 | isec1(8)=0.995 ! lowest sigma level |
---|
227 | elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.104)) then ! V sigma 0 |
---|
228 | isec1(6)=34 ! indicatorOfParameter |
---|
229 | isec1(7)=107 ! indicatorOfTypeOfLevel |
---|
230 | isec1(8)=0.995 ! lowest sigma level |
---|
231 | elseif ((parCat.eq.3).and.(parNum.eq.5).and.(typSurf.eq.1)) then ! TOPO |
---|
232 | isec1(6)=7 ! indicatorOfParameter |
---|
233 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
234 | isec1(8)=0 |
---|
235 | elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.1) & |
---|
236 | .and.(discipl.eq.2)) then ! LSM |
---|
237 | isec1(6)=81 ! indicatorOfParameter |
---|
238 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
239 | isec1(8)=0 |
---|
240 | elseif ((parCat.eq.3).and.(parNum.eq.196).and.(typSurf.eq.1)) then ! BLH |
---|
241 | isec1(6)=221 ! indicatorOfParameter |
---|
242 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
243 | isec1(8)=0 |
---|
244 | elseif ((parCat.eq.1).and.(parNum.eq.7).and.(typSurf.eq.1)) then ! LSP/TP |
---|
245 | isec1(6)=62 ! indicatorOfParameter |
---|
246 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
247 | isec1(8)=0 |
---|
248 | elseif ((parCat.eq.1).and.(parNum.eq.196).and.(typSurf.eq.1)) then ! CP |
---|
249 | isec1(6)=63 ! indicatorOfParameter |
---|
250 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
251 | isec1(8)=0 |
---|
252 | endif |
---|
253 | |
---|
254 | endif ! gribVer |
---|
255 | |
---|
256 | if (isec1(6).ne.-1) then |
---|
257 | ! get the size and data of the values array |
---|
258 | call grib_get_real4_array(igrib,'values',zsec4,iret) |
---|
259 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
260 | endif |
---|
261 | |
---|
262 | if(ifield.eq.1) then |
---|
263 | |
---|
264 | !get the required fields from section 2 |
---|
265 | !store compatible to gribex input |
---|
266 | call grib_get_int(igrib,'numberOfPointsAlongAParallel', & |
---|
267 | isec2(2),iret) |
---|
268 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
269 | call grib_get_int(igrib,'numberOfPointsAlongAMeridian', & |
---|
270 | isec2(3),iret) |
---|
271 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
272 | call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', & |
---|
273 | xauxin,iret) |
---|
274 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
275 | call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', & |
---|
276 | yauxin,iret) |
---|
277 | ! call grib_check(iret,gribFunction,gribErrorMsg) |
---|
278 | xaux=xauxin+real(nxshift)*dx |
---|
279 | yaux=yauxin |
---|
280 | |
---|
281 | ! CHECK GRID SPECIFICATIONS |
---|
282 | |
---|
283 | if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT' |
---|
284 | if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT' |
---|
285 | if(xaux.eq.0.) xaux=-179.0 ! NCEP DATA |
---|
286 | xaux0=xlon0 |
---|
287 | yaux0=ylat0 |
---|
288 | if(xaux.lt.0.) xaux=xaux+360. |
---|
289 | if(yaux.lt.0.) yaux=yaux+360. |
---|
290 | if(xaux0.lt.0.) xaux0=xaux0+360. |
---|
291 | if(yaux0.lt.0.) yaux0=yaux0+360. |
---|
292 | if(abs(xaux-xaux0).gt.eps) & |
---|
293 | stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT' |
---|
294 | if(abs(yaux-yaux0).gt.eps) & |
---|
295 | stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT' |
---|
296 | endif |
---|
297 | !HSO end of edits |
---|
298 | |
---|
299 | i179=nint(179./dx) |
---|
300 | if (dx.lt.0.7) then |
---|
301 | i180=nint(180./dx)+1 ! 0.5 deg data |
---|
302 | else |
---|
303 | i180=nint(179./dx)+1 ! 1 deg data |
---|
304 | endif |
---|
305 | i181=i180+1 |
---|
306 | |
---|
307 | if (isec1(6).ne.-1) then |
---|
308 | |
---|
309 | do j=0,nymin1 |
---|
310 | do i=0,nxfield-1 |
---|
311 | if((isec1(6).eq.011).and.(isec1(7).eq.100)) then |
---|
312 | ! TEMPERATURE |
---|
313 | if((i.eq.0).and.(j.eq.0)) then |
---|
314 | do ii=1,nuvz |
---|
315 | if ((isec1(8)*100.0).eq.akz(ii)) numpt=ii |
---|
316 | end do |
---|
317 | endif |
---|
318 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
319 | if(i.le.i180) then |
---|
320 | tth(i179+i,j,numpt,n)=help |
---|
321 | else |
---|
322 | tth(i-i181,j,numpt,n)=help |
---|
323 | endif |
---|
324 | endif |
---|
325 | if((isec1(6).eq.033).and.(isec1(7).eq.100)) then |
---|
326 | ! U VELOCITY |
---|
327 | if((i.eq.0).and.(j.eq.0)) then |
---|
328 | do ii=1,nuvz |
---|
329 | if ((isec1(8)*100.0).eq.akz(ii)) numpu=ii |
---|
330 | end do |
---|
331 | endif |
---|
332 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
333 | if(i.le.i180) then |
---|
334 | uuh(i179+i,j,numpu)=help |
---|
335 | else |
---|
336 | uuh(i-i181,j,numpu)=help |
---|
337 | endif |
---|
338 | endif |
---|
339 | if((isec1(6).eq.034).and.(isec1(7).eq.100)) then |
---|
340 | ! V VELOCITY |
---|
341 | if((i.eq.0).and.(j.eq.0)) then |
---|
342 | do ii=1,nuvz |
---|
343 | if ((isec1(8)*100.0).eq.akz(ii)) numpv=ii |
---|
344 | end do |
---|
345 | endif |
---|
346 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
347 | if(i.le.i180) then |
---|
348 | vvh(i179+i,j,numpv)=help |
---|
349 | else |
---|
350 | vvh(i-i181,j,numpv)=help |
---|
351 | endif |
---|
352 | endif |
---|
353 | if((isec1(6).eq.052).and.(isec1(7).eq.100)) then |
---|
354 | ! RELATIVE HUMIDITY -> CONVERT TO SPECIFIC HUMIDITY LATER |
---|
355 | if((i.eq.0).and.(j.eq.0)) then |
---|
356 | do ii=1,nuvz |
---|
357 | if ((isec1(8)*100.0).eq.akz(ii)) numprh=ii |
---|
358 | end do |
---|
359 | endif |
---|
360 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
361 | if(i.le.i180) then |
---|
362 | qvh(i179+i,j,numprh,n)=help |
---|
363 | else |
---|
364 | qvh(i-i181,j,numprh,n)=help |
---|
365 | endif |
---|
366 | endif |
---|
367 | if((isec1(6).eq.001).and.(isec1(7).eq.001)) then |
---|
368 | ! SURFACE PRESSURE |
---|
369 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
370 | if(i.le.i180) then |
---|
371 | ps(i179+i,j,1,n)=help |
---|
372 | else |
---|
373 | ps(i-i181,j,1,n)=help |
---|
374 | endif |
---|
375 | endif |
---|
376 | if((isec1(6).eq.039).and.(isec1(7).eq.100)) then |
---|
377 | ! W VELOCITY |
---|
378 | if((i.eq.0).and.(j.eq.0)) then |
---|
379 | do ii=1,nuvz |
---|
380 | if ((isec1(8)*100.0).eq.akz(ii)) numpw=ii |
---|
381 | end do |
---|
382 | endif |
---|
383 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
384 | if(i.le.i180) then |
---|
385 | wwh(i179+i,j,numpw)=help |
---|
386 | else |
---|
387 | wwh(i-i181,j,numpw)=help |
---|
388 | endif |
---|
389 | endif |
---|
390 | if((isec1(6).eq.066).and.(isec1(7).eq.001)) then |
---|
391 | ! SNOW DEPTH |
---|
392 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
393 | if(i.le.i180) then |
---|
394 | sd(i179+i,j,1,n)=help |
---|
395 | else |
---|
396 | sd(i-i181,j,1,n)=help |
---|
397 | endif |
---|
398 | endif |
---|
399 | if((isec1(6).eq.002).and.(isec1(7).eq.102)) then |
---|
400 | ! MEAN SEA LEVEL PRESSURE |
---|
401 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
402 | if(i.le.i180) then |
---|
403 | msl(i179+i,j,1,n)=help |
---|
404 | else |
---|
405 | msl(i-i181,j,1,n)=help |
---|
406 | endif |
---|
407 | endif |
---|
408 | if((isec1(6).eq.071).and.(isec1(7).eq.244)) then |
---|
409 | ! TOTAL CLOUD COVER |
---|
410 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
411 | if(i.le.i180) then |
---|
412 | tcc(i179+i,j,1,n)=help |
---|
413 | else |
---|
414 | tcc(i-i181,j,1,n)=help |
---|
415 | endif |
---|
416 | endif |
---|
417 | if((isec1(6).eq.033).and.(isec1(7).eq.105).and. & |
---|
418 | (isec1(8).eq.10)) then |
---|
419 | ! 10 M U VELOCITY |
---|
420 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
421 | if(i.le.i180) then |
---|
422 | u10(i179+i,j,1,n)=help |
---|
423 | else |
---|
424 | u10(i-i181,j,1,n)=help |
---|
425 | endif |
---|
426 | endif |
---|
427 | if((isec1(6).eq.034).and.(isec1(7).eq.105).and. & |
---|
428 | (isec1(8).eq.10)) then |
---|
429 | ! 10 M V VELOCITY |
---|
430 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
431 | if(i.le.i180) then |
---|
432 | v10(i179+i,j,1,n)=help |
---|
433 | else |
---|
434 | v10(i-i181,j,1,n)=help |
---|
435 | endif |
---|
436 | endif |
---|
437 | if((isec1(6).eq.011).and.(isec1(7).eq.105).and. & |
---|
438 | (isec1(8).eq.02)) then |
---|
439 | ! 2 M TEMPERATURE |
---|
440 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
441 | if(i.le.i180) then |
---|
442 | tt2(i179+i,j,1,n)=help |
---|
443 | else |
---|
444 | tt2(i-i181,j,1,n)=help |
---|
445 | endif |
---|
446 | endif |
---|
447 | if((isec1(6).eq.017).and.(isec1(7).eq.105).and. & |
---|
448 | (isec1(8).eq.02)) then |
---|
449 | ! 2 M DEW POINT TEMPERATURE |
---|
450 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
451 | if(i.le.i180) then |
---|
452 | td2(i179+i,j,1,n)=help |
---|
453 | else |
---|
454 | td2(i-i181,j,1,n)=help |
---|
455 | endif |
---|
456 | endif |
---|
457 | if((isec1(6).eq.062).and.(isec1(7).eq.001)) then |
---|
458 | ! LARGE SCALE PREC. |
---|
459 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
460 | if(i.le.i180) then |
---|
461 | lsprec(i179+i,j,1,n)=help |
---|
462 | else |
---|
463 | lsprec(i-i181,j,1,n)=help |
---|
464 | endif |
---|
465 | endif |
---|
466 | if((isec1(6).eq.063).and.(isec1(7).eq.001)) then |
---|
467 | ! CONVECTIVE PREC. |
---|
468 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
469 | if(i.le.i180) then |
---|
470 | convprec(i179+i,j,1,n)=help |
---|
471 | else |
---|
472 | convprec(i-i181,j,1,n)=help |
---|
473 | endif |
---|
474 | endif |
---|
475 | if((isec1(6).eq.007).and.(isec1(7).eq.001)) then |
---|
476 | ! TOPOGRAPHY |
---|
477 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
478 | if(i.le.i180) then |
---|
479 | oro(i179+i,j)=help |
---|
480 | excessoro(i179+i,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED |
---|
481 | else |
---|
482 | oro(i-i181,j)=help |
---|
483 | excessoro(i-i181,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED |
---|
484 | endif |
---|
485 | endif |
---|
486 | if((isec1(6).eq.081).and.(isec1(7).eq.001)) then |
---|
487 | ! LAND SEA MASK |
---|
488 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
489 | if(i.le.i180) then |
---|
490 | lsm(i179+i,j)=help |
---|
491 | else |
---|
492 | lsm(i-i181,j)=help |
---|
493 | endif |
---|
494 | endif |
---|
495 | if((isec1(6).eq.221).and.(isec1(7).eq.001)) then |
---|
496 | ! MIXING HEIGHT |
---|
497 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
498 | if(i.le.i180) then |
---|
499 | hmix(i179+i,j,1,n)=help |
---|
500 | else |
---|
501 | hmix(i-i181,j,1,n)=help |
---|
502 | endif |
---|
503 | endif |
---|
504 | if((isec1(6).eq.052).and.(isec1(7).eq.105).and. & |
---|
505 | (isec1(8).eq.02)) then |
---|
506 | ! 2 M RELATIVE HUMIDITY |
---|
507 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
508 | if(i.le.i180) then |
---|
509 | qvh2(i179+i,j)=help |
---|
510 | else |
---|
511 | qvh2(i-i181,j)=help |
---|
512 | endif |
---|
513 | endif |
---|
514 | if((isec1(6).eq.011).and.(isec1(7).eq.107)) then |
---|
515 | ! TEMPERATURE LOWEST SIGMA LEVEL |
---|
516 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
517 | if(i.le.i180) then |
---|
518 | tlev1(i179+i,j)=help |
---|
519 | else |
---|
520 | tlev1(i-i181,j)=help |
---|
521 | endif |
---|
522 | endif |
---|
523 | if((isec1(6).eq.033).and.(isec1(7).eq.107)) then |
---|
524 | ! U VELOCITY LOWEST SIGMA LEVEL |
---|
525 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
526 | if(i.le.i180) then |
---|
527 | ulev1(i179+i,j)=help |
---|
528 | else |
---|
529 | ulev1(i-i181,j)=help |
---|
530 | endif |
---|
531 | endif |
---|
532 | if((isec1(6).eq.034).and.(isec1(7).eq.107)) then |
---|
533 | ! V VELOCITY LOWEST SIGMA LEVEL |
---|
534 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
535 | if(i.le.i180) then |
---|
536 | vlev1(i179+i,j)=help |
---|
537 | else |
---|
538 | vlev1(i-i181,j)=help |
---|
539 | endif |
---|
540 | endif |
---|
541 | ! SEC & IP 12/2018 read GFS clouds |
---|
542 | if((isec1(6).eq.153).and.(isec1(7).eq.100)) then !! CLWCR Cloud liquid water content [kg/kg] |
---|
543 | if((i.eq.0).and.(j.eq.0)) then |
---|
544 | do ii=1,nuvz |
---|
545 | if ((isec1(8)*100.0).eq.akz(ii)) numpclwch=ii |
---|
546 | end do |
---|
547 | endif |
---|
548 | help=zsec4(nxfield*(ny-j-1)+i+1) |
---|
549 | if(i.le.i180) then |
---|
550 | clwch(i179+i,j,numpclwch,n)=help |
---|
551 | else |
---|
552 | clwch(i-i181,j,numpclwch,n)=help |
---|
553 | endif |
---|
554 | readclouds=.true. |
---|
555 | sumclouds=.true. |
---|
556 | ! readclouds=.false. |
---|
557 | ! sumclouds=.false. |
---|
558 | endif |
---|
559 | |
---|
560 | |
---|
561 | end do |
---|
562 | end do |
---|
563 | |
---|
564 | endif |
---|
565 | |
---|
566 | if((isec1(6).eq.33).and.(isec1(7).eq.100)) then |
---|
567 | ! NCEP ISOBARIC LEVELS |
---|
568 | iumax=iumax+1 |
---|
569 | endif |
---|
570 | |
---|
571 | call grib_release(igrib) |
---|
572 | goto 10 !! READ NEXT LEVEL OR PARAMETER |
---|
573 | ! |
---|
574 | ! CLOSING OF INPUT DATA FILE |
---|
575 | ! |
---|
576 | |
---|
577 | !HSO close grib file |
---|
578 | 50 continue |
---|
579 | call grib_close_file(ifile) |
---|
580 | |
---|
581 | ! SENS. HEAT FLUX |
---|
582 | sshf(:,:,1,n)=0.0 ! not available from gfs.tccz.pgrbfxx files |
---|
583 | hflswitch=.false. ! Heat flux not available |
---|
584 | ! SOLAR RADIATIVE FLUXES |
---|
585 | ssr(:,:,1,n)=0.0 ! not available from gfs.tccz.pgrbfxx files |
---|
586 | ! EW SURFACE STRESS |
---|
587 | ewss=0.0 ! not available from gfs.tccz.pgrbfxx files |
---|
588 | ! NS SURFACE STRESS |
---|
589 | nsss=0.0 ! not available from gfs.tccz.pgrbfxx files |
---|
590 | strswitch=.false. ! stress not available |
---|
591 | |
---|
592 | ! CONVERT TP TO LSP (GRIB2 only) |
---|
593 | if (gribVer.eq.2) then |
---|
594 | do j=0,nymin1 |
---|
595 | do i=0,nxfield-1 |
---|
596 | if(i.le.i180) then |
---|
597 | if (convprec(i179+i,j,1,n).lt.lsprec(i179+i,j,1,n)) then ! neg precip would occur |
---|
598 | lsprec(i179+i,j,1,n)= & |
---|
599 | lsprec(i179+i,j,1,n)-convprec(i179+i,j,1,n) |
---|
600 | else |
---|
601 | lsprec(i179+i,j,1,n)=0 |
---|
602 | endif |
---|
603 | else |
---|
604 | if (convprec(i-i181,j,1,n).lt.lsprec(i-i181,j,1,n)) then |
---|
605 | lsprec(i-i181,j,1,n)= & |
---|
606 | lsprec(i-i181,j,1,n)-convprec(i-i181,j,1,n) |
---|
607 | else |
---|
608 | lsprec(i-i181,j,1,n)=0 |
---|
609 | endif |
---|
610 | endif |
---|
611 | enddo |
---|
612 | enddo |
---|
613 | endif |
---|
614 | !HSO end edits |
---|
615 | |
---|
616 | |
---|
617 | ! TRANSFORM RH TO SPECIFIC HUMIDITY |
---|
618 | |
---|
619 | do j=0,ny-1 |
---|
620 | do i=0,nxfield-1 |
---|
621 | do k=1,nuvz |
---|
622 | help=qvh(i,j,k,n) |
---|
623 | temp=tth(i,j,k,n) |
---|
624 | plev1=akm(k)+bkm(k)*ps(i,j,1,n) |
---|
625 | elev=ew(temp)*help/100.0 |
---|
626 | qvh(i,j,k,n)=xmwml*(elev/(plev1-((1.0-xmwml)*elev))) |
---|
627 | end do |
---|
628 | end do |
---|
629 | end do |
---|
630 | |
---|
631 | ! CALCULATE 2 M DEW POINT FROM 2 M RELATIVE HUMIDITY |
---|
632 | ! USING BOLTON'S (1980) FORMULA |
---|
633 | ! BECAUSE td2 IS NOT AVAILABLE FROM NCEP GFS DATA |
---|
634 | |
---|
635 | do j=0,ny-1 |
---|
636 | do i=0,nxfield-1 |
---|
637 | help=qvh2(i,j) |
---|
638 | temp=tt2(i,j,1,n) |
---|
639 | elev=ew(temp)/100.*help/100. !vapour pressure in hPa |
---|
640 | td2(i,j,1,n)=243.5/(17.67/log(elev/6.112)-1)+273. |
---|
641 | if (help.le.0.) td2(i,j,1,n)=tt2(i,j,1,n) |
---|
642 | end do |
---|
643 | end do |
---|
644 | |
---|
645 | if(levdiff2.eq.0) then |
---|
646 | iwmax=nlev_ec+1 |
---|
647 | do i=0,nxmin1 |
---|
648 | do j=0,nymin1 |
---|
649 | wwh(i,j,nlev_ec+1)=0. |
---|
650 | end do |
---|
651 | end do |
---|
652 | endif |
---|
653 | |
---|
654 | |
---|
655 | ! For global fields, assign the leftmost data column also to the rightmost |
---|
656 | ! data column; if required, shift whole grid by nxshift grid points |
---|
657 | !************************************************************************* |
---|
658 | |
---|
659 | if (xglobal) then |
---|
660 | call shift_field_0(ewss,nxfield,ny) |
---|
661 | call shift_field_0(nsss,nxfield,ny) |
---|
662 | call shift_field_0(oro,nxfield,ny) |
---|
663 | call shift_field_0(excessoro,nxfield,ny) |
---|
664 | call shift_field_0(lsm,nxfield,ny) |
---|
665 | call shift_field_0(ulev1,nxfield,ny) |
---|
666 | call shift_field_0(vlev1,nxfield,ny) |
---|
667 | call shift_field_0(tlev1,nxfield,ny) |
---|
668 | call shift_field_0(qvh2,nxfield,ny) |
---|
669 | call shift_field(ps,nxfield,ny,1,1,2,n) |
---|
670 | call shift_field(sd,nxfield,ny,1,1,2,n) |
---|
671 | call shift_field(msl,nxfield,ny,1,1,2,n) |
---|
672 | call shift_field(tcc,nxfield,ny,1,1,2,n) |
---|
673 | call shift_field(u10,nxfield,ny,1,1,2,n) |
---|
674 | call shift_field(v10,nxfield,ny,1,1,2,n) |
---|
675 | call shift_field(tt2,nxfield,ny,1,1,2,n) |
---|
676 | call shift_field(td2,nxfield,ny,1,1,2,n) |
---|
677 | call shift_field(lsprec,nxfield,ny,1,1,2,n) |
---|
678 | call shift_field(convprec,nxfield,ny,1,1,2,n) |
---|
679 | call shift_field(sshf,nxfield,ny,1,1,2,n) |
---|
680 | call shift_field(ssr,nxfield,ny,1,1,2,n) |
---|
681 | call shift_field(hmix,nxfield,ny,1,1,2,n) |
---|
682 | call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n) |
---|
683 | call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n) |
---|
684 | call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1) |
---|
685 | call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1) |
---|
686 | call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1) |
---|
687 | ! IP & SEC adding GFS Clouds 20181205 |
---|
688 | call shift_field(clwch,nxfield,ny,nuvzmax,nuvz,2,n) |
---|
689 | endif |
---|
690 | |
---|
691 | do i=0,nxmin1 |
---|
692 | do j=0,nymin1 |
---|
693 | ! Convert precip. from mm/s -> mm/hour |
---|
694 | convprec(i,j,1,n)=convprec(i,j,1,n)*3600. |
---|
695 | lsprec(i,j,1,n)=lsprec(i,j,1,n)*3600. |
---|
696 | surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2) |
---|
697 | end do |
---|
698 | end do |
---|
699 | |
---|
700 | if ((.not.hflswitch).or.(.not.strswitch)) then |
---|
701 | ! write(*,*) 'WARNING: No flux data contained in GRIB file ', |
---|
702 | ! + wfname(indj) |
---|
703 | |
---|
704 | ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD |
---|
705 | !*************************************************************************** |
---|
706 | |
---|
707 | do i=0,nxmin1 |
---|
708 | do j=0,nymin1 |
---|
709 | hlev1=30.0 ! HEIGHT OF FIRST MODEL SIGMA LAYER |
---|
710 | ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2) |
---|
711 | fflev1=sqrt(ulev1(i,j)**2+vlev1(i,j)**2) |
---|
712 | call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, & |
---|
713 | tt2(i,j,1,n),tlev1(i,j),ff10m,fflev1, & |
---|
714 | surfstr(i,j,1,n),sshf(i,j,1,n)) |
---|
715 | if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200. |
---|
716 | if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400. |
---|
717 | end do |
---|
718 | end do |
---|
719 | endif |
---|
720 | |
---|
721 | if(iumax.ne.nuvz) stop 'READWIND: NUVZ NOT CONSISTENT' |
---|
722 | if(iumax.ne.nwz) stop 'READWIND: NWZ NOT CONSISTENT' |
---|
723 | |
---|
724 | return |
---|
725 | 888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### ' |
---|
726 | write(*,*) ' #### ',wfname(indj),' #### ' |
---|
727 | write(*,*) ' #### IS NOT GRIB FORMAT !!! #### ' |
---|
728 | stop 'Execution terminated' |
---|
729 | 999 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### ' |
---|
730 | write(*,*) ' #### ',wfname(indj),' #### ' |
---|
731 | write(*,*) ' #### CANNOT BE OPENED !!! #### ' |
---|
732 | stop 'Execution terminated' |
---|
733 | |
---|
734 | end subroutine readwind_gfs |
---|