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 gridcheck |
---|
23 | |
---|
24 | !********************************************************************** |
---|
25 | ! * |
---|
26 | ! FLEXPART MODEL SUBROUTINE GRIDCHECK * |
---|
27 | ! * |
---|
28 | !********************************************************************** |
---|
29 | ! * |
---|
30 | ! AUTHOR: G. WOTAWA * |
---|
31 | ! DATE: 1997-08-06 * |
---|
32 | ! LAST UPDATE: 1997-10-10 * |
---|
33 | ! * |
---|
34 | ! Update: 1999-02-08, global fields allowed, A. Stohl* |
---|
35 | ! CHANGE: 17/11/2005, Caroline Forster, GFS data * |
---|
36 | ! * |
---|
37 | !********************************************************************** |
---|
38 | ! * |
---|
39 | ! DESCRIPTION: * |
---|
40 | ! * |
---|
41 | ! THIS SUBROUTINE DETERMINES THE GRID SPECIFICATIONS (LOWER LEFT * |
---|
42 | ! LONGITUDE, LOWER LEFT LATITUDE, NUMBER OF GRID POINTS, GRID DIST- * |
---|
43 | ! ANCE AND VERTICAL DISCRETIZATION OF THE ECMWF MODEL) FROM THE * |
---|
44 | ! GRIB HEADER OF THE FIRST INPUT FILE. THE CONSISTANCY (NO CHANGES * |
---|
45 | ! WITHIN ONE FLEXPART RUN) IS CHECKED IN THE ROUTINE "READWIND" AT * |
---|
46 | ! ANY CALL. * |
---|
47 | ! * |
---|
48 | ! XLON0 geographical longitude of lower left gridpoint * |
---|
49 | ! YLAT0 geographical latitude of lower left gridpoint * |
---|
50 | ! NX number of grid points x-direction * |
---|
51 | ! NY number of grid points y-direction * |
---|
52 | ! DX grid distance x-direction * |
---|
53 | ! DY grid distance y-direction * |
---|
54 | ! NUVZ number of grid points for horizontal wind * |
---|
55 | ! components in z direction * |
---|
56 | ! NWZ number of grid points for vertical wind * |
---|
57 | ! component in z direction * |
---|
58 | ! sizesouth, sizenorth give the map scale (i.e. number of virtual grid* |
---|
59 | ! points of the polar stereographic grid): * |
---|
60 | ! used to check the CFL criterion * |
---|
61 | ! UVHEIGHT(1)- heights of gridpoints where u and v are * |
---|
62 | ! UVHEIGHT(NUVZ) given * |
---|
63 | ! WHEIGHT(1)- heights of gridpoints where w is given * |
---|
64 | ! WHEIGHT(NWZ) * |
---|
65 | ! * |
---|
66 | !********************************************************************** |
---|
67 | |
---|
68 | use par_mod |
---|
69 | use com_mod |
---|
70 | use conv_mod |
---|
71 | |
---|
72 | implicit none |
---|
73 | |
---|
74 | integer :: ix,jy,i,ifn,ifield,j,k,iumax,iwmax,numskip |
---|
75 | real :: xaux1,xaux2,yaux1,yaux2,sizesouth,sizenorth,xauxa,pint |
---|
76 | real :: akm_usort(nwzmax) |
---|
77 | |
---|
78 | ! NCEP GFS |
---|
79 | real :: pres(nwzmax), help |
---|
80 | |
---|
81 | ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING |
---|
82 | |
---|
83 | ! dimension of isec2 at least (22+n), where n is the number of parallels or |
---|
84 | ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid |
---|
85 | |
---|
86 | ! dimension of zsec2 at least (10+nn), where nn is the number of vertical |
---|
87 | ! coordinate parameters |
---|
88 | |
---|
89 | integer :: isec0(2),isec1(56),isec2(22+nxmax+nymax),isec3(2) |
---|
90 | integer :: isec4(64),inbuff(jpack),ilen,iswap,ierr,iword,lunit |
---|
91 | real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp) |
---|
92 | character(len=1) :: opt, yoper = 'D' |
---|
93 | |
---|
94 | |
---|
95 | iumax=0 |
---|
96 | iwmax=0 |
---|
97 | |
---|
98 | if(ideltas.gt.0) then |
---|
99 | ifn=1 |
---|
100 | else |
---|
101 | ifn=numbwf |
---|
102 | endif |
---|
103 | ! |
---|
104 | ! OPENING OF DATA FILE (GRIB CODE) |
---|
105 | ! |
---|
106 | 5 call pbopen(lunit,path(3)(1:length(3))//wfname(ifn),'r',ierr) |
---|
107 | if(ierr.lt.0) goto 999 |
---|
108 | |
---|
109 | ifield=0 |
---|
110 | 10 ifield=ifield+1 |
---|
111 | ! |
---|
112 | ! GET NEXT FIELDS |
---|
113 | ! |
---|
114 | call pbgrib(lunit,inbuff,jpack,ilen,ierr) |
---|
115 | if(ierr.eq.-1) goto 30 ! EOF DETECTED |
---|
116 | if(ierr.lt.-1) goto 999 ! ERROR DETECTED |
---|
117 | |
---|
118 | ierr=1 |
---|
119 | |
---|
120 | ! Check whether we are on a little endian or on a big endian computer |
---|
121 | !******************************************************************** |
---|
122 | |
---|
123 | !if (inbuff(1).eq.1112101447) then ! little endian, swap bytes |
---|
124 | ! iswap=1+ilen/4 |
---|
125 | ! call swap32(inbuff,iswap) |
---|
126 | !else if (inbuff(1).ne.1196575042) then ! big endian |
---|
127 | ! stop 'subroutine gridcheck: corrupt GRIB data' |
---|
128 | !endif |
---|
129 | |
---|
130 | call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, & |
---|
131 | zsec4,jpunp,inbuff,jpack,iword,yoper,ierr) |
---|
132 | if (ierr.ne.0) goto 10 ! ERROR DETECTED |
---|
133 | |
---|
134 | if(ifield.eq.1) then |
---|
135 | nxfield=isec2(2) |
---|
136 | ny=isec2(3) |
---|
137 | xaux1=real(isec2(5))/1000. |
---|
138 | xaux2=real(isec2(8))/1000. |
---|
139 | if((xaux1.eq.0.).and.(xaux2.eq.-1.0)) then ! NCEP DATA FROM 0 TO |
---|
140 | xaux1=-179.0 ! 359 DEG EAST -> |
---|
141 | xaux2=180.0 ! TRANSFORMED TO -179 |
---|
142 | endif ! TO 180 DEG EAST |
---|
143 | if (xaux1.gt.180) xaux1=xaux1-360.0 |
---|
144 | if (xaux2.gt.180) xaux2=xaux2-360.0 |
---|
145 | if (xaux1.lt.-180) xaux1=xaux1+360.0 |
---|
146 | if (xaux2.lt.-180) xaux2=xaux2+360.0 |
---|
147 | if (xaux2.lt.xaux1) xaux2=xaux2+360. |
---|
148 | yaux1=real(isec2(7))/1000. |
---|
149 | yaux2=real(isec2(4))/1000. |
---|
150 | xlon0=xaux1 |
---|
151 | ylat0=yaux1 |
---|
152 | dx=(xaux2-xaux1)/real(nxfield-1) |
---|
153 | dy=(yaux2-yaux1)/real(ny-1) |
---|
154 | dxconst=180./(dx*r_earth*pi) |
---|
155 | dyconst=180./(dy*r_earth*pi) |
---|
156 | |
---|
157 | |
---|
158 | ! Check whether fields are global |
---|
159 | ! If they contain the poles, specify polar stereographic map |
---|
160 | ! projections using the stlmbr- and stcm2p-calls |
---|
161 | !*********************************************************** |
---|
162 | |
---|
163 | xauxa=abs(xaux2+dx-360.-xaux1) |
---|
164 | if (xauxa.lt.0.001) then |
---|
165 | nx=nxfield+1 ! field is cyclic |
---|
166 | xglobal=.true. |
---|
167 | if (abs(nxshift).ge.nx) & |
---|
168 | stop 'nxshift in file par_mod is too large' |
---|
169 | xlon0=xlon0+real(nxshift)*dx |
---|
170 | else |
---|
171 | nx=nxfield |
---|
172 | xglobal=.false. |
---|
173 | if (nxshift.ne.0) & |
---|
174 | stop 'nxshift (par_mod) must be zero for non-global domain' |
---|
175 | endif |
---|
176 | nxmin1=nx-1 |
---|
177 | nymin1=ny-1 |
---|
178 | if (xlon0.gt.180.) xlon0=xlon0-360. |
---|
179 | xauxa=abs(yaux1+90.) |
---|
180 | if (xglobal.and.xauxa.lt.0.001) then |
---|
181 | sglobal=.true. ! field contains south pole |
---|
182 | ! Enhance the map scale by factor 3 (*2=6) compared to north-south |
---|
183 | ! map scale |
---|
184 | sizesouth=6.*(switchsouth+90.)/dy |
---|
185 | call stlmbr(southpolemap,-90.,0.) |
---|
186 | call stcm2p(southpolemap,0.,0.,switchsouth,0.,sizesouth, & |
---|
187 | sizesouth,switchsouth,180.) |
---|
188 | switchsouthg=(switchsouth-ylat0)/dy |
---|
189 | else |
---|
190 | sglobal=.false. |
---|
191 | switchsouthg=999999. |
---|
192 | endif |
---|
193 | xauxa=abs(yaux2-90.) |
---|
194 | if (xglobal.and.xauxa.lt.0.001) then |
---|
195 | nglobal=.true. ! field contains north pole |
---|
196 | ! Enhance the map scale by factor 3 (*2=6) compared to north-south |
---|
197 | ! map scale |
---|
198 | sizenorth=6.*(90.-switchnorth)/dy |
---|
199 | call stlmbr(northpolemap,90.,0.) |
---|
200 | call stcm2p(northpolemap,0.,0.,switchnorth,0.,sizenorth, & |
---|
201 | sizenorth,switchnorth,180.) |
---|
202 | switchnorthg=(switchnorth-ylat0)/dy |
---|
203 | else |
---|
204 | nglobal=.false. |
---|
205 | switchnorthg=999999. |
---|
206 | endif |
---|
207 | endif |
---|
208 | |
---|
209 | if (nxshift.lt.0) stop 'nxshift (par_mod) must not be negative' |
---|
210 | if (nxshift.ge.nxfield) stop 'nxshift (par_mod) too large' |
---|
211 | |
---|
212 | ! NCEP ISOBARIC LEVELS |
---|
213 | !********************* |
---|
214 | |
---|
215 | if((isec1(6).eq.33).and.(isec1(7).eq.100)) then |
---|
216 | iumax=iumax+1 |
---|
217 | pres(iumax)=real(isec1(8))*100.0 |
---|
218 | endif |
---|
219 | |
---|
220 | ! NCEP TERRAIN |
---|
221 | !************* |
---|
222 | |
---|
223 | if((isec1(6).eq.007).and.(isec1(7).eq.001)) then |
---|
224 | do jy=0,ny-1 |
---|
225 | do ix=0,nxfield-1 |
---|
226 | help=zsec4(nxfield*(ny-jy-1)+ix+1) |
---|
227 | if(ix.le.180) then |
---|
228 | oro(179+ix,jy)=help |
---|
229 | excessoro(179+ix,jy)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED |
---|
230 | else |
---|
231 | oro(ix-181,jy)=help |
---|
232 | excessoro(ix-181,jy)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED |
---|
233 | endif |
---|
234 | end do |
---|
235 | end do |
---|
236 | endif |
---|
237 | |
---|
238 | ! NCEP LAND SEA MASK |
---|
239 | !******************* |
---|
240 | |
---|
241 | if((isec1(6).eq.081).and.(isec1(7).eq.001)) then |
---|
242 | do jy=0,ny-1 |
---|
243 | do ix=0,nxfield-1 |
---|
244 | help=zsec4(nxfield*(ny-jy-1)+ix+1) |
---|
245 | if(ix.le.180) then |
---|
246 | lsm(179+ix,jy)=help |
---|
247 | else |
---|
248 | lsm(ix-181,jy)=help |
---|
249 | endif |
---|
250 | end do |
---|
251 | end do |
---|
252 | endif |
---|
253 | |
---|
254 | goto 10 !! READ NEXT LEVEL OR PARAMETER |
---|
255 | ! |
---|
256 | ! CLOSING OF INPUT DATA FILE |
---|
257 | ! |
---|
258 | 30 call pbclose(lunit,ierr) !! FINISHED READING / CLOSING GRIB FILE |
---|
259 | |
---|
260 | nuvz=iumax |
---|
261 | nwz =iumax |
---|
262 | nlev_ec=iumax |
---|
263 | |
---|
264 | if (nx.gt.nxmax) then |
---|
265 | write(*,*) 'FLEXPART error: Too many grid points in x direction.' |
---|
266 | write(*,*) 'Reduce resolution of wind fields.' |
---|
267 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
268 | write(*,*) nx,nxmax |
---|
269 | stop |
---|
270 | endif |
---|
271 | |
---|
272 | if (ny.gt.nymax) then |
---|
273 | write(*,*) 'FLEXPART error: Too many grid points in y direction.' |
---|
274 | write(*,*) 'Reduce resolution of wind fields.' |
---|
275 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
276 | write(*,*) ny,nymax |
---|
277 | stop |
---|
278 | endif |
---|
279 | |
---|
280 | if (nuvz.gt.nuvzmax) then |
---|
281 | write(*,*) 'FLEXPART error: Too many u,v grid points in z '// & |
---|
282 | 'direction.' |
---|
283 | write(*,*) 'Reduce resolution of wind fields.' |
---|
284 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
285 | write(*,*) nuvz,nuvzmax |
---|
286 | stop |
---|
287 | endif |
---|
288 | |
---|
289 | if (nwz.gt.nwzmax) then |
---|
290 | write(*,*) 'FLEXPART error: Too many w grid points in z '// & |
---|
291 | 'direction.' |
---|
292 | write(*,*) 'Reduce resolution of wind fields.' |
---|
293 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
294 | write(*,*) nwz,nwzmax |
---|
295 | stop |
---|
296 | endif |
---|
297 | |
---|
298 | ! If desired, shift all grids by nxshift grid cells |
---|
299 | !************************************************** |
---|
300 | |
---|
301 | if (xglobal) then |
---|
302 | call shift_field_0(oro,nxfield,ny) |
---|
303 | call shift_field_0(lsm,nxfield,ny) |
---|
304 | call shift_field_0(excessoro,nxfield,ny) |
---|
305 | endif |
---|
306 | |
---|
307 | ! Output of grid info |
---|
308 | !******************** |
---|
309 | |
---|
310 | write(*,*) |
---|
311 | write(*,*) |
---|
312 | write(*,'(a,2i7)') '# of vertical levels in ECMWF data: ', & |
---|
313 | nuvz,nwz |
---|
314 | write(*,*) |
---|
315 | write(*,'(a)') 'Mother domain:' |
---|
316 | write(*,'(a,f10.1,a1,f10.1,a,f10.1)') ' Longitude range: ', & |
---|
317 | xlon0,' to ',xlon0+(nx-1)*dx,' Grid distance: ',dx |
---|
318 | write(*,'(a,f10.1,a1,f10.1,a,f10.1)') ' Latitude range: ', & |
---|
319 | ylat0,' to ',ylat0+(ny-1)*dy,' Grid distance: ',dy |
---|
320 | write(*,*) |
---|
321 | |
---|
322 | |
---|
323 | ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL |
---|
324 | ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM |
---|
325 | |
---|
326 | numskip=nlev_ec-nuvz ! number of ecmwf model layers not used |
---|
327 | ! by trajectory model |
---|
328 | do i=1,nwz |
---|
329 | j=numskip+i |
---|
330 | k=nlev_ec+1+numskip+i |
---|
331 | akm_usort(nwz-i+1)=pres(nwz-i+1) |
---|
332 | bkm(nwz-i+1)=0.0 |
---|
333 | end do |
---|
334 | |
---|
335 | !****************************** |
---|
336 | ! change Sabine Eckhardt: akm should always be in descending order ... readwind adapted! |
---|
337 | !****************************** |
---|
338 | do i=1,nwz |
---|
339 | if (akm_usort(1).gt.akm_usort(2)) then |
---|
340 | akm(i)=akm_usort(i) |
---|
341 | else |
---|
342 | akm(i)=akm_usort(nwz-i+1) |
---|
343 | endif |
---|
344 | end do |
---|
345 | |
---|
346 | ! |
---|
347 | ! CALCULATION OF AKZ, BKZ |
---|
348 | ! AKZ,BKZ: model discretization parameters at the center of each model |
---|
349 | ! layer |
---|
350 | ! |
---|
351 | ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0, |
---|
352 | ! i.e. ground level |
---|
353 | !***************************************************************************** |
---|
354 | |
---|
355 | do i=1,nuvz |
---|
356 | akz(i)=akm(i) |
---|
357 | bkz(i)=bkm(i) |
---|
358 | end do |
---|
359 | |
---|
360 | ! NOTE: In FLEXPART versions up to 4.0, the number of model levels was doubled |
---|
361 | ! upon the transformation to z levels. In order to save computer memory, this is |
---|
362 | ! not done anymore in the standard version. However, this option can still be |
---|
363 | ! switched on by replacing the following lines with those below, that are |
---|
364 | ! currently commented out. For this, similar changes are necessary in |
---|
365 | ! verttransform.f and verttranform_nests.f |
---|
366 | !***************************************************************************** |
---|
367 | |
---|
368 | nz=nuvz |
---|
369 | if (nz.gt.nzmax) stop 'nzmax too small' |
---|
370 | do i=1,nuvz |
---|
371 | aknew(i)=akz(i) |
---|
372 | bknew(i)=bkz(i) |
---|
373 | end do |
---|
374 | |
---|
375 | ! Switch on following lines to use doubled vertical resolution |
---|
376 | !************************************************************* |
---|
377 | !nz=nuvz+nwz-1 |
---|
378 | !if (nz.gt.nzmax) stop 'nzmax too small' |
---|
379 | !do 100 i=1,nwz |
---|
380 | ! aknew(2*(i-1)+1)=akm(i) |
---|
381 | !00 bknew(2*(i-1)+1)=bkm(i) |
---|
382 | !do 110 i=2,nuvz |
---|
383 | ! aknew(2*(i-1))=akz(i) |
---|
384 | !10 bknew(2*(i-1))=bkz(i) |
---|
385 | ! End doubled vertical resolution |
---|
386 | |
---|
387 | |
---|
388 | ! Determine the uppermost level for which the convection scheme shall be applied |
---|
389 | ! by assuming that there is no convection above 50 hPa (for standard SLP) |
---|
390 | !***************************************************************************** |
---|
391 | |
---|
392 | do i=1,nuvz-2 |
---|
393 | pint=akz(i)+bkz(i)*101325. |
---|
394 | if (pint.lt.5000.) goto 96 |
---|
395 | end do |
---|
396 | 96 nconvlev=i |
---|
397 | if (nconvlev.gt.nconvlevmax-1) then |
---|
398 | nconvlev=nconvlevmax-1 |
---|
399 | write(*,*) 'Attention, convection only calculated up to ', & |
---|
400 | akz(nconvlev)+bkz(nconvlev)*1013.25,' hPa' |
---|
401 | endif |
---|
402 | |
---|
403 | return |
---|
404 | |
---|
405 | 999 write(*,*) |
---|
406 | write(*,*) ' ###########################################'// & |
---|
407 | '###### ' |
---|
408 | write(*,*) ' TRAJECTORY MODEL SUBROUTINE GRIDCHECK:' |
---|
409 | write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfname(ifn) |
---|
410 | write(*,*) ' ###########################################'// & |
---|
411 | '###### ' |
---|
412 | write(*,*) |
---|
413 | write(*,'(a)') '!!! PLEASE INSERT A NEW CD-ROM AND !!!' |
---|
414 | write(*,'(a)') '!!! PRESS ANY KEY TO CONTINUE... !!!' |
---|
415 | write(*,'(a)') '!!! ...OR TERMINATE FLEXPART PRESSING!!!' |
---|
416 | write(*,'(a)') '!!! THE "X" KEY... !!!' |
---|
417 | write(*,*) |
---|
418 | read(*,'(a)') opt |
---|
419 | if(opt.eq.'X') then |
---|
420 | stop |
---|
421 | else |
---|
422 | goto 5 |
---|
423 | endif |
---|
424 | |
---|
425 | end subroutine gridcheck |
---|