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_nests |
---|
23 | |
---|
24 | !***************************************************************************** |
---|
25 | ! * |
---|
26 | ! This routine checks the grid specification for the nested model * |
---|
27 | ! domains. It is similar to subroutine gridcheck, which checks the * |
---|
28 | ! mother domain. * |
---|
29 | ! * |
---|
30 | ! Authors: A. Stohl, G. Wotawa * |
---|
31 | ! * |
---|
32 | ! 8 February 1999 * |
---|
33 | ! * |
---|
34 | !***************************************************************************** |
---|
35 | ! CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with ECMWF grib_api * |
---|
36 | ! CHANGE: 03/12/2008, Harald Sodemann, change to f90 grib_api * |
---|
37 | !***************************************************************************** |
---|
38 | |
---|
39 | use grib_api |
---|
40 | use par_mod |
---|
41 | use com_mod |
---|
42 | use class_vtable |
---|
43 | |
---|
44 | implicit none |
---|
45 | |
---|
46 | !HSO parameters for grib_api |
---|
47 | integer :: ifile |
---|
48 | integer :: iret |
---|
49 | integer :: igrib |
---|
50 | integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl |
---|
51 | #if defined CTBTO |
---|
52 | integer :: parId |
---|
53 | #endif |
---|
54 | integer :: gotGrib |
---|
55 | !HSO end |
---|
56 | integer :: i,j,k,l,ifn,ifield,iumax,iwmax,numskip,nlev_ecn |
---|
57 | integer :: nuvzn,nwzn |
---|
58 | real :: akmn(nwzmax),bkmn(nwzmax),akzn(nuvzmax),bkzn(nuvzmax) |
---|
59 | real(kind=4) :: xaux1,xaux2,yaux1,yaux2 |
---|
60 | real(kind=8) :: xaux1in,xaux2in,yaux1in,yaux2in |
---|
61 | |
---|
62 | ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING |
---|
63 | |
---|
64 | ! dimension of isec2 at least (22+n), where n is the number of parallels or |
---|
65 | ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid |
---|
66 | |
---|
67 | ! dimension of zsec2 at least (10+nn), where nn is the number of vertical |
---|
68 | ! coordinate parameters |
---|
69 | |
---|
70 | integer :: isec1(56),isec2(22+nxmaxn+nymaxn) |
---|
71 | real(kind=4) :: zsec2(60+2*nuvzmax),zsec4(jpunp) |
---|
72 | #if defined CTBTO |
---|
73 | real(kind=4) ::conversion_factor |
---|
74 | #endif |
---|
75 | |
---|
76 | !HSO grib api error messages |
---|
77 | character(len=24) :: gribErrorMsg = 'Error reading grib file' |
---|
78 | character(len=20) :: gribFunction = 'gridcheck_nests' |
---|
79 | |
---|
80 | |
---|
81 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
82 | !!!! Vtable related variables |
---|
83 | |
---|
84 | ! Paths to Vtables (current implementation will assume they are in the cwd |
---|
85 | ! This is a crappy place for these parameters. Need to move them. |
---|
86 | character(LEN=255), parameter :: VTABLE_ECMWF_GRIB1_PATH = & |
---|
87 | "Vtable_ecmwf_grib1", & |
---|
88 | VTABLE_ECMWF_GRIB2_PATH = & |
---|
89 | "Vtable_ecmwf_grib2", & |
---|
90 | VTABLE_ECMWF_GRIB1_2_PATH = & |
---|
91 | "Vtable_ecmwf_grib1_2" |
---|
92 | |
---|
93 | integer :: gribfile_type |
---|
94 | integer :: current_grib_level ! This "was" isec1(8) in previous version |
---|
95 | character(len=255) :: gribfile_name |
---|
96 | character(len=255) :: vtable_path |
---|
97 | character(len=15) :: fpname |
---|
98 | |
---|
99 | type(Vtable) :: my_vtable ! unallocated |
---|
100 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
101 | |
---|
102 | |
---|
103 | |
---|
104 | |
---|
105 | xresoln(0)=1. ! resolution enhancement for mother grid |
---|
106 | yresoln(0)=1. ! resolution enhancement for mother grid |
---|
107 | |
---|
108 | ! Loop about all nesting levels |
---|
109 | !****************************** |
---|
110 | |
---|
111 | do l=1,numbnests |
---|
112 | |
---|
113 | iumax=0 |
---|
114 | iwmax=0 |
---|
115 | |
---|
116 | if(ideltas.gt.0) then |
---|
117 | ifn=1 |
---|
118 | else |
---|
119 | ifn=numbwf |
---|
120 | endif |
---|
121 | |
---|
122 | |
---|
123 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
124 | !!!!!!!!!!!!!!!!!!! VTABLE code |
---|
125 | !!!!!!! Vtable choice |
---|
126 | gribfile_name = path(3)(1:length(3))//trim(wfname(ifn)) |
---|
127 | print *, 'gribfile_name: ', gribfile_name |
---|
128 | |
---|
129 | gribfile_type = vtable_detect_gribfile_type( gribfile_name ) |
---|
130 | |
---|
131 | print *, 'gribfile_type: ', gribfile_type |
---|
132 | |
---|
133 | if (gribfile_type .eq. VTABLE_GRIBFILE_TYPE_ECMWF_GRIB1) then |
---|
134 | vtable_path = VTABLE_ECMWF_GRIB1_PATH |
---|
135 | else if (gribfile_type .eq. VTABLE_GRIBFILE_TYPE_ECMWF_GRIB2) then |
---|
136 | vtable_path = VTABLE_ECMWF_GRIB2_PATH |
---|
137 | else if (gribfile_type .eq. VTABLE_GRIBFILE_TYPE_ECMWF_GRIB1_2) then |
---|
138 | vtable_path = VTABLE_ECMWF_GRIB1_2_PATH |
---|
139 | else |
---|
140 | print *, 'Unsupported gribfile_type: ', gribfile_type |
---|
141 | stop |
---|
142 | endif |
---|
143 | |
---|
144 | |
---|
145 | ! Load the Vtable into 'my_vtable' |
---|
146 | print *, 'Loading Vtable: ', vtable_path |
---|
147 | call vtable_load_by_name(vtable_path, my_vtable) |
---|
148 | print *, 'Vtable Initialized: ', my_vtable%initialized |
---|
149 | print *, 'Vtable num_entries: ', my_vtable%num_entries |
---|
150 | !!!!!!!!!!!!!!!!!!! VTABLE code |
---|
151 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
152 | |
---|
153 | |
---|
154 | |
---|
155 | |
---|
156 | |
---|
157 | ! |
---|
158 | ! OPENING OF DATA FILE (GRIB CODE) |
---|
159 | ! |
---|
160 | ifile=0 |
---|
161 | igrib=0 |
---|
162 | iret=0 |
---|
163 | |
---|
164 | 5 call grib_open_file(ifile,path(numpath+2*(l-1)+1) & |
---|
165 | (1:length(numpath+2*(l-1)+1))//trim(wfnamen(l,ifn)),'r',iret) |
---|
166 | if (iret.ne.GRIB_SUCCESS) then |
---|
167 | goto 999 ! ERROR DETECTED |
---|
168 | endif |
---|
169 | !turn on support for multi fields messages |
---|
170 | !call grib_multi_support_on |
---|
171 | |
---|
172 | gotGrib=0 |
---|
173 | ifield=0 |
---|
174 | 10 ifield=ifield+1 |
---|
175 | |
---|
176 | ! |
---|
177 | ! GET NEXT FIELDS |
---|
178 | ! |
---|
179 | call grib_new_from_file(ifile,igrib,iret) |
---|
180 | if (iret.eq.GRIB_END_OF_FILE) then |
---|
181 | goto 30 ! EOF DETECTED |
---|
182 | elseif (iret.ne.GRIB_SUCCESS) then |
---|
183 | goto 999 ! ERROR DETECTED |
---|
184 | endif |
---|
185 | |
---|
186 | |
---|
187 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
188 | !!!!!!!!!!!!!!!!!!! VTABLE code |
---|
189 | ! Get the fpname |
---|
190 | fpname = vtable_get_fpname(igrib, my_vtable) |
---|
191 | !print *, 'fpname: ', trim(fpname) |
---|
192 | |
---|
193 | |
---|
194 | !!!!!!!!!!!!!!!!!!! VTABLE code |
---|
195 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
196 | |
---|
197 | |
---|
198 | |
---|
199 | |
---|
200 | |
---|
201 | |
---|
202 | #if defined CTBTO |
---|
203 | conversion_factor=1.0 |
---|
204 | #endif |
---|
205 | |
---|
206 | !first see if we read GRIB1 or GRIB2 |
---|
207 | call grib_get_int(igrib,'editionNumber',gribVer,iret) |
---|
208 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
209 | |
---|
210 | if (gribVer.eq.1) then ! GRIB Edition 1 |
---|
211 | |
---|
212 | !print*,'GRiB Edition 1' |
---|
213 | !read the grib2 identifiers |
---|
214 | !call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret) |
---|
215 | !call grib_check(iret,gribFunction,gribErrorMsg) |
---|
216 | call grib_get_int(igrib,'level',current_grib_level,iret) |
---|
217 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
218 | |
---|
219 | else |
---|
220 | |
---|
221 | #if defined CTBTO |
---|
222 | call grib2check(igrib, fpname, conversion_factor) |
---|
223 | #endif |
---|
224 | |
---|
225 | call grib_get_int(igrib, 'level', current_grib_level, iret) |
---|
226 | call grib_check(iret, gribFunction, gribErrorMsg) |
---|
227 | |
---|
228 | endif |
---|
229 | |
---|
230 | |
---|
231 | |
---|
232 | !get the size and data of the values array |
---|
233 | if (trim(fpname) .ne. 'None') then |
---|
234 | call grib_get_real4_array(igrib,'values',zsec4,iret) |
---|
235 | #if defined CTBTO |
---|
236 | zsec4=zsec4/conversion_factor |
---|
237 | #endif |
---|
238 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
239 | endif |
---|
240 | |
---|
241 | !HSO get the required fields from section 2 in a gribex compatible manner |
---|
242 | if (ifield.eq.1) then |
---|
243 | call grib_get_int(igrib,'numberOfPointsAlongAParallel', & |
---|
244 | isec2(2),iret) |
---|
245 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
246 | call grib_get_int(igrib,'numberOfPointsAlongAMeridian', & |
---|
247 | isec2(3),iret) |
---|
248 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
249 | call grib_get_int(igrib,'numberOfVerticalCoordinateValues', & |
---|
250 | isec2(12),iret) |
---|
251 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
252 | !HSO get the size and data of the vertical coordinate array |
---|
253 | call grib_get_real4_array(igrib,'pv',zsec2,iret) |
---|
254 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
255 | |
---|
256 | nxn(l)=isec2(2) |
---|
257 | nyn(l)=isec2(3) |
---|
258 | nlev_ecn=isec2(12)/2-1 |
---|
259 | endif ! ifield |
---|
260 | |
---|
261 | !HSO get the second part of the grid dimensions only from GRiB1 messages |
---|
262 | #if defined CTBTO |
---|
263 | if (gotGrib.eq.0) then |
---|
264 | #else |
---|
265 | if ((gribVer.eq.1).and.(gotGrib.eq.0)) then |
---|
266 | #endif |
---|
267 | call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', & |
---|
268 | xaux1in,iret) |
---|
269 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
270 | call grib_get_real8(igrib,'longitudeOfLastGridPointInDegrees', & |
---|
271 | xaux2in,iret) |
---|
272 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
273 | call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', & |
---|
274 | yaux1in,iret) |
---|
275 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
276 | call grib_get_real8(igrib,'latitudeOfFirstGridPointInDegrees', & |
---|
277 | yaux2in,iret) |
---|
278 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
279 | xaux1=xaux1in |
---|
280 | xaux2=xaux2in |
---|
281 | yaux1=yaux1in |
---|
282 | yaux2=yaux2in |
---|
283 | if(xaux1.gt.180) xaux1=xaux1-360.0 |
---|
284 | if(xaux2.gt.180) xaux2=xaux2-360.0 |
---|
285 | if(xaux1.lt.-180) xaux1=xaux1+360.0 |
---|
286 | if(xaux2.lt.-180) xaux2=xaux2+360.0 |
---|
287 | if (xaux2.lt.xaux1) xaux2=xaux2+360. |
---|
288 | xlon0n(l)=xaux1 |
---|
289 | ylat0n(l)=yaux1 |
---|
290 | dxn(l)=(xaux2-xaux1)/real(nxn(l)-1) |
---|
291 | dyn(l)=(yaux2-yaux1)/real(nyn(l)-1) |
---|
292 | gotGrib=1 |
---|
293 | endif ! ifield.eq.1 |
---|
294 | |
---|
295 | k=current_grib_level |
---|
296 | if(trim(fpname) .eq. 'UU') iumax=max(iumax,nlev_ec-k+1) |
---|
297 | if(trim(fpname) .eq. 'ETADOT') iwmax=max(iwmax,nlev_ec-k+1) |
---|
298 | |
---|
299 | if(trim(fpname) .eq. 'ORO') then |
---|
300 | do j=0,nyn(l)-1 |
---|
301 | do i=0,nxn(l)-1 |
---|
302 | oron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga |
---|
303 | end do |
---|
304 | end do |
---|
305 | endif |
---|
306 | if(trim(fpname) .eq. 'LSM') then |
---|
307 | do j=0,nyn(l)-1 |
---|
308 | do i=0,nxn(l)-1 |
---|
309 | lsmn(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga |
---|
310 | end do |
---|
311 | end do |
---|
312 | endif |
---|
313 | if(trim(fpname) .eq. 'EXCESSORO') then |
---|
314 | do j=0,nyn(l)-1 |
---|
315 | do i=0,nxn(l)-1 |
---|
316 | excessoron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga |
---|
317 | end do |
---|
318 | end do |
---|
319 | endif |
---|
320 | |
---|
321 | call grib_release(igrib) |
---|
322 | goto 10 !! READ NEXT LEVEL OR PARAMETER |
---|
323 | ! |
---|
324 | ! CLOSING OF INPUT DATA FILE |
---|
325 | ! |
---|
326 | |
---|
327 | 30 call grib_close_file(ifile) |
---|
328 | |
---|
329 | !error message if no fields found with correct first longitude in it |
---|
330 | if (gotGrib.eq.0) then |
---|
331 | print*,'***ERROR: input file needs to contain GRiB1 formatted'// & |
---|
332 | 'messages' |
---|
333 | stop |
---|
334 | endif |
---|
335 | |
---|
336 | nuvzn=iumax |
---|
337 | nwzn=iwmax |
---|
338 | if(nuvzn.eq.nlev_ec) nwzn=nlev_ecn+1 |
---|
339 | |
---|
340 | if (nxn(l).gt.nxmaxn) then |
---|
341 | write(*,*) 'FLEXPART error: Too many grid points in x direction.' |
---|
342 | write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)' |
---|
343 | write(*,*) 'for nesting level ',l |
---|
344 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
345 | write(*,*) nxn(l),nxmaxn |
---|
346 | stop |
---|
347 | endif |
---|
348 | |
---|
349 | if (nyn(l).gt.nymaxn) then |
---|
350 | write(*,*) 'FLEXPART error: Too many grid points in y direction.' |
---|
351 | write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)' |
---|
352 | write(*,*) 'for nesting level ',l |
---|
353 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
354 | write(*,*) nyn(l),nymaxn |
---|
355 | stop |
---|
356 | endif |
---|
357 | |
---|
358 | if ((nuvzn.gt.nuvzmax).or.(nwzn.gt.nwzmax)) then |
---|
359 | write(*,*) 'FLEXPART error: Nested wind fields have too many'// & |
---|
360 | 'vertical levels.' |
---|
361 | write(*,*) 'Problem was encountered for nesting level ',l |
---|
362 | stop |
---|
363 | endif |
---|
364 | |
---|
365 | |
---|
366 | ! Output of grid info |
---|
367 | !******************** |
---|
368 | |
---|
369 | write(*,'(a,i2)') 'Nested domain #: ',l |
---|
370 | write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Longitude range: ', & |
---|
371 | xlon0n(l),' to ',xlon0n(l)+(nxn(l)-1)*dxn(l), & |
---|
372 | ' Grid distance: ',dxn(l) |
---|
373 | write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Latitude range: ', & |
---|
374 | ylat0n(l),' to ',ylat0n(l)+(nyn(l)-1)*dyn(l), & |
---|
375 | ' Grid distance: ',dyn(l) |
---|
376 | write(*,*) |
---|
377 | |
---|
378 | ! Determine, how much the resolutions in the nests are enhanced as |
---|
379 | ! compared to the mother grid |
---|
380 | !***************************************************************** |
---|
381 | |
---|
382 | xresoln(l)=dx/dxn(l) |
---|
383 | yresoln(l)=dy/dyn(l) |
---|
384 | |
---|
385 | ! Determine the mother grid coordinates of the corner points of the |
---|
386 | ! nested grids |
---|
387 | ! Convert first to geographical coordinates, then to grid coordinates |
---|
388 | !******************************************************************** |
---|
389 | |
---|
390 | xaux1=xlon0n(l) |
---|
391 | xaux2=xlon0n(l)+real(nxn(l)-1)*dxn(l) |
---|
392 | yaux1=ylat0n(l) |
---|
393 | yaux2=ylat0n(l)+real(nyn(l)-1)*dyn(l) |
---|
394 | |
---|
395 | xln(l)=(xaux1-xlon0)/dx |
---|
396 | xrn(l)=(xaux2-xlon0)/dx |
---|
397 | yln(l)=(yaux1-ylat0)/dy |
---|
398 | yrn(l)=(yaux2-ylat0)/dy |
---|
399 | |
---|
400 | |
---|
401 | if ((xln(l).lt.0.).or.(yln(l).lt.0.).or. & |
---|
402 | (xrn(l).gt.real(nxmin1)).or.(yrn(l).gt.real(nymin1))) then |
---|
403 | write(*,*) 'Nested domain does not fit into mother domain' |
---|
404 | write(*,*) 'For global mother domain fields, you can shift' |
---|
405 | write(*,*) 'shift the mother domain into x-direction' |
---|
406 | write(*,*) 'by setting nxshift (file par_mod) to a' |
---|
407 | write(*,*) 'positive value. Execution is terminated.' |
---|
408 | stop |
---|
409 | endif |
---|
410 | |
---|
411 | |
---|
412 | ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL |
---|
413 | ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM |
---|
414 | |
---|
415 | numskip=nlev_ecn-nuvzn ! number of ecmwf model layers not used by FLEXPART |
---|
416 | do i=1,nwzn |
---|
417 | j=numskip+i |
---|
418 | k=nlev_ecn+1+numskip+i |
---|
419 | akmn(nwzn-i+1)=zsec2(j) |
---|
420 | bkmn(nwzn-i+1)=zsec2(k) |
---|
421 | end do |
---|
422 | |
---|
423 | ! |
---|
424 | ! CALCULATION OF AKZ, BKZ |
---|
425 | ! AKZ,BKZ: model discretization parameters at the center of each model |
---|
426 | ! layer |
---|
427 | ! |
---|
428 | ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0, |
---|
429 | ! i.e. ground level |
---|
430 | !***************************************************************************** |
---|
431 | |
---|
432 | akzn(1)=0. |
---|
433 | bkzn(1)=1.0 |
---|
434 | do i=1,nuvzn |
---|
435 | akzn(i+1)=0.5*(akmn(i+1)+akmn(i)) |
---|
436 | bkzn(i+1)=0.5*(bkmn(i+1)+bkmn(i)) |
---|
437 | end do |
---|
438 | nuvzn=nuvzn+1 |
---|
439 | |
---|
440 | ! Check, whether the heights of the model levels of the nested |
---|
441 | ! wind fields are consistent with those of the mother domain. |
---|
442 | ! If not, terminate model run. |
---|
443 | !************************************************************* |
---|
444 | |
---|
445 | do i=1,nuvz |
---|
446 | if ((akzn(i).ne.akz(i)).or.(bkzn(i).ne.bkz(i))) then |
---|
447 | write(*,*) 'FLEXPART error: The wind fields of nesting level',l |
---|
448 | write(*,*) 'are not consistent with the mother domain:' |
---|
449 | write(*,*) 'Differences in vertical levels detected.' |
---|
450 | stop |
---|
451 | endif |
---|
452 | end do |
---|
453 | |
---|
454 | do i=1,nwz |
---|
455 | if ((akmn(i).ne.akm(i)).or.(bkmn(i).ne.bkm(i))) then |
---|
456 | write(*,*) 'FLEXPART error: The wind fields of nesting level',l |
---|
457 | write(*,*) 'are not consistent with the mother domain:' |
---|
458 | write(*,*) 'Differences in vertical levels detected.' |
---|
459 | stop |
---|
460 | endif |
---|
461 | end do |
---|
462 | |
---|
463 | end do |
---|
464 | |
---|
465 | return |
---|
466 | |
---|
467 | 999 write(*,*) |
---|
468 | write(*,*) ' ###########################################'// & |
---|
469 | '###### ' |
---|
470 | write(*,*) ' FLEXPART SUBROUTINE GRIDCHECK:' |
---|
471 | write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfnamen(l,ifn) |
---|
472 | write(*,*) ' FOR NESTING LEVEL ',k |
---|
473 | write(*,*) ' ###########################################'// & |
---|
474 | '###### ' |
---|
475 | stop |
---|
476 | |
---|
477 | end subroutine gridcheck_nests |
---|