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