CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
gridcheck.F90
Go to the documentation of this file.
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_ecmwf
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: 11/01/2008, Harald Sodemann, GRIB1/2 input with *
36  ! ECMWF grib_api *
37  ! CHANGE: 03/12/2008, Harald Sodemann, update to f90 with *
38  ! ECMWF grib_api *
39  ! *
40  !**********************************************************************
41  ! *
42  ! DESCRIPTION: *
43  ! *
44  ! THIS SUBROUTINE DETERMINES THE GRID SPECIFICATIONS (LOWER LEFT *
45  ! LONGITUDE, LOWER LEFT LATITUDE, NUMBER OF GRID POINTS, GRID DIST- *
46  ! ANCE AND VERTICAL DISCRETIZATION OF THE ECMWF MODEL) FROM THE *
47  ! GRIB HEADER OF THE FIRST INPUT FILE. THE CONSISTANCY (NO CHANGES *
48  ! WITHIN ONE FLEXPART RUN) IS CHECKED IN THE ROUTINE "READWIND" AT *
49  ! ANY CALL. *
50  ! *
51  ! XLON0 geographical longitude of lower left gridpoint *
52  ! YLAT0 geographical latitude of lower left gridpoint *
53  ! NX number of grid points x-direction *
54  ! NY number of grid points y-direction *
55  ! DX grid distance x-direction *
56  ! DY grid distance y-direction *
57  ! NUVZ number of grid points for horizontal wind *
58  ! components in z direction *
59  ! NWZ number of grid points for vertical wind *
60  ! component in z direction *
61  ! sizesouth, sizenorth give the map scale (i.e. number of virtual grid*
62  ! points of the polar stereographic grid): *
63  ! used to check the CFL criterion *
64  ! UVHEIGHT(1)- heights of gridpoints where u and v are *
65  ! UVHEIGHT(NUVZ) given *
66  ! WHEIGHT(1)- heights of gridpoints where w is given *
67  ! WHEIGHT(NWZ) *
68  ! *
69  !**********************************************************************
70 
71  use grib_api
72  use par_mod
73  use com_mod
74  use conv_mod
75  use cmapf_mod, only: stlmbr,stcm2p
76 
77  implicit none
78 
79  !HSO parameters for grib_api
80  integer :: ifile
81  integer :: iret
82  integer :: igrib
83  integer :: gotgrid
84  real(kind=4) :: xaux1,xaux2,yaux1,yaux2
85  real(kind=8) :: xaux1in,xaux2in,yaux1in,yaux2in
86  integer :: gribver,parcat,parnum,typsurf,valsurf,discipl
87 #if defined WITH_CTBTO_PATCHES
88  integer :: parid
89 #endif
90  !HSO end
91  integer :: ix,jy,i,ifn,ifield,j,k,iumax,iwmax,numskip
92  real :: sizesouth,sizenorth,xauxa,pint
93 
94  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
95 
96  ! dimension of isec2 at least (22+n), where n is the number of parallels or
97  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
98 
99  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
100  ! coordinate parameters
101 
102  integer :: isec1(56),isec2(22+nxmax+nymax)
103  real(kind=4) :: zsec2(60+2*nuvzmax),zsec4(jpunp)
104 #if defined WITH_CTBTO_PATCHES
105  real(kind=4) :: conversion_factor
106 #endif
107  character(len=1) :: opt
108 
109  !HSO grib api error messages
110  character(len=24) :: griberrormsg = 'Error reading grib file'
111  character(len=20) :: gribfunction = 'gridcheck'
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  ! OPENING OF DATA FILE (GRIB CODE)
123  !
124 5 call grib_open_file(ifile,path(3)(1:length(3)) &
125  //trim(wfname(ifn)),'r',iret)
126  if (iret.ne.grib_success) then
127  goto 999 ! ERROR DETECTED
128  endif
129  !turn on support for multi fields messages
130  !call grib_multi_support_on
131 
132  gotgrid=0
133  ifield=0
134 10 ifield=ifield+1
135 
136  !
137  ! GET NEXT FIELDS
138  !
139  call grib_new_from_file(ifile,igrib,iret)
140  if (iret.eq.grib_end_of_file ) then
141  goto 30 ! EOF DETECTED
142  elseif (iret.ne.grib_success) then
143  goto 999 ! ERROR DETECTED
144  endif
145 
146 #if defined WITH_CTBTO_PATCHES
147  conversion_factor=1.0
148 #endif
149 
150  !first see if we read GRIB1 or GRIB2
151  call grib_get_int(igrib,'editionNumber',gribver,iret)
152  call grib_check(iret,gribfunction,griberrormsg)
153 
154  if (gribver.eq.1) then ! GRIB Edition 1
155 
156  !print*,'GRiB Edition 1'
157  !read the grib2 identifiers
158  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
159  call grib_check(iret,gribfunction,griberrormsg)
160  call grib_get_int(igrib,'level',isec1(8),iret)
161  call grib_check(iret,gribfunction,griberrormsg)
162 
163  !change code for etadot to code for omega
164 #if defined WITH_CTBTO_PATCHES
165  if (isec1(6).eq.77.or.isec1(6).eq.87) then
166 #else
167  if (isec1(6).eq.77) then
168 #endif
169  isec1(6)=135
170  endif
171 
172  !print*,isec1(6),isec1(8)
173 
174  else
175 
176 #if defined WITH_CTBTO_PATCHES
177  call grib2check(igrib,isec1,conversion_factor)
178 #else
179  !print*,'GRiB Edition 2'
180  !read the grib2 identifiers
181  call grib_get_int(igrib,'discipline',discipl,iret)
182  call grib_check(iret,gribfunction,griberrormsg)
183  call grib_get_int(igrib,'parameterCategory',parcat,iret)
184  call grib_check(iret,gribfunction,griberrormsg)
185  call grib_get_int(igrib,'parameterNumber',parnum,iret)
186  call grib_check(iret,gribfunction,griberrormsg)
187  call grib_get_int(igrib,'typeOfFirstFixedSurface',typsurf,iret)
188  call grib_check(iret,gribfunction,griberrormsg)
189  call grib_get_int(igrib,'level',valsurf,iret)
190  call grib_check(iret,gribfunction,griberrormsg)
191 
192  !print*,discipl,parCat,parNum,typSurf,valSurf
193 
194  !convert to grib1 identifiers
195  isec1(6)=-1
196  isec1(7)=-1
197  isec1(8)=-1
198  isec1(8)=valsurf ! level
199  if ((parcat.eq.0).and.(parnum.eq.0).and.(typsurf.eq.105)) then ! T
200  isec1(6)=130 ! indicatorOfParameter
201  elseif ((parcat.eq.2).and.(parnum.eq.2).and.(typsurf.eq.105)) then ! U
202  isec1(6)=131 ! indicatorOfParameter
203  elseif ((parcat.eq.2).and.(parnum.eq.3).and.(typsurf.eq.105)) then ! V
204  isec1(6)=132 ! indicatorOfParameter
205  elseif ((parcat.eq.1).and.(parnum.eq.0).and.(typsurf.eq.105)) then ! Q
206  isec1(6)=133 ! indicatorOfParameter
207  elseif ((parcat.eq.3).and.(parnum.eq.0).and.(typsurf.eq.1)) then !SP
208  isec1(6)=134 ! indicatorOfParameter
209  elseif ((parcat.eq.2).and.(parnum.eq.32)) then ! W, actually eta dot
210  isec1(6)=135 ! indicatorOfParameter
211  elseif ((parcat.eq.3).and.(parnum.eq.0).and.(typsurf.eq.101)) then !SLP
212  isec1(6)=151 ! indicatorOfParameter
213  elseif ((parcat.eq.2).and.(parnum.eq.2).and.(typsurf.eq.103)) then ! 10U
214  isec1(6)=165 ! indicatorOfParameter
215  elseif ((parcat.eq.2).and.(parnum.eq.3).and.(typsurf.eq.103)) then ! 10V
216  isec1(6)=166 ! indicatorOfParameter
217  elseif ((parcat.eq.0).and.(parnum.eq.0).and.(typsurf.eq.103)) then ! 2T
218  isec1(6)=167 ! indicatorOfParameter
219  elseif ((parcat.eq.0).and.(parnum.eq.6).and.(typsurf.eq.103)) then ! 2D
220  isec1(6)=168 ! indicatorOfParameter
221  elseif ((parcat.eq.1).and.(parnum.eq.11).and.(typsurf.eq.1)) then ! SD
222  isec1(6)=141 ! indicatorOfParameter
223  elseif ((parcat.eq.6).and.(parnum.eq.1)) then ! CC
224  isec1(6)=164 ! indicatorOfParameter
225  elseif ((parcat.eq.1).and.(parnum.eq.9)) then ! LSP
226  isec1(6)=142 ! indicatorOfParameter
227  elseif ((parcat.eq.1).and.(parnum.eq.10)) then ! CP
228  isec1(6)=143 ! indicatorOfParameter
229  elseif ((parcat.eq.0).and.(parnum.eq.11).and.(typsurf.eq.1)) then ! SHF
230  isec1(6)=146 ! indicatorOfParameter
231  elseif ((parcat.eq.4).and.(parnum.eq.9).and.(typsurf.eq.1)) then ! SR
232  isec1(6)=176 ! indicatorOfParameter
233  elseif ((parcat.eq.2).and.(parnum.eq.17)) then ! EWSS
234  isec1(6)=180 ! indicatorOfParameter
235  elseif ((parcat.eq.2).and.(parnum.eq.18)) then ! NSSS
236  isec1(6)=181 ! indicatorOfParameter
237  elseif ((parcat.eq.3).and.(parnum.eq.4)) then ! ORO
238  isec1(6)=129 ! indicatorOfParameter
239  elseif ((parcat.eq.3).and.(parnum.eq.7)) then ! SDO
240  isec1(6)=160 ! indicatorOfParameter
241  elseif ((discipl.eq.2).and.(parcat.eq.0).and.(parnum.eq.0).and. &
242  (typsurf.eq.1)) then ! LSM
243  isec1(6)=172 ! indicatorOfParameter
244  else
245  print*,'***ERROR: undefined GRiB2 message found!',discipl, &
246  parcat,parnum,typsurf
247  endif
248 #endif
249 
250  endif
251 
252  !get the size and data of the values array
253  if (isec1(6).ne.-1) then
254  call grib_get_real4_array(igrib,'values',zsec4,iret)
255 #if defined WITH_CTBTO_PATCHES
256  zsec4=zsec4/conversion_factor
257 #endif
258  call grib_check(iret,gribfunction,griberrormsg)
259  endif
260 
261  if (ifield.eq.1) then
262 
263  !HSO get the required fields from section 2 in a gribex compatible manner
264  call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
265  isec2(2),iret)
266  call grib_check(iret,gribfunction,griberrormsg)
267  call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
268  isec2(3),iret)
269  call grib_check(iret,gribfunction,griberrormsg)
270  call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
271  xaux1in,iret)
272  call grib_check(iret,gribfunction,griberrormsg)
273  call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
274  isec2(12),iret)
275  call grib_check(iret,gribfunction,griberrormsg)
276 
277  ! get the size and data of the vertical coordinate array
278  call grib_get_real4_array(igrib,'pv',zsec2,iret)
279  call grib_check(iret,gribfunction,griberrormsg)
280 
281  nxfield=isec2(2)
282  ny=isec2(3)
283  nlev_ec=isec2(12)/2-1
284  endif
285 
286  !HSO get the second part of the grid dimensions only from GRiB1 messages
287 #if defined WITH_CTBTO_PATCHES
288  if (gotgrid.eq.0) then
289 #else
290  if ((gribver.eq.1).and.(gotgrid.eq.0)) then
291 #endif
292  call grib_get_real8(igrib,'longitudeOfLastGridPointInDegrees', &
293  xaux2in,iret)
294  call grib_check(iret,gribfunction,griberrormsg)
295  call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
296  yaux1in,iret)
297  call grib_check(iret,gribfunction,griberrormsg)
298  call grib_get_real8(igrib,'latitudeOfFirstGridPointInDegrees', &
299  yaux2in,iret)
300  call grib_check(iret,gribfunction,griberrormsg)
301  xaux1=xaux1in
302  xaux2=xaux2in
303  yaux1=yaux1in
304  yaux2=yaux2in
305 #if defined WITH_CTBTO_PATCHES
306  if (xaux1.ge.180.) xaux1=xaux1-360.0
307 #else
308  if (xaux1.gt.180.) xaux1=xaux1-360.0
309 #endif
310  if (xaux2.gt.180.) xaux2=xaux2-360.0
311  if (xaux1.lt.-180.) xaux1=xaux1+360.0
312  if (xaux2.lt.-180.) xaux2=xaux2+360.0
313  if (xaux2.lt.xaux1) xaux2=xaux2+360.0
314  xlon0=xaux1
315  ylat0=yaux1
316  dx=(xaux2-xaux1)/real(nxfield-1)
317  dy=(yaux2-yaux1)/real(ny-1)
318  dxconst=180./(dx*r_earth*pi)
319  dyconst=180./(dy*r_earth*pi)
320  gotgrid=1
321 
322  ! Check whether fields are global
323  ! If they contain the poles, specify polar stereographic map
324  ! projections using the stlmbr- and stcm2p-calls
325  !***********************************************************
326 
327  xauxa=abs(xaux2+dx-360.-xaux1)
328  if (xauxa.lt.0.001) then
329  nx=nxfield+1 ! field is cyclic
330  xglobal=.true.
331  if (abs(nxshift).ge.nx) &
332  stop 'nxshift in file par_mod is too large'
333  xlon0=xlon0+real(nxshift)*dx
334  else
335  nx=nxfield
336  xglobal=.false.
337  if (nxshift.ne.0) &
338  stop 'nxshift (par_mod) must be zero for non-global domain'
339  endif
340  nxmin1=nx-1
341  nymin1=ny-1
342  if (xlon0.gt.180.) xlon0=xlon0-360.
343  xauxa=abs(yaux1+90.)
344  if (xglobal.and.xauxa.lt.0.001) then
345  sglobal=.true. ! field contains south pole
346  ! Enhance the map scale by factor 3 (*2=6) compared to north-south
347  ! map scale
348  sizesouth=6.*(switchsouth+90.)/dy
349  call stlmbr(southpolemap,-90.,0.)
350  call stcm2p(southpolemap,0.,0.,switchsouth,0.,sizesouth, &
351  sizesouth,switchsouth,180.)
352  switchsouthg=(switchsouth-ylat0)/dy
353  else
354  sglobal=.false.
355  switchsouthg=999999.
356  endif
357  xauxa=abs(yaux2-90.)
358  if (xglobal.and.xauxa.lt.0.001) then
359  nglobal=.true. ! field contains north pole
360  ! Enhance the map scale by factor 3 (*2=6) compared to north-south
361  ! map scale
362  sizenorth=6.*(90.-switchnorth)/dy
363  call stlmbr(northpolemap,90.,0.)
364  call stcm2p(northpolemap,0.,0.,switchnorth,0.,sizenorth, &
365  sizenorth,switchnorth,180.)
366  switchnorthg=(switchnorth-ylat0)/dy
367  else
368  nglobal=.false.
369  switchnorthg=999999.
370  endif
371  if (nxshift.lt.0) &
372  stop 'nxshift (par_mod) must not be negative'
373  if (nxshift.ge.nxfield) stop 'nxshift (par_mod) too large'
374  endif ! gotGrid
375 
376  k=isec1(8)
377  if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
378  if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
379 
380  if(isec1(6).eq.129) then
381  do jy=0,ny-1
382  do ix=0,nxfield-1
383  oro(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)/ga
384  end do
385  end do
386  endif
387  if(isec1(6).eq.172) then
388  do jy=0,ny-1
389  do ix=0,nxfield-1
390  lsm(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)
391  end do
392  end do
393  endif
394  if(isec1(6).eq.160) then
395  do jy=0,ny-1
396  do ix=0,nxfield-1
397  excessoro(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)
398  end do
399  end do
400  endif
401 
402  call grib_release(igrib)
403  goto 10 !! READ NEXT LEVEL OR PARAMETER
404  !
405  ! CLOSING OF INPUT DATA FILE
406  !
407 
408 30 call grib_close_file(ifile)
409 
410  !error message if no fields found with correct first longitude in it
411  if (gotgrid.eq.0) then
412  print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
413  'messages'
414  stop
415  endif
416 
417  nuvz=iumax
418  nwz =iwmax
419  if(nuvz.eq.nlev_ec) nwz=nlev_ec+1
420 
421  if (nx.gt.nxmax) then
422  write(*,*) 'FLEXPART error: Too many grid points in x direction.'
423  write(*,*) 'Reduce resolution of wind fields.'
424  write(*,*) 'Or change parameter settings in file par_mod.'
425  write(*,*) nx,nxmax
426  stop
427  endif
428 
429  if (ny.gt.nymax) then
430  write(*,*) 'FLEXPART error: Too many grid points in y direction.'
431  write(*,*) 'Reduce resolution of wind fields.'
432  write(*,*) 'Or change parameter settings in file par_mod.'
433  write(*,*) ny,nymax
434  stop
435  endif
436 
437  if (nuvz+1.gt.nuvzmax) then
438  write(*,*) 'FLEXPART error: Too many u,v grid points in z '// &
439  'direction.'
440  write(*,*) 'Reduce resolution of wind fields.'
441  write(*,*) 'Or change parameter settings in file par_mod.'
442  write(*,*) nuvz+1,nuvzmax
443  stop
444  endif
445 
446  if (nwz.gt.nwzmax) then
447  write(*,*) 'FLEXPART error: Too many w grid points in z '// &
448  'direction.'
449  write(*,*) 'Reduce resolution of wind fields.'
450  write(*,*) 'Or change parameter settings in file par_mod.'
451  write(*,*) nwz,nwzmax
452  stop
453  endif
454 
455  ! If desired, shift all grids by nxshift grid cells
456  !**************************************************
457 
458  if (xglobal) then
459  call shift_field_0(oro,nxfield,ny)
460  call shift_field_0(lsm,nxfield,ny)
461  call shift_field_0(excessoro,nxfield,ny)
462  endif
463 
464  ! Output of grid info
465  !********************
466 
467  write(*,*)
468  write(*,*)
469  write(*,'(a,2i7)') '# of vertical levels in ECMWF data: ', &
470  nuvz+1,nwz
471  write(*,*)
472  write(*,'(a)') 'Mother domain:'
473  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Longitude range: ', &
474  xlon0,' to ',xlon0+(nx-1)*dx,' Grid distance: ',dx
475  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Latitude range: ', &
476  ylat0,' to ',ylat0+(ny-1)*dy,' Grid distance: ',dy
477  write(*,*)
478 
479 
480  ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL
481  ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM
482 
483  numskip=nlev_ec-nuvz ! number of ecmwf model layers not used
484  ! by trajectory model
485  !do 8940 i=1,244
486  ! write (*,*) 'zsec2:',i,ifield,zsec2(i),numskip
487  !940 continue
488  ! stop
489  ! SEC SEC SEC
490  ! for unknown reason zsec 1 to 10 is filled in this version
491  ! compared to the old one
492  ! SEC SEC SE
493  do i=1,nwz
494  j=numskip+i
495  k=nlev_ec+1+numskip+i
496  akm(nwz-i+1)=zsec2(j)
497  ! write (*,*) 'ifield:',ifield,k,j,zsec2(10+j)
498  bkm(nwz-i+1)=zsec2(k)
499  end do
500 
501  !
502  ! CALCULATION OF AKZ, BKZ
503  ! AKZ,BKZ: model discretization parameters at the center of each model
504  ! layer
505  !
506  ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0,
507  ! i.e. ground level
508  !*****************************************************************************
509 
510  akz(1)=0.
511  bkz(1)=1.0
512  do i=1,nuvz
513  akz(i+1)=0.5*(akm(i+1)+akm(i))
514  bkz(i+1)=0.5*(bkm(i+1)+bkm(i))
515  end do
516  nuvz=nuvz+1
517 
518  ! NOTE: In FLEXPART versions up to 4.0, the number of model levels was doubled
519  ! upon the transformation to z levels. In order to save computer memory, this is
520  ! not done anymore in the standard version. However, this option can still be
521  ! switched on by replacing the following lines with those below, that are
522  ! currently commented out. For this, similar changes are necessary in
523  ! verttransform.f and verttranform_nests.f
524  !*****************************************************************************
525 
526  nz=nuvz
527  if (nz.gt.nzmax) stop 'nzmax too small'
528  do i=1,nuvz
529  aknew(i)=akz(i)
530  bknew(i)=bkz(i)
531  end do
532 
533  ! Switch on following lines to use doubled vertical resolution
534  !*************************************************************
535  !nz=nuvz+nwz-1
536  !if (nz.gt.nzmax) stop 'nzmax too small'
537  !do 100 i=1,nwz
538  ! aknew(2*(i-1)+1)=akm(i)
539  !00 bknew(2*(i-1)+1)=bkm(i)
540  !do 110 i=2,nuvz
541  ! aknew(2*(i-1))=akz(i)
542  !10 bknew(2*(i-1))=bkz(i)
543  ! End doubled vertical resolution
544 
545 
546  ! Determine the uppermost level for which the convection scheme shall be applied
547  ! by assuming that there is no convection above 50 hPa (for standard SLP)
548  !*****************************************************************************
549 
550  do i=1,nuvz-2
551  pint=akz(i)+bkz(i)*101325.
552  if (pint.lt.5000.) goto 96
553  end do
554 96 nconvlev=i
555  if (nconvlev.gt.nconvlevmax-1) then
556  nconvlev=nconvlevmax-1
557  write(*,*) 'Attention, convection only calculated up to ', &
558  akz(nconvlev)+bkz(nconvlev)*1013.25,' hPa'
559  endif
560 
561  return
562 
563 999 write(*,*)
564  write(*,*) ' ###########################################'// &
565  '###### '
566  write(*,*) ' TRAJECTORY MODEL SUBROUTINE GRIDCHECK:'
567  write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfname(ifn)
568  write(*,*) ' ###########################################'// &
569  '###### '
570  write(*,*)
571  write(*,'(a)') '!!! PLEASE INSERT A NEW CD-ROM AND !!!'
572  write(*,'(a)') '!!! PRESS ANY KEY TO CONTINUE... !!!'
573  write(*,'(a)') '!!! ...OR TERMINATE FLEXPART PRESSING!!!'
574  write(*,'(a)') '!!! THE "X" KEY... !!!'
575  write(*,*)
576  read(*,'(a)') opt
577  if(opt.eq.'X') then
578  stop
579  else
580  goto 5
581  endif
582 
583 end subroutine gridcheck_ecmwf
subroutine, public stcm2p(strcmp, x1, y1, xlat1, xlong1, x2, y2, xlat2, xlong2)
Definition: cmapf_mod.f90:621
subroutine gridcheck_ecmwf
Definition: gridcheck.F90:22
subroutine, public stlmbr(strcmp, tnglat, xlong)
Definition: cmapf_mod.f90:802
subroutine shift_field_0(field, nxf, nyf)
subroutine grib2check(igrib, isec1, conversion_factor)
Definition: grib2check.F90:1