28   real(dp)                     :: plat,plon,pazi=0.0
    31   namelist /regional_grid_nml/ plat,plon,pazi,delx,dely,lx,ly
    33   real(dp),
parameter           :: re=6371200.0
    34   real(dp),
parameter           :: lam=0.8
    36   integer                      :: nxh,nyh, nx,ny, nxm,nym
    39   real(dp),
dimension(:,:),
allocatable:: glat,glon
    40   real(dp),
dimension(:,:),
allocatable:: garea
    41   real(dp),
dimension(:,:),
allocatable:: dx,dy,angle_dx,angle_dy
    43   character(len=256)           :: nml_fname
    47   integer                      :: string_dimid, nxp_dimid, nyp_dimid, nx_dimid, ny_dimid
    48   integer                      :: tile_varid, x_varid, y_varid, area_varid
    49   integer                      :: dx_varid, dy_varid, angle_dx_varid, angle_dy_varid
    50   integer, 
dimension(2)        :: dimids
    52   real(dp)                     :: a,k,m_arcx,m_arcy,q
    53   real(dp)                     :: m_delx,m_dely, delxre,delyre
    58   if (command_argument_count() == 1) 
then    59     call get_command_argument(1, nml_fname)
    61     nml_fname = 
"regional_grid.nml"    64   open(10,file=trim(nml_fname),status=
"old",action=
"read")
    65   read(10,nml=regional_grid_nml)
    75   allocate(glat(0:nx,0:ny))
    76   allocate(glon(0:nx,0:ny))
    77   allocate(garea(0:nxm,0:nym))
    79   allocate(dx(0:nxm,0:ny))
    80   allocate(dy(0:nx,0:nym))
    81   allocate(angle_dx(0:nx,0:ny))
    82   allocate(angle_dy(0:nx,0:ny))
    86   print
'("arcx, arcy ",f8.4,f8.4)',arcx,arcy
    87   call bestesg_geo(lam,arcx,arcy, a,k,m_arcx,m_arcy,q,ff)
    88   if(ff)stop 
'Failure flag returned from get_bestesg'    89   print
'("For lam=",f8.2," the best [smallest possible]")',lam
    90   print
'("optimality criterion, Q, for this domain: ",e13.6 )',q
    91   print
'("The corresponding optimal A and K:",2(1x,f8.4))',a,k
    92   print
'("The corresponding m_arcx,y:",2(1x,e20.13))',m_arcx,m_arcy
    96   print
'("x and y central grid resolution in map units:",2(1x,e12.5))',m_delx,m_dely
    98   print
'("Get additional diagnostics from hgrid_ak_rr")'   103   call hgrid_ak(lx,ly, nx,ny, a,k, plat*
dtor,plon*
dtor,pazi*
dtor, re,delxre,delyre, &
   104                 glat,glon,garea, dx,dy,angle_dx,angle_dy, ff)
   105   if(ff)stop 
'Failure flag raised in hgrid routine'   109   where (glon < 0.0) glon = glon + 360.0
   111   call check( nf90_create(
"regional_grid.nc", nf90_64bit_offset, ncid) )
   112   call check( nf90_def_dim(ncid, 
"string", 255, string_dimid) )
   113   call check( nf90_def_dim(ncid, 
"nx", nx, nx_dimid) )
   114   call check( nf90_def_dim(ncid, 
"ny", ny, ny_dimid) )
   115   call check( nf90_def_dim(ncid, 
"nxp", nx+1, nxp_dimid) )
   116   call check( nf90_def_dim(ncid, 
"nyp", ny+1, nyp_dimid) )
   118   call check( nf90_def_var(ncid, 
"tile", nf90_char, [string_dimid], tile_varid) )
   119   call check( nf90_put_att(ncid, tile_varid, 
"standard_name", 
"grid_tile_spec") )
   121   dimids = (/ nxp_dimid, nyp_dimid /)
   122   call check( nf90_def_var(ncid, 
"x", nf90_double, dimids, x_varid) )
   123   call check( nf90_put_att(ncid, x_varid, 
"standard_name", 
"geographic_longitude") )
   124   call check( nf90_put_att(ncid, x_varid, 
"units", 
"degree_east") )
   125   call check( nf90_put_att(ncid, x_varid, 
"hstagger", 
"C") )
   126   call check( nf90_def_var(ncid, 
"y", nf90_double, dimids, y_varid) )
   127   call check( nf90_put_att(ncid, y_varid, 
"standard_name", 
"geographic_latitude") )
   128   call check( nf90_put_att(ncid, y_varid, 
"units", 
"degree_north") )
   129   call check( nf90_put_att(ncid, y_varid, 
"hstagger", 
"C") )
   130   dimids = (/ nx_dimid, ny_dimid /)
   131   call check( nf90_def_var(ncid, 
"area", nf90_double, dimids, area_varid) )
   132   call check( nf90_put_att(ncid, area_varid, 
"standard_name", 
"grid_cell_area") )
   133   call check( nf90_put_att(ncid, area_varid, 
"units", 
"m2") )
   134   call check( nf90_put_att(ncid, area_varid, 
"hstagger", 
"H") )
   136   dimids = (/ nx_dimid, nyp_dimid /)
   137   call check( nf90_def_var(ncid, 
"dx", nf90_double, dimids, dx_varid) )
   138   call check( nf90_put_att(ncid, dx_varid, 
"standard_name", 
"dx") )
   139   call check( nf90_put_att(ncid, dx_varid, 
"units", 
"m") )
   140   call check( nf90_put_att(ncid, dx_varid, 
"hstagger", 
"H") )
   142   dimids = (/ nxp_dimid, ny_dimid /)
   143   call check( nf90_def_var(ncid, 
"dy", nf90_double, dimids, dy_varid) )
   144   call check( nf90_put_att(ncid, dy_varid, 
"standard_name", 
"dy") )
   145   call check( nf90_put_att(ncid, dy_varid, 
"units", 
"m") )
   146   call check( nf90_put_att(ncid, dy_varid, 
"hstagger", 
"H") )
   148   dimids = (/ nxp_dimid, nyp_dimid /)
   149   call check( nf90_def_var(ncid, 
"angle_dx", nf90_double, dimids, angle_dx_varid) )
   150   call check( nf90_put_att(ncid, angle_dx_varid, 
"standard_name", 
"angle_dx") )
   151   call check( nf90_put_att(ncid, angle_dx_varid, 
"units", 
"deg") )
   152   call check( nf90_put_att(ncid, angle_dx_varid, 
"hstagger", 
"C") )
   153   call check( nf90_def_var(ncid, 
"angle_dy", nf90_double, dimids, angle_dy_varid) )
   154   call check( nf90_put_att(ncid, angle_dy_varid, 
"standard_name", 
"angle_dy") )
   155   call check( nf90_put_att(ncid, angle_dy_varid, 
"units", 
"deg") )
   156   call check( nf90_put_att(ncid, angle_dy_varid, 
"hstagger", 
"C") )
   158   call check( nf90_put_att(ncid, nf90_global, 
"history", 
"gnomonic_ed") )
   159   call check( nf90_put_att(ncid, nf90_global, 
"source", 
"FV3GFS") )
   160   call check( nf90_put_att(ncid, nf90_global, 
"grid", 
"akappa") )
   161   call check( nf90_put_att(ncid, nf90_global, 
"plat", plat) )
   162   call check( nf90_put_att(ncid, nf90_global, 
"plon", plon) )
   163   call check( nf90_put_att(ncid, nf90_global, 
"pazi", pazi) )
   164   call check( nf90_put_att(ncid, nf90_global, 
"delx", m_delx*
rtod) )
   165   call check( nf90_put_att(ncid, nf90_global, 
"dely", m_dely*
rtod) )
   166   call check( nf90_put_att(ncid, nf90_global, 
"lx", lx) )
   167   call check( nf90_put_att(ncid, nf90_global, 
"ly", ly) )
   168   call check( nf90_put_att(ncid, nf90_global, 
"a", a) )
   169   call check( nf90_put_att(ncid, nf90_global, 
"k", k) )
   171   call check( nf90_enddef(ncid) )
   173   call check( nf90_put_var(ncid, tile_varid, 
"tile7") )
   174   call check( nf90_put_var(ncid, x_varid, glon) )
   175   call check( nf90_put_var(ncid, y_varid, glat) )
   176   call check( nf90_put_var(ncid, area_varid, garea) )
   177   call check( nf90_put_var(ncid, dx_varid, dx) )
   178   call check( nf90_put_var(ncid, dy_varid, dy) )
   179   call check( nf90_put_var(ncid, angle_dx_varid, angle_dx) )
   180   call check( nf90_put_var(ncid, angle_dy_varid, angle_dy) )
   182   call check( nf90_close(ncid) )
   190 subroutine check(status)
   192 integer,
intent(in) :: status
   194 if(status /= nf90_noerr) 
then   195   write(0,*)
' check netcdf status=',status
   196   write(0,
'("error ", a)')trim(nf90_strerror(status))
 integer, parameter dp
Double precision real kind. 
 
Standard integer, real, and complex single and double precision kinds. 
 
Suite of routines to perform the 2-parameter family of Extended Schmidt Gnomonic (ESG) regional grid ...
 
real(dp), parameter dtor
Degrees to radians. 
 
real(dp), parameter rtod
radians to degrees 
 
Some of the commonly used constants (pi etc) mainly for double-precision subroutines. 
 
program regional_grid
Driver routine to compute geo-referencing parameters for the Extended Schmidt Gnomonic (ESG) regional...
 
subroutine check(status)
Check results of netCDF call.