orog_mask_tools  1.8.0
mtnlm7_oclsm.f File Reference

Terrain maker for global spectral model. More...

Go to the source code of this file.

Functions/Subroutines

program __mtnlm7_oclsm.f__
 This program creates 7 terrain-related files computed from the navy 10-minute terrain dataset. More...
 
subroutine get_index (IMN, JMN, npts, lonO, latO, DELXN, jst, jen, ilist, numx)
 Determine the location of a cubed-sphere point within the high-resolution orography data. More...
 
real function get_lat_angle (dy, DEGRAD)
 Convert the 'y' direction distance of a cubed-sphere grid point to the corresponding distance in latitude. More...
 
real function get_lon_angle (dx, lat, DEGRAD)
 Convert the 'x' direction distance of a cubed-sphere grid point to the corresponding distance in longitude. More...
 
subroutine get_mismatch_index (im_in, jm_in, geolon_in, geolat_in, bitmap_in, num_out, lon_out, lat_out, iindx, jindx)
 For unmapped land points, find the nearest land point on the input data and pass back its i/j index. More...
 
real function get_xnsum (lon1, lat1, lon2, lat2, IMN, JMN, glat, zavg, zslm, delxn)
 Count the number of high-resolution orography points that are higher than the model grid box average orography height. More...
 
subroutine get_xnsum2 (lon1, lat1, lon2, lat2, IMN, JMN, glat, zavg, zslm, delxn, xnsum1, xnsum2, HC)
 Count the number of high-resolution orography points that are higher than a critical value inside a model grid box (or a portion of a model grid box). More...
 
subroutine get_xnsum3 (lon1, lat1, lon2, lat2, IMN, JMN, glat, zavg, zslm, delxn, xnsum1, xnsum2, HC)
 Count the number of high-resolution orography points that are higher than a critical value inside a model grid box (or a portion of a model grid box). More...
 
subroutine gg2rg (im, jm, numi, a)
 Convert from a full grid to a reduced grid. More...
 
logical function inside_a_polygon (lon1, lat1, npts, lon2, lat2)
 Check if a point is inside a polygon. More...
 
subroutine interpolate_mismatch (im_in, jm_in, data_in, num_out, data_out, iindx, jindx)
 Replace unmapped model land points with the nearest land point on the input grid. More...
 
subroutine latlon2xyz (siz, lon, lat, x, y, z)
 Convert from latitude and longitude to x,y,z coordinates. More...
 
subroutine MAKEMT (ZAVG, ZSLM, ORO, SLM, VAR, VAR4, GLAT, IST, IEN, JST, JEN, IM, JM, IMN, JMN, XLAT, numi)
 Create the orography, land-mask, standard deviation of orography and the convexity on a model gaussian grid. More...
 
subroutine MAKEMT2 (ZAVG, ZSLM, ORO, SLM, land_frac, VAR, VAR4, GLAT, IM, JM, IMN, JMN, lon_c, lat_c)
 Create the orography, land-mask, land fraction, standard deviation of orography and the convexity on a model cubed-sphere tile. More...
 
subroutine MAKEOA (ZAVG, VAR, GLAT, OA4, OL, IOA4, ELVMAX, ORO, oro1, XNSUM, XNSUM1, XNSUM2, XNSUM3, XNSUM4, IST, IEN, JST, JEN, IM, JM, IMN, JMN, XLAT, numi)
 Create orographic asymmetry and orographic length scale on the model grid. More...
 
subroutine MAKEOA2 (ZAVG, zslm, VAR, GLAT, OA4, OL, IOA4, ELVMAX, ORO, oro1, XNSUM, XNSUM1, XNSUM2, XNSUM3, XNSUM4, IM, JM, IMN, JMN, lon_c, lat_c, lon_t, lat_t, dx, dy, is_south_pole, is_north_pole)
 Create orographic asymmetry and orographic length scale on the model grid. More...
 
subroutine MAKEOA3 (ZAVG, zslm, VAR, GLAT, OA4, OL, IOA4, ELVMAX, ORO, SLM, oro1, XNSUM, XNSUM1, XNSUM2, XNSUM3, XNSUM4, IM, JM, IMN, JMN, lon_c, lat_c, lon_t, lat_t, is_south_pole, is_north_pole, IMI, JMI, OA_IN, OL_IN, slm_in, lon_in, lat_in)
 Create orographic asymmetry and orographic length scale on the model grid. More...
 
subroutine MAKEPC (ZAVG, ZSLM, THETA, GAMMA, SIGMA, GLAT, IST, IEN, JST, JEN, IM, JM, IMN, JMN, XLAT, numi)
 Make the principle coordinates - slope of orography, anisotropy, angle of mountain range with respect to east. More...
 
subroutine MAKEPC2 (ZAVG, ZSLM, THETA, GAMMA, SIGMA, GLAT, IM, JM, IMN, JMN, lon_c, lat_c)
 Make the principle coordinates - slope of orography, anisotropy, angle of mountain range with respect to east. More...
 
subroutine maxmin (ia, len, tile)
 Print the maximum, mininum, mean and standard deviation of an array. More...
 
subroutine minmaxj (IM, JM, A, title)
 Print out the maximum and minimum values of an array and their i/j location. More...
 
subroutine minmxj (IM, JM, A, title)
 Print out the maximum and minimum values of an array. More...
 
subroutine mnmxja (IM, JM, A, imax, jmax, title)
 Print out the maximum and minimum values of an array. More...
 
subroutine nanc (a, l, c)
 Report NaNS and NaNQ within an address range for 8-byte real words. More...
 
subroutine read_g (glob, ITOPO)
 Read input global 30-arc second orography data. More...
 
subroutine REVERS (IM, JM, numi, F, WRK)
 Reverse the east-west and north-south axes in a two-dimensional array. More...
 
subroutine rg2gg (im, jm, numi, a)
 Convert from a reduced grid to a full grid. More...
 
subroutine SPFFT1 (IMAX, INCW, INCG, KMAX, W, G, IDIR)
 Perform multiple fast fourier transforms. More...
 
real function spherical_angle (v1, v2, v3)
 Compute spherical angle. More...
 
real function spherical_distance (theta1, phi1, theta2, phi2)
 Compute a great circle distance between two points. More...
 
subroutine TERSUB (IMN, JMN, IM, JM, NM, NR, NF0, NF1, NW, EFAC, BLAT, OUTGRID, INPUTOROG)
 Driver routine to compute terrain. More...
 
real function timef ()
 Get the date/time for the system clock. More...
 

Detailed Description

Terrain maker for global spectral model.

Author
Mark Iredell
Date
92-04-16

Definition in file mtnlm7_oclsm.f.

Function/Subroutine Documentation

◆ __mtnlm7_oclsm.f__()

program __mtnlm7_oclsm.f__ ( )

This program creates 7 terrain-related files computed from the navy 10-minute terrain dataset.

The model physics grid parameters and spectral truncation and filter parameters are read by this program as input.

The 7 files produced are:

  1. sea-land mask on model physics grid
  2. gridded orography on model physics grid
  3. mountain std dev on model physics grid
  4. spectral orography in spectral domain
  5. unfiltered gridded orography on model physics grid
  6. grib sea-land mask on model physics grid
  7. grib gridded orography on model physics grid

The orography is only filtered for wavenumbers greater than nf0. For wavenumbers n between nf0 and nf1, the orography is filtered by the factor 1-((n-nf0)/(nf1-nf0))**2. The filtered orography will not have information beyond wavenumber nf1.

PROGRAM HISTORY LOG:

  • 92-04-16 IREDELL
  • 98-02-02 IREDELL FILTER
  • 98-05-31 HONG Modified for subgrid orography used in Kim's scheme
  • 98-12-31 HONG Modified for high-resolution GTOPO orography
  • 99-05-31 HONG Modified for getting OL4 (mountain fraction)
    • 00-02-10 Moorthi's modifications
  • 00-04-11 HONG Modified for reduced grids
  • 00-04-12 Iredell Modified for reduced grids
  • 02-01-07 (j) modified for principal axes of orography There are now 14 files, 4 additional for lm mb
    • 04-04-04 (j) re-Test on IST/ilen calc for sea-land mask(j)
    • 04-09-04 minus sign here in MAKEOA IST and IEN as in MAKEMT!
    • 05-09-05 if test on HK and HLPRIM for GAMMA SQRT
    • 07-08-07 replace 8' with 30" incl GICE, conintue w/ S-Y. lake slm - 08-08-07 All input 30", UMD option, and filter as described below Quadratic filter applied by default. NF0 is normally set to an even value beyond the previous truncation, for example, for jcap=382, NF0=254+2 NF1 is set as jcap+2 (and/or nearest even), eg., for t382, NF1=382+2=384 if no filter is desired then NF1=NF0=0 and ORF=ORO but if no filter but spectral to grid (with gibbs) then NF1=jcap+2, and NF1=jcap+1

      INPUT FILES:

    • UNIT5 - PHYSICS LONGITUDES (IM), PHYSICS LATITUDES (JM), SPECTRAL TRUNCATION (NM), RHOMBOIDAL FLAG (NR), AND FIRST AND SECOND FILTER PARAMETERS (NF0,NF1). RESPECTIVELY READ IN FREE FORMAT.
    • UNIT235 - GTOPO 30" AVR for ZAVG elevation - UNIT10 - 30" UMD land (lake) cover mask see MSKSRC switch
    • XUNIT11 - GTOPO AVR
    • XUNIT12 - GTOPO STD DEV
    • XUNIT13 - GTOPO MAX
    • UNIT14 - GTOPO SLM (10' NAVY if switched to get lakes
    • UNIT15 - GICE Grumbine 30" RAMP Antarctica orog IMNx3616
    • UNIT25 - Ocean land-sea mask on gaussian grid
      OUTPUT FILES:
    • UNIT51 - SEA-LAND MASK (IM,JM)
    • UNIT52 - GRIDDED OROGRAPHY (IM,JM)
    • UNIT54 - SPECTRAL OROGRAPHY ((NM+1)*((NR+1)*NM+2))
    • UNIT55 - UNFILTERED GRIDDED OROGRAPHY (IM,JM)
    • UNIT57 - GRIB GRIDDED OROGRAPHY (IM,JM)

      SUBPROGRAMS CALLED:

    • UNIQUE:
    • TERSUB - MAIN SUBPROGRAM
    • SPLAT - COMPUTE GAUSSIAN LATITUDES OR EQUALLY-SPACED LATITUDES
    • LIBRARY:
    • SPTEZ - SPHERICAL TRANSFORM
    • GBYTES - UNPACK BITS
Returns
0 for success, error code otherwise.

Definition at line 78 of file mtnlm7_oclsm.f.

References netcdf_err().

◆ get_index()

subroutine get_index ( integer, intent(in)  IMN,
integer, intent(in)  JMN,
integer  npts,
real, dimension(npts), intent(in)  lonO,
real, dimension(npts), intent(in)  latO,
real, intent(in)  DELXN,
integer, intent(out)  jst,
integer, intent(out)  jen,
integer, dimension(imn), intent(out)  ilist,
integer, intent(out)  numx 
)

Determine the location of a cubed-sphere point within the high-resolution orography data.

The location is described by the range of i/j indices on the high-res grid.

Parameters
[in]imn'i' dimension of the high-resolution orography data set.
[in]jmn'j' dimension of the high-resolution orography data set.
[in]nptsNumber of vertices to describe the cubed-sphere point.
[in]lonOThe longitudes of the cubed-sphere vertices.
[in]latOThe latitudes of the cubed-sphere vertices.
[in]delxnResolution of the high-resolution orography data set.
[out]jstStarting 'j' index on the high-resolution grid.
[out]jenEnding 'j' index on the high-resolution grid.
[out]ilistList of 'i' indices on the high-resolution grid.
[out]numxThe number of 'i' indices on the high-resolution grid.
Author
GFDL programmer

Definition at line 1790 of file mtnlm7_oclsm.f.

Referenced by MAKEMT2(), MAKEOA2(), MAKEOA3(), and MAKEPC2().

◆ get_lat_angle()

real function get_lat_angle ( real  dy,
real  DEGRAD 
)

Convert the 'y' direction distance of a cubed-sphere grid point to the corresponding distance in latitude.

Parameters
[in]dyDistance along the 'y' direction of a cubed-sphere point.
[in]degradConversion from radians to degrees.
Returns
get_lat_angle Corresponding distance in latitude.
Author
GFDL programmer

Definition at line 2938 of file mtnlm7_oclsm.f.

◆ get_lon_angle()

real function get_lon_angle ( real  dx,
real  lat,
real  DEGRAD 
)

Convert the 'x' direction distance of a cubed-sphere grid point to the corresponding distance in longitude.

Parameters
[in]dxDistance along the 'x' direction of a cubed-sphere grid point.
[in]latLatitude of the cubed-sphere point.
[in]degradConversion from radians to degrees.
Returns
get_lon_angle Corresponding distance in longitude.
Author
GFDL programmer

Definition at line 2919 of file mtnlm7_oclsm.f.

◆ get_mismatch_index()

subroutine get_mismatch_index ( integer, intent(in)  im_in,
integer, intent(in)  jm_in,
real, dimension(im_in,jm_in), intent(in)  geolon_in,
real, dimension(im_in,jm_in), intent(in)  geolat_in,
logical*1, dimension(im_in,jm_in), intent(in)  bitmap_in,
integer, intent(in)  num_out,
real, dimension(num_out), intent(in)  lon_out,
real, dimension(num_out), intent(in)  lat_out,
integer, dimension(num_out), intent(out)  iindx,
integer, dimension(num_out), intent(out)  jindx 
)

For unmapped land points, find the nearest land point on the input data and pass back its i/j index.

Parameters
[in]im_in'i' dimension of input data.
[in]jm_in'j' dimension of input data.
[in]geolon_inLongitude of input data.
[in]geolat_inLatitude of input data.
[in]bitmap_inBitmap (mask) of input data.
[in]num_outNumber of unmapped points.
[in]lon_outLongitude of unmapped points.
[in]lat_outLatitude of unmapped points.
[out]iindx'i' indices of nearest land points on the input data.
[out]jindx'j' indices of nearest land points on the input data.
Author
GFDL progammer

Definition at line 3415 of file mtnlm7_oclsm.f.

References spherical_distance().

Referenced by MAKEOA3().

◆ get_xnsum()

real function get_xnsum ( real  lon1,
real  lat1,
real  lon2,
real  lat2,
integer  IMN,
integer  JMN,
real, dimension(jmn)  glat,
integer, dimension(imn,jmn)  zavg,
integer, dimension(imn,jmn)  zslm,
real  delxn 
)

Count the number of high-resolution orography points that are higher than the model grid box average orography height.

Parameters
[in]lon1Longitude of corner point 1 of the model grid box.
[in]lat1Latitude of corner point 1 of the model grid box.
[in]lon2Longitude of corner point 2 of the model grid box.
[in]lat2Latitude of corner point 2 of the model grid box.
[in]imn'i' dimension of the high-resolution orography data.
[in]jmn'j' dimension of the high-resolution orography data.
[in]glatLatitude of each row of the high-resolution orography data.
[in]zavgThe high-resolution orography.
[in]zslmThe high-resolution land mask.
[in]delxnResolution of the high-res orography data.
Returns
get_xnsum The number of high-res points above the mean orography.
Author
GFDL Programmer

Definition at line 4512 of file mtnlm7_oclsm.f.

◆ get_xnsum2()

subroutine get_xnsum2 ( real  lon1,
real  lat1,
real  lon2,
real  lat2,
integer  IMN,
integer  JMN,
real, dimension(jmn)  glat,
integer, dimension(imn,jmn)  zavg,
integer, dimension(imn,jmn)  zslm,
real  delxn,
real, intent(out)  xnsum1,
real, intent(out)  xnsum2,
real, intent(out)  HC 
)

Count the number of high-resolution orography points that are higher than a critical value inside a model grid box (or a portion of a model grid box).

The critical value is a function of the standard deviation of orography.

Parameters
[in]lon1Longitude of corner point 1 of the model grid box.
[in]lat1Latitude of corner point 1 of the model grid box.
[in]lon2Longitude of corner point 2 of the model grid box.
[in]lat2Latitude of corner point 2 of the model grid box.
[in]imn'i' dimension of the high-resolution orography data.
[in]jmn'j' dimension of the high-resolution orography data.
[in]glatLatitude of each row of the high-resolution orography data.
[in]zavgThe high-resolution orography.
[in]zslmThe high-resolution land mask.
[in]delxnResolution of the high-res orography data.
[out]xnsum1The number of high-resolution orography above the critical value inside a model grid box.
[out]xnsum2The number of high-resolution orography points inside a model grid box.
[out]hcCritical height.
Author
GFDL Programmer

Definition at line 4619 of file mtnlm7_oclsm.f.

Referenced by MAKEOA2().

◆ get_xnsum3()

subroutine get_xnsum3 ( real  lon1,
real  lat1,
real  lon2,
real  lat2,
integer  IMN,
integer  JMN,
real, dimension(jmn)  glat,
integer, dimension(imn,jmn)  zavg,
integer, dimension(imn,jmn)  zslm,
real  delxn,
real, intent(out)  xnsum1,
real, intent(out)  xnsum2,
real  HC 
)

Count the number of high-resolution orography points that are higher than a critical value inside a model grid box (or a portion of a model grid box).

Unlike routine get_xnsum2(), this routine does not compute the critical value. Rather, it is passed in.

Parameters
[in]lon1Longitude of corner point 1 of the model grid box.
[in]lat1Latitude of corner point 1 of the model grid box.
[in]lon2Longitude of corner point 2 of the model grid box.
[in]lat2Latitude of corner point 2 of the model grid box.
[in]imn'i' dimension of the high-resolution orography data.
[in]jmn'j' dimension of the high-resolution orography data.
[in]glatLatitude of each row of the high-resolution orography data.
[in]zavgThe high-resolution orography.
[in]zslmThe high-resolution land mask.
[in]delxnResolution of the high-res orography data.
[out]xnsum1The number of high-resolution orography above the critical value inside a model grid box.
[out]xnsum2The number of high-resolution orography points inside a model grid box.
[in]hcCritical height.
Author
GFDL Programmer

Definition at line 4715 of file mtnlm7_oclsm.f.

Referenced by MAKEOA2().

◆ gg2rg()

subroutine gg2rg ( integer, intent(in)  im,
integer, intent(in)  jm,
integer, dimension(jm), intent(in)  numi,
real, dimension(im,jm), intent(inout)  a 
)

Convert from a full grid to a reduced grid.

Parameters
[in]im'i' dimension of the full grid.
[in]jm'j' dimension of the full grid.
[in]numiNumber of 'i' points for each row of the reduced grid.
[in,out]aThe data to be converted.
Author
Jordan Alpert NOAA/EMC

Definition at line 4010 of file mtnlm7_oclsm.f.

◆ inside_a_polygon()

logical function inside_a_polygon ( real  lon1,
real  lat1,
integer  npts,
real, dimension(npts)  lon2,
real, dimension(npts)  lat2 
)

Check if a point is inside a polygon.

Parameters
[in]lon1Longitude of the point to check.
[in]lat1Latitude of the point to check.
[in]nptsNumber of polygon vertices.
[in]lon2Longitude of the polygon vertices.
[in]lat2Latitude of the polygon vertices.
Returns
inside_a_polygon When true, point is within the polygon.
Author
GFDL programmer

Definition at line 4416 of file mtnlm7_oclsm.f.

References latlon2xyz(), and spherical_angle().

◆ interpolate_mismatch()

subroutine interpolate_mismatch ( integer, intent(in)  im_in,
integer, intent(in)  jm_in,
real, dimension(im_in,jm_in), intent(in)  data_in,
integer, intent(in)  num_out,
real, dimension(num_out), intent(out)  data_out,
integer, dimension(num_out), intent(in)  iindx,
integer, dimension(num_out), intent(in)  jindx 
)

Replace unmapped model land points with the nearest land point on the input grid.

Parameters
[in]im_in'i' dimension of input grid.
[in]jm_in'j' dimension of input grid.
[in]data_inInput grid data.
[in]num_outNumber of unmapped model points.
[out]data_outData on the model tile.
[in]iindx'i' indices of the nearest land points on the input grid.
[in]jindx'j' indices of the nearest land points on the input grid.
Author
GFDL programmer

Definition at line 3496 of file mtnlm7_oclsm.f.

Referenced by MAKEOA3().

◆ latlon2xyz()

subroutine latlon2xyz ( integer, intent(in)  siz,
real, dimension(siz), intent(in)  lon,
real, dimension(siz), intent(in)  lat,
real, dimension(siz), intent(out)  x,
real, dimension(siz), intent(out)  y,
real, dimension(siz), intent(out)  z 
)

Convert from latitude and longitude to x,y,z coordinates.

Parameters
[in]sizNumber of points to convert.
[in]lonLongitude of points to convert.
[in]latLatitude of points to convert.
[out]x'x' coordinate of the converted points.
[out]y'y' coordinate of the converted points.
[out]z'z' coordinate of the converted points.
Author
GFDL programmer

Definition at line 4344 of file mtnlm7_oclsm.f.

Referenced by inside_a_polygon().

◆ MAKEMT()

subroutine MAKEMT ( integer, dimension(imn,jmn)  ZAVG,
integer, dimension(imn,jmn)  ZSLM,
dimension(im,jm)  ORO,
dimension(im,jm)  SLM,
dimension(im,jm)  VAR,
dimension(im,jm)  VAR4,
dimension(jmn)  GLAT,
dimension(im,jm)  IST,
dimension(im,jm)  IEN,
dimension(jm)  JST,
dimension(jm)  JEN,
  IM,
  JM,
  IMN,
  JMN,
dimension(jm)  XLAT,
dimension(jm)  numi 
)

Create the orography, land-mask, standard deviation of orography and the convexity on a model gaussian grid.

This routine was used for the spectral GFS model.

Parameters
[in]zavgThe high-resolution input orography dataset.
[in]zslmThe high-resolution input land-mask dataset.
[out]oroOrography on the model grid.
[out]slmLand-mask on the model grid.
[out]varStandard deviation of orography on the model grid.
[out]var4Convexity on the model grid.
[out]glatLatitude of each row of the high-resolution orography and land-mask datasets.
[out]istThis is the 'i' index of high-resolution data set at the east edge of the model grid cell. the high-resolution dataset with respect to the 'east' edge
[out]ienThis is the 'i' index of high-resolution data set at the west edge of the model grid cell.
[out]jstThis is the 'j' index of high-resolution data set at the south edge of the model grid cell.
[out]jenThis is the 'j' index of high-resolution data set at the north edge of the model grid cell.
[in]im"i" dimension of the model grid.
[in]jm"j" dimension of the model grid.
[in]imn"i" dimension of the hi-res input orog/mask dataset.
[in]jmn"j" dimension of the hi-res input orog/mask dataset.
[in]xlatThe latitude of each row of the model grid.
[in]numiFor reduced gaussian grids, the number of 'i' points for each 'j' row.
Author
Jordan Alpert NOAA/EMC

Definition at line 1624 of file mtnlm7_oclsm.f.

◆ MAKEMT2()

subroutine MAKEMT2 ( integer, dimension(imn,jmn)  ZAVG,
integer, dimension(imn,jmn)  ZSLM,
real, dimension(im,jm)  ORO,
real, dimension(im,jm)  SLM,
real, dimension(im,jm)  land_frac,
real, dimension(im,jm)  VAR,
real, dimension(im,jm)  VAR4,
real, dimension(jmn)  GLAT,
integer  IM,
integer  JM,
integer  IMN,
integer  JMN,
real, dimension(im+1,jm+1)  lon_c,
real, dimension(im+1,jm+1)  lat_c 
)

Create the orography, land-mask, land fraction, standard deviation of orography and the convexity on a model cubed-sphere tile.

This routine is used for the FV3GFS model.

Parameters
[in]zavgThe high-resolution input orography dataset.
[in]zslmThe high-resolution input land-mask dataset.
[out]oroOrography on the model tile.
[out]slmLand-mask on the model tile.
[out]land_fracLand fraction on the model tile.
[out]varStandard deviation of orography on the model tile.
[out]var4Convexity on the model tile.
[out]glatLatitude of each row of the high-resolution orography and land-mask datasets.
[in]im"i" dimension of the model grid.
[in]jm"j" dimension of the model grid.
[in]imn"i" dimension of the hi-res input orog/mask datasets.
[in]jmn"j" dimension of the hi-res input orog/mask datasets.
[in]lon_cLongitude of the model grid corner points.
[in]lat_cLatitude on the model grid corner points.
Author
GFDL Programmer

Definition at line 1884 of file mtnlm7_oclsm.f.

References get_index().

◆ MAKEOA()

subroutine MAKEOA ( integer, dimension(imn,jmn)  ZAVG,
dimension(im,jm)  VAR,
dimension(jmn)  GLAT,
dimension(im,jm,4)  OA4,
dimension(im,jm,4)  OL,
dimension(im,jm,4)  IOA4,
dimension(im,jm)  ELVMAX,
dimension(im,jm)  ORO,
dimension(im,jm)  oro1,
dimension(im,jm)  XNSUM,
dimension(im,jm)  XNSUM1,
dimension(im,jm)  XNSUM2,
dimension(im,jm)  XNSUM3,
dimension(im,jm)  XNSUM4,
dimension(im,jm)  IST,
dimension(im,jm)  IEN,
dimension(jm)  JST,
dimension(jm)  JEN,
  IM,
  JM,
  IMN,
  JMN,
dimension(jm)  XLAT,
dimension(jm)  numi 
)

Create orographic asymmetry and orographic length scale on the model grid.

This routine is used for the spectral GFS gaussian grid.

Parameters
[in]zavgThe high-resolution input orography dataset.
[in]varStandard deviation of orography on the model grid.
[out]glatLatitude of each row of input terrain dataset.
[out]oa4Orographic asymmetry on the model grid. Four directional components - W/S/SW/NW
[out]olOrographic length scale on the model grid. Four directional components - W/S/SW/NW
[out]ioa4Count of oa4 values between certain thresholds.
[out]elvmaxMaximum elevation on the model grid.
[in]oroOrography on the model grid.
[out]oro1Save array for model grid orography.
[out]xnsumNumber of high-resolution orography points higher than the model grid box average.
[out]xnsum1Number of high-resolution orography points higher than the critical height.
[out]xnsum2Total number of high-resolution orography points within a model grid box.
[out]xnsum3Same as xnsum1, except shifted by half a model grid box.
[out]xnsum4Same as xnsum2, except shifted by half a model grid box.
[out]istThis is the 'i' index of high-resolution data set at the east edge of the model grid cell.
[out]ienThis is the 'i' index of high-resolution data set at the west edge of the model grid cell.
[out]jstThis is the 'j' index of high-resolution data set at the south edge of the model grid cell.
[out]jenThis is the 'j' index of high-resolution data set at the north edge of the model grid cell.
[in]im"i" dimension of the model grid.
[in]jm"j" dimension of the model grid.
[in]imn"i" dimension of the input terrain dataset.
[in]jmn"j" dimension of the input terrain dataset.
[in]xlatThe latitude of each row of the model grid.
[in]numiFor reduced gaussian grids, the number of 'i' points for each 'j' row.
Author
Jordan Alpert NOAA/EMC

Definition at line 2593 of file mtnlm7_oclsm.f.

◆ MAKEOA2()

subroutine MAKEOA2 ( integer, dimension(imn,jmn)  ZAVG,
integer, dimension(imn,jmn)  zslm,
real, dimension(im,jm)  VAR,
real, dimension(jmn)  GLAT,
real, dimension(im,jm,4)  OA4,
real, dimension(im,jm,4)  OL,
integer, dimension(im,jm,4)  IOA4,
real, dimension(im,jm)  ELVMAX,
real, dimension(im,jm)  ORO,
real, dimension(im,jm)  oro1,
real, dimension(im,jm)  XNSUM,
real, dimension(im,jm)  XNSUM1,
real, dimension(im,jm)  XNSUM2,
real, dimension(im,jm)  XNSUM3,
real, dimension(im,jm)  XNSUM4,
integer  IM,
integer  JM,
integer  IMN,
integer  JMN,
real, dimension(im+1,jm+1)  lon_c,
real, dimension(im+1,jm+1)  lat_c,
real, dimension(im,jm)  lon_t,
real, dimension(im,jm)  lat_t,
real, dimension(im,jm)  dx,
real, dimension(im,jm)  dy,
logical, dimension(im,jm)  is_south_pole,
logical, dimension(im,jm)  is_north_pole 
)

Create orographic asymmetry and orographic length scale on the model grid.

This routine is used for the cubed-sphere grid.

Parameters
[in]zavgHigh-resolution orography data.
[in]zslmHigh-resolution land-mask data.
[in]varStandard deviation of orography on the model grid.
[out]glatLatitude of each row of input terrain dataset.
[out]oa4Orographic asymmetry on the model grid. Four directional components - W/S/SW/NW
[out]olOrographic length scale on the model grid. Four directional components - W/S/SW/NW
[out]ioa4Count of oa4 values between certain thresholds.
[out]elvmaxMaximum elevation within a model grid box.
[in]oroOrography on the model grid.
[out]oro1Save array for model grid orography.
[out]xnsumNot used.
[out]xnsum1Not used.
[out]xnsum2Not used.
[out]xnsum3Not used.
[out]xnsum4Not used.
[in]im"i" dimension of the model grid tile.
[in]jm"j" dimension of the model grid tile.
[in]imn"i" dimension of the high-resolution orography and mask data.
[in]jmn"j" dimension of the high-resolution orography and mask data.
[in]lon_cCorner point longitudes of the model grid points.
[in]lat_cCorner point latitudes of the model grid points.
[in]lon_tCenter point longitudes of the model grid points.
[in]lat_tCenter point latitudes of the model grid points.
[in]dxLength of model grid points in the 'x' direction.
[in]dyLength of model grid points in the 'y' direction.
[in]is_south_poleIs the model point at the south pole?
[in]is_north_poleis the model point at the north pole?
Author
GFDL Programmer

Definition at line 2988 of file mtnlm7_oclsm.f.

References get_index(), get_xnsum2(), and get_xnsum3().

◆ MAKEOA3()

subroutine MAKEOA3 ( integer, dimension(imn,jmn)  ZAVG,
integer, dimension(imn,jmn)  zslm,
real, dimension(im,jm)  VAR,
real, dimension(jmn)  GLAT,
real, dimension(im,jm,4)  OA4,
real, dimension(im,jm,4)  OL,
integer, dimension(im,jm,4)  IOA4,
real, dimension(im,jm)  ELVMAX,
real, dimension(im,jm)  ORO,
real, dimension(im,jm)  SLM,
real, dimension(im,jm)  oro1,
real, dimension(im,jm)  XNSUM,
real, dimension(im,jm)  XNSUM1,
real, dimension(im,jm)  XNSUM2,
real, dimension(im,jm)  XNSUM3,
real, dimension(im,jm)  XNSUM4,
integer  IM,
integer  JM,
integer  IMN,
integer  JMN,
real, dimension(im+1,jm+1)  lon_c,
real, dimension(im+1,jm+1)  lat_c,
real, dimension(im,jm)  lon_t,
real, dimension(im,jm)  lat_t,
logical, dimension(im,jm)  is_south_pole,
logical, dimension(im,jm)  is_north_pole,
integer  IMI,
integer  JMI,
real, dimension(imi,jmi,4)  OA_IN,
real, dimension(imi,jmi,4)  OL_IN,
real, dimension(imi,jmi)  slm_in,
real, dimension(imi,jmi)  lon_in,
real, dimension(imi,jmi)  lat_in 
)

Create orographic asymmetry and orographic length scale on the model grid.

This routine is used for the cubed-sphere grid. The asymmetry and length scales are interpolated from an existing gfs orography file. The maximum elevation is computed from the high-resolution orography data.

Parameters
[in]zavgHigh-resolution orography data.
[in]zslmHigh-resolution land-mask data. Not used.
[in]varStandard deviation of orography on the model grid.
[out]glatLatitude of each row of input terrain dataset.
[out]oa4Orographic asymmetry on the model grid. Four directional components - W/S/SW/NW
[out]olOrographic length scale on the model grid. Four directional components - W/S/SW/NW
[out]ioa4Count of oa4 values between certain thresholds.
[out]elvmaxMaximum elevation within a model grid box.
[in]slmLand-mask on model grid.
[in]oroOrography on the model grid.
[out]oro1Save array for model grid orography.
[in]xnsumNot used.
[in]xnsum1Not used.
[in]xnsum2Not used.
[in]xnsum3Not used.
[in]xnsum4Not used.
[in]im"i" dimension of the model grid tile.
[in]jm"j" dimension of the model grid tile.
[in]imn"i" dimension of the high-resolution orography and mask data.
[in]jmn"j" dimension of the high-resolution orography and mask data.
[in]lon_cCorner point longitudes of the model grid points.
[in]lat_cCorner point latitudes of the model grid points.
[in]lon_tCenter point longitudes of the model grid points.
[in]lat_tCenter point latitudes of the model grid points.
[in]is_south_poleNot used.
[in]is_north_poleNot used.
[in]imi'i' dimension of input gfs orography data.
[in]jmi'j' dimension of input gfs orography data.
[in]oa_inAsymmetry on the input gfs orography data.
[in]ol_inLength scales on the input gfs orography data.
[in]slm_inLand-sea mask on the input gfs orography data.
[in]lon_inLongitude on the input gfs orography data.
[in]lat_inLatitude on the input gfs orography data.
Author
Jordan Alpert NOAA/EMC

Definition at line 3556 of file mtnlm7_oclsm.f.

References get_index(), get_mismatch_index(), and interpolate_mismatch().

◆ MAKEPC()

subroutine MAKEPC ( integer, dimension(imn,jmn)  ZAVG,
integer, dimension(imn,jmn)  ZSLM,
dimension(im,jm)  THETA,
dimension(im,jm)  GAMMA,
dimension(im,jm)  SIGMA,
dimension(jmn)  GLAT,
dimension(im,jm)  IST,
dimension(im,jm)  IEN,
dimension(jm)  JST,
dimension(jm)  JEN,
  IM,
  JM,
  IMN,
  JMN,
dimension(jm)  XLAT,
dimension(jm)  numi 
)

Make the principle coordinates - slope of orography, anisotropy, angle of mountain range with respect to east.

This routine is used for spectral GFS gaussian grids.

Parameters
[in]zavgThe high-resolution input orography dataset.
[in]zslmThe high-resolution input land-mask dataset.
[out]thetaAngle of mountain range with respect to east for each model point.
[out]gammaAnisotropy for each model point.
[out]sigmaSlope of orography for each model point.
[out]glatLatitude of each row of the high-resolution orography and land-mask datasets.
[out]istThis is the 'i' index of high-resolution data set at the east edge of the model grid cell.
[out]ienThis is the 'i' index of high-resolution data set at the west edge of the model grid cell.
[out]jstThis is the 'j' index of high-resolution data set at the south edge of the model grid cell.
[out]jenThis is the 'j' index of high-resolution data set at the north edge of the model grid cell.
[in]im"i" dimension of the model grid tile.
[in]jm"j" dimension of the model grid tile.
[in]imn"i" dimension of the hi-res input orog/mask datasets.
[in]jmn"j" dimension of the hi-res input orog/mask datasets.
[in]xlatThe latitude of each row of the model grid.
[in]numiFor reduced gaussian grids, the number of 'i' points for each 'j' row.
Author
Jordan Alpert NOAA/EMC

Definition at line 2044 of file mtnlm7_oclsm.f.

◆ MAKEPC2()

subroutine MAKEPC2 ( integer, dimension(imn,jmn)  ZAVG,
integer, dimension(imn,jmn)  ZSLM,
real, dimension(im,jm)  THETA,
real, dimension(im,jm)  GAMMA,
real, dimension(im,jm)  SIGMA,
real, dimension(jmn)  GLAT,
integer  IM,
integer  JM,
integer  IMN,
integer  JMN,
real, dimension(im+1,jm+1)  lon_c,
real, dimension(im+1,jm+1)  lat_c 
)

Make the principle coordinates - slope of orography, anisotropy, angle of mountain range with respect to east.

This routine is used for the FV3GFS cubed-sphere grid.

Parameters
[in]zavgThe high-resolution input orography dataset.
[in]zslmThe high-resolution input land-mask dataset.
[out]thetaAngle of mountain range with respect to east for each model point.
[out]gammaAnisotropy for each model point.
[out]sigmaSlope of orography for each model point.
[out]glatLatitude of each row of the high-resolution orography and land-mask datasets.
[in]im"i" dimension of the model grid tile.
[in]jm"j" dimension of the model grid tile.
[in]imn"i" dimension of the hi-res input orog/mask datasets.
[in]jmn"j" dimension of the hi-res input orog/mask datasets.
[in]lon_cLongitude of model grid corner points.
[in]lat_cLatitude of the model grid corner points.
Author
GFDL Programmer

Definition at line 2318 of file mtnlm7_oclsm.f.

References get_index().

◆ maxmin()

subroutine maxmin ( integer*2, dimension(len)  ia,
integer  len,
character*7  tile 
)

Print the maximum, mininum, mean and standard deviation of an array.

Parameters
[in]iaThe array to be checked.
[in]lenThe number of points to be checked.
[in]tileA name associated with the array.
Author
Jordan Alpert NOAA/EMC

Definition at line 4234 of file mtnlm7_oclsm.f.

Referenced by read_g().

◆ minmaxj()

subroutine minmaxj ( integer  IM,
integer  JM,
real(kind=4), dimension(im,jm)  A,
character*8  title 
)

Print out the maximum and minimum values of an array and their i/j location.

Also print out the number of undefined points.

Parameters
[in]imThe 'i' dimension of the array.
[in]jmThe 'i' dimension of the array.
[in]aThe array to check.
[in]titleName of the data to be checked.
Author
Jordan Alpert NOAA/EMC

Definition at line 4292 of file mtnlm7_oclsm.f.

◆ minmxj()

subroutine minmxj ( integer  IM,
integer  JM,
real, dimension(im,jm)  A,
character*8  title 
)

Print out the maximum and minimum values of an array.

Parameters
[in]imThe 'i' dimension of the array.
[in]jmThe 'i' dimension of the array.
[in]aThe array to check.
[in]titleName of the data to be checked.
Author
Jordan Alpert NOAA/EMC

Definition at line 4036 of file mtnlm7_oclsm.f.

Referenced by TERSUB().

◆ mnmxja()

subroutine mnmxja ( integer  IM,
integer  JM,
real, dimension(im,jm)  A,
integer  imax,
integer  jmax,
character*8  title 
)

Print out the maximum and minimum values of an array.

Pass back the i/j location of the maximum value.

Parameters
[in]imThe 'i' dimension of the array.
[in]jmThe 'i' dimension of the array.
[in]aThe array to check.
[out]imax'i' location of maximum
[out]jmax'j' location of maximum
[in]titleName of the data to be checked.
Author
Jordan Alpert NOAA/EMC

Definition at line 4071 of file mtnlm7_oclsm.f.

Referenced by TERSUB().

◆ nanc()

subroutine nanc (   a,
  l,
character*(*)  c 
)

Report NaNS and NaNQ within an address range for 8-byte real words.

This routine prints a single line for each call and prints a message and returns to the caller on detection of the FIRST NaN in the range. If no NaN values are found it returns silently.

Parameters
[in]AReal*8 variable or array
[in]LNumber of words to scan (length of array)
[in]CDistinctive message set in caller to indicate where the routine was called.
Author
Jordan Alpert NOAA/EMC

Definition at line 4778 of file mtnlm7_oclsm.f.

◆ read_g()

subroutine read_g ( integer*2, dimension(360*120,180*120)  glob,
integer  ITOPO 
)

Read input global 30-arc second orography data.

Parameters
[out]globThe orography data.
[in]itopoNot used.
Author
Jordan Alpert NOAA/EMC

Definition at line 4183 of file mtnlm7_oclsm.f.

References maxmin().

Referenced by TERSUB().

◆ REVERS()

subroutine REVERS (   IM,
  JM,
integer, dimension(jm)  numi,
real, dimension(im,jm)  F,
real, dimension(im,jm)  WRK 
)

Reverse the east-west and north-south axes in a two-dimensional array.

Parameters
[in]im'i' dimension of the 2-d array.
[in]jm'j' dimension of the 2-d array.
[in]numiNot used.
[in,out]fThe two-dimensional array to be processed.
[out]wrkTwo-dimensional work array.
Author
Jordan Alpert NOAA/EMC

Definition at line 3947 of file mtnlm7_oclsm.f.

◆ rg2gg()

subroutine rg2gg ( integer, intent(in)  im,
integer, intent(in)  jm,
integer, dimension(jm), intent(in)  numi,
real, dimension(im,jm), intent(inout)  a 
)

Convert from a reduced grid to a full grid.

Parameters
[in]im'i' dimension of the full grid.
[in]jm'j' dimension of the full grid.
[in]numiNumber of 'i' points for each row of the reduced grid.
[in,out]aThe data to be converted.
Author
Jordan Alpert NOAA/EMC

Definition at line 3984 of file mtnlm7_oclsm.f.

Referenced by TERSUB().

◆ SPFFT1()

subroutine SPFFT1 ( integer, intent(in)  IMAX,
integer, intent(in)  INCW,
integer, intent(in)  INCG,
integer, intent(in)  KMAX,
complex, dimension(incw,kmax), intent(inout)  W,
real, dimension(incg,kmax), intent(inout)  G,
integer, intent(in)  IDIR 
)

Perform multiple fast fourier transforms.

This subprogram performs multiple fast fourier transforms between complex amplitudes in fourier space and real values in cyclic physical space.

Subprograms called (NCEPLIB SP Library):

  • scrft Complex to real fourier transform
  • dcrft Complex to real fourier transform
  • srcft Real to complex fourier transform
  • drcft Real to complex fourier transform

Program history log: 1998-12-18 Mark Iredell

Parameters
[in]imaxInteger number of values in the cyclic physical space. See limitations on imax in remarks below.
[in]incwInteger first dimension of the complex amplitude array. (incw >= imax/2+1).
[in]incgInteger first dimension of the real value array. (incg >= imax).
[in]kmaxInteger number of transforms to perform.
[in]wComplex amplitudes on input if idir>0, and on output if idir<0.
[in]gReal values on input if idir<0, and on output if idir>0.
[in]idirInteger direction flag. idir>0 to transform from fourier to physical space. idir<0 to transform from physical to fourier space.
Note
The restrictions on imax are that it must be a multiple of 1 to 25 factors of two, up to 2 factors of three, and up to 1 factor of five, seven and eleven.
Author
Mark Iredell ORG: W/NMC23
Date
96-02-20

Definition at line 4133 of file mtnlm7_oclsm.f.

◆ spherical_angle()

real function spherical_angle ( real, dimension(3)  v1,
real, dimension(3)  v2,
real, dimension(3)  v3 
)

Compute spherical angle.

Parameters
[in]v1Vector 1.
[in]v2Vector 2.
[in]v3Vector 3.
Returns
spherical_angle Spherical Angle.
Author
GFDL programmer

Definition at line 4366 of file mtnlm7_oclsm.f.

Referenced by inside_a_polygon().

◆ spherical_distance()

real function spherical_distance ( real, intent(in)  theta1,
real, intent(in)  phi1,
real, intent(in)  theta2,
real, intent(in)  phi2 
)

Compute a great circle distance between two points.

Parameters
[in]theta1Longitude of point 1.
[in]phi1Latitude of point 1.
[in]theta2Longitude of point 2.
[in]phi2Latitude of point2.
Returns
spherical_distance Great circle distance.
Author
GFDL programmer

Definition at line 3379 of file mtnlm7_oclsm.f.

Referenced by get_mismatch_index().

◆ TERSUB()

subroutine TERSUB ( integer  IMN,
integer  JMN,
integer  IM,
integer  JM,
integer  NM,
integer  NR,
integer  NF0,
integer  NF1,
integer  NW,
integer  EFAC,
integer  BLAT,
character(len=*), intent(in)  OUTGRID,
character(len=*), intent(in)  INPUTOROG 
)

Driver routine to compute terrain.

Parameters
[in]IMN"i" dimension of the input terrain dataset.
[in]JMN"j" dimension of the input terrain dataset.
[in]IM"i" dimension of the model grid tile.
[in]JM"j" dimension of the model grid tile.
[in]NMSpectral truncation.
[in]NRRhomboidal flag.
[in]NF0First order spectral filter parameters.
[in]NF1Second order spectral filter parameters.
[in]NWNumber of waves.
[in]EFACFactor to adjust orography by its variance.
[in]BLATWhen less than zero, reverse latitude/ longitude for output.
[in]OUTGRIDThe 'grid' file for the model tile.
[in]INPUTOROGInput orography/GWD file on gaussian grid. When specified, will be interpolated to model tile. When not specified, program will create fields from raw high-resolution topography data.
Author
Jordan Alpert NOAA/EMC

Definition at line 184 of file mtnlm7_oclsm.f.

References minmxj(), mnmxja(), netcdf_err(), read_g(), rg2gg(), timef(), and write_netcdf().

◆ timef()

real function timef ( )

Get the date/time for the system clock.

Author
Mark Iredell
Returns
timef

Definition at line 4827 of file mtnlm7_oclsm.f.

Referenced by TERSUB().