MODULE MODMPASFIO New May-Nov 2017 for I/O API-3.2!!Unstructured MPAS-grid descriptions; grid related utility routines; and I/O for MPAS-format netCDF.
- Notes and Limitations
- Environment Variables
- Declarations
PARAMETERs and publicMODULE-variables describing the current MPAS grid.
INITMPGRID()andINITREARTH()- module initialization;
SPHEREDIST()earth-radius reset.SHUTMPGRID()- module shutdown, flushing and closing all files.
OPENMPAS()- open existing MPAS-format input or output file
CREATEMPAS()- create new MPAS-format output file
DESCMPAS()- read header and return description for MPAS-format input file
READMPSTEPS()- Read MPAS-format ASCII timestep-value list from an MPAS-format file and convert to I/O API conventions.
WRITEMPSTEP()- Write timestep-value record to an MPAS-format file as MPAS-format ASCII
READMPAS()- read variable or time step of variable from MPAS-format input file
WRITEMPAS()- write variable or time step of variable to MPAS-format output file
FINDCELL()andFINDVRTX()- find subscript for MPAS-cell or MPAS-dual-cell containing the point
(ALAT,ALON)ARC2MPAS()- compute 2-D and 3-D weights for arc-segment (e.g., for SMOKE)
MPCELLMATX()- Create a sparse incidence matrix (or index-array) for MPAS-cell subscripts for set of points or a rectangular grid, for cell-based interpolation.
MPBARYMATX()- Create a 3-band sparse matrix for barycentric-linear interpolation to a set of points or a 2-D rectangular grid of points.
MPBARYMULT()- Multiply a 3-band barycentric-linear sparse matrix by an MPAS variable, to perform interpolation.
MPINTERP()- Interpolate an MPAS-grid variable to a point, set of points, or 2-D grid of points, using barycentric-linear interpolation for
REALvariables, or cell-based interpolation forINTEGERvariables.SPHEREDIST()- compute spherical-Earth distance from
(LAT1,LON1)to(LAT2,LON2)or from MPAS Earth-centered Cartesian(X1,Y1,Z1)to(X2,Y2,Z2)MPSTR2DT()andMPDT2STR()- convert MPAS-format ASCII date&time to/from I/O API convention date&time
JDATE:JTIME
Neither MPAS nor netCDF has internal data structures
adequate to detect such gaps. It is your
responsibility to detect such gaps in input data, and to
take appropriate measures if such gaps exist.
(NetCDF will "fill" the gaps with special values
like NF_FILL_INT=-2147483647,
NF_FILL_REAL=9.9692099683868690e+36, etc.)
setenv MPAS_CHKFILL Y to have such checks
automatically performed by routine
READMPAS().
MODULE supports only Earth-sphere based
modeling and its files.
MODULE saves a
copy of these mandatory variables from the
INITMPGRID() initialization-file as
PUBLIC, PROTECTED module-variables: see
the next section, below.
REARTH, which is mandated
to be DOUBLE PRECISION, the MPAS specification
does not specify the types for the mandatory variables.
MODULE MPASFIO assumes that these
variables are either default REAL or
DOUBLE PRECISION for coordinate related
variables (converting the former to the latter internally,
to match REARTH), or default
INTEGER for index and ID variables.
-π/2 ≤ lat ≤ π/2
0 ≤ lon ≤ 2 π
Everyone else uses either degrees, minutes, and
seconds or real-number degrees, and (except for the WMO) in
a fashion compliant with ISO Standard 6709. This
MODULE uses real-number degrees in the range
-90 ≤ lat ≤ 90
0 ≤ lon ≤ 360
internally, and converts back and forth between that and
the on-file MPAS radians usage. Arguments to the routines
are converted to this convention, as necessary.
(X,Y,Z) coordinate system is relative
to 3-D Cartesian axes with origin X=Y=Z=0 at
the center of the Earth, with units measured in meters.
FINDCELL(), FINDVRTX(),
MPBARYMATX(), MPBARYMULT(),
and MPINTERP() when the domain is not
whole-Earth (therefore having boundary-cells) and the
requested operation is for a boundary cell or intersects
a boundary cell or if the MPAS domain is not convex (e.g.,
"C-shaped")..
MODULE
MODNCFIO also provides facilities for describing
and reading "raw" netCDF. However, the paradigm
in that module is open-action-close within each routine,
whereas the paradigm here is to separate the actions from
the opening and the closing, with one action per routine
(the file remaining open until a SHUTMPGRID()
call). These two paradigms do not mix well: attempting to
open or close an already-open netCDF file can lead to
strange-looking (e.g., "double-free") failures.
- https://mpas-dev.github.io/files/documents/MPAS-MeshSpec.pdf,
- for MPAS grid and netCDF-file specifications
MODULE MODNCFIO,- for netCDF declarations and "raw netCDF" file routines
MODULE MODWRFIO- for WRF-related file routines
- m3tools program mpasdiff
- for MPAS-format-netCDF variable and layer-range-of-variable statistical comparisons
- m3tools program mpasstat
- for MPAS-format-netCDF variable and layer-range-of-variable statistics
- m3tools program mpastom3
- for interpolating variables from MPAS-format-netCDF files to GRIDDED I/O API files.
- m3toolssample program mpaswtest.
- Example program for
ARC2MPAS()MPAS-grid allocation of emissions line-sources.
Back to Contents
Back to Contents
MPAS_VERBOSEY N [N]- Generate verbose diagnostic log output?
MPAS_CHKFILLY N [Y]- Check data read in by
READMPAS()against netCDF fill-values, to find potential "holes"?
- Logical file names
- Following the usual I/O API logical file name conventions for the files used by routines
INITMPGRID(), OPENMPAS(), CREATEMPAS(), DESCMPAS(), READMPAS(), WRITEMPAS(), READMPSTEPS(), WRITEMPSTEP():
setenv <name> <file-path>
PARAMETERsPARAMETERs:
INTEGER, PUBLIC, PARAMETER :: MPSTRLEN = 64- String length for use in MPAS formatted netCDF files (e.g., for MPAS formatted ASCII version and date-and-time strings)
INTEGER, PUBLIC, PARAMETER :: NMPASDIMS = 11- Number of "standard" MPAS-file dimensions for this module (which augments the MPAS-mandated list by adding standard time and vertical dimensions).
Note that
OPENMPAS()andCREATEMPAS()have provisions for defining additional dimensions on a per-file basis. (Having a potentially open-ended set of dimensions and their names adds substantial complexity ;-( ).
INTEGER, PUBLIC, PARAMETER :: NMPASVARS = 36- Number of "standard" MPAS-file variables for this module.
INTEGER, PUBLIC, PARAMETER :: MXMPASVARS = 2048- Maximum number of MPAS-file variables for this module.
INTEGER, PUBLIC, PARAMETER :: MXMPASFILE = 64- Maximum number of MPAS-files for this module.
CHARACTER(LEN=16), PUBLIC, PARAMETER :: MPASDIMNAMES( NMPASDIMS )- Names of the "standard" MPAS-file dimensions:
Time !! 1 TWO = 2 !! 2 StrLen = 64 !! 3 nCells !! 4 nEdges !! 5 nVertices !! 6 vertexDegree !! 7 maxEdges !! 8 maxEdges2 !! 9 nVertLevels !! 10 nVertLevelsP1 !! 11
CHARACTER(LEN=MPSTRLEN), PUBLIC, PARAMETER :: MPASVARNAMES( NMPASVARS )- Names of the "standard" MPAS-file variables, as specified by the MPAS Mesh Specification document:
indexToCellID !! 1 indexToEdgeID !! 2 indexToVertexID !! 3 nEdgesOnCell !! 4 nEdgesOnEdge !! 5 cellsOnCell !! 6 edgesOnCell !! 7 verticesOnCell !! 8 cellsOnEdge !! 9 edgesOnEdge !! 10 verticesOnEdge !! 11 cellsOnVertex !! 12 edgesOnVertex !! 13 latCell !! 14 lonCell !! 15 latEdge !! 16 lonEdge !! 17 latVertex !! 18 lonVertex !! 19 xCell !! 20 yCell !! 21 zCell !! 22 xEdge !! 23 yEdge !! 24 zEdge !! 25 xVertex !! 26 yVertex !! 27 zVertex !! 28 weightsOnEdge !! 29 dvEdge !! 30 dcEdge !! 31 angleEdge !! 32 areaCell !! 33 areaTriangle !! 34 kiteAreasOnVertex !! 35 meshDensity !! 36
PUBLIC Variables:
Note that these arePROTECTEDvariables: as far as codes thatUSEthis module are concerned, they are read-only, being set within the module by routineINITMPGRID().
We give here both the names of theMODULE-variables, and the (clumsier, more-verbose) on-file names as specified by the MPAS-grid specification-document.
CHARACTER*64 ONSPHERE- (MPAS-spec name on_a_sphere):
YESfor on-sphere MPAS grids
CHARACTER*64 MESH_ID- (mesh_id): ASCII string for grid name
CHARACTER*64 MESH_SPEC- (mesh_spec): MPAS grid specification version
INTEGER MPSTEPS- maximum timestep number in the file used by
INITMPGRID()
INTEGER MPCELLS- (nCells): number of primary cells in the mesh for the file used by
INITMPGRID()
INTEGER MPEDGES- (nEdges): number of edges in the mesh
INTEGER MPVRTXS- (nVertices): number of vertices in the mesh
INTEGER MPLVLS- (nVertLevels): number of vertical levels in the input file used by INITMPGRID()
INTEGER MPVORDR- (vertexDegree): maximum number of cells (or edges) per vertex
INTEGER MPBNDYC- (maxEdges): maximum number of boundary-edges or vertices for any cell
REAL*8 REARTH- (sphere_radius): Earth-radius (Meters)
In the absence of initialization, the default for use bySPHEREDIST()is 6370,000 Meters.
INTEGER MPDATES(MPSTEPS)- (from xtime): Julian dates (coded
YYYYDDD) for the time steps in the input file used by INITMPGRID()
INTEGER MPTIMES(MPSTEPS)- (from xtime): times (coded
HHMMSS) for the time steps in the input file used by INITMPGRID()
INTEGER CELLID(MPCELLS)- (indexToCellID): global cell-ID
INTEGER EDGEID(MPEDGES)- (indexToEdgeID): global edge-ID
INTEGER VRTXID(MPVRTXS)- (indexToVertexID): global vertex-ID
INTEGER NBNDYE(MPCELLS)- (nEdgesOnCel): Number of boundary-edges on a given cell
INTEGER NEDGEE(MPEDGES)- (nEdgesOnEdge): Number of edges on a given edge. Used to reconstruct tangential velocities.
INTEGER BNDYCELL(MPBNDYC, MPCELLS)- (cellsOnCell): indices of cells that are on the boundary of a given cell
INTEGER BNDYEDGE(MPBNDYC, MPCELLS)- (edgesOnCell): indices of edges that are on the boundary of a given cell
INTEGER BNDYVRTX(MPBNDYC, MPCELLS)- (verticesOnCell): boundary-vertex indices for this cell
INTEGER ECELLS(2, MPEDGES)- (cellsOnEdge): cell-indices for cells bounding this edge
INTEGER EVRTXS(2, MPEDGES)- (verticesOnEdge): vertex-indices for vertices bounding this edge
INTEGER VCELLS(MPVORDR, MPVRTXS)- (verticesOnCell): indices of cells that radiate from a given vertex
INTEGER VEDGES(MPVORDR, MPVRTXS)- (edgesOnVertex): indices of edges that radiate from a given vertex
INTEGER EEDGES(2*MPBNDYC, MPEDGES)- (edgesOnEdge): Edge indices used to reconstruct tangential velocities.
REAL*8 EWGHTS(2*MPBNDYC, MPEDGES)- (weightsOnEdge): Weights used to reconstruct tangential velocities.
REAL*8 ALATC(MPCELLS)- (latCell): latitudes of cell-centers
REAL*8 ALONC(MPCELLS)- (lonCell): longitudes of cell-centers, normalized so that
0 ≤ lon < 360
REAL*8 ALATE(MPEDGES)- (latEdge): latitudes of edge-centers
REAL*8 ALONE(MPEDGES)- (lonEdge): longitudes of edge-centers, normalized so that
0 ≤ lon < 360
REAL*8 ALATV(MPVRTXS)- (latVertex): latitudes of vertices
REAL*8 ALONV(MPVRTXS)- (lonVertex): longitudes of vertices, normalized so that
0 ≤ lon < 360
REAL*8 XCELL(MPCELLS)- (xCell):
X-coordinates of cell-centers
REAL*8 YCELL(MPCELLS)- (yCell):
Y-coordinates of cell-centers
REAL*8 ZCELL(MPCELLS)- (zCell):
Z-coordinates of cell-centers
REAL*8 XEDGE(MPEDGES)- (xEdge):
X-coordinates of edge-centers
REAL*8 YEDGE(MPEDGES)- (yEdge):
Y-coordinates of edge-centers
REAL*8 ZEDGE(MPEDGES)- (zEdge):
Z-coordinates of edge-centers
REAL*8 XVRTX(MPVRTXS)- (xVertex):
X-coordinates of vertices
REAL*8 YVRTX(MPVRTXS)- (yVertex):
Y-coordinates of vertices
REAL*8 ZVRTX(MPVRTXS)- (zVertex):
Z-coordinates of vertices
REAL*8 DVEDGE(MPEDGES)- (dvEdge): Distance in meters between the vertices that bound a given edge.
REAL*8 DCEDGE(MPEDGES)- (dcEdge): Distance in meters between centers of the cells that bound a given edge
REAL*8 EANGLE(MPEDGES)- (angleEdge): Angle in radians an edge's normal vector makes with the local eastward direction.
REAL*8 CAREAS(MPCELLS)- (areaCell): Areas in square meters for a given cell
REAL*8 VAREAS(MPVRTXS)- (areaTriangle): Area in square meters for the triangle of the dual mesh centered at the given vertex.
REAL*8 KAREAS(:MPVORDR,MPVRTXS)- (kiteAreasOnVertex): area of the intersection of the dual-mesh triangle centered at the given vertex with the original-mesh cells that intersect it.
REAL*8 MSHDEN(MPCELLS)- (meshDensity): Value of the generating density function at each cell center (used for irregular meshes)
Back to Contents
LOGICAL FUNCTIONs INITMPGRID() and
INITREARTH()This is a generic function with two currently implemented forms, and an explicit variant:
LOGICAL FUNCTION INITREARTH( RADEARTH )
REAL(8), INTENT( IN ) :: RADEARTH !! Earth-radius (Meters)
LOGICAL FUNCTION INITMPGRID( FNAME )
CHARACTER*(*), INTENT( IN ) :: FNAME !! logical file name
LOGICAL FUNCTION INITMPGRID( NCELLS, &
NEDGES, &
NVRTXS, &
NVLVLS, &
NVORDR, &
NBNDYC, &
NBNDY2, &
CELLID1, &
EDGEID1, &
VRTXID1, &
NBNDYE1, &
NEDGEE1, &
BYCELL1, &
BYEDGE1, &
BYVRTX1, &
ECELLS1, &
EVRTXS1, &
VCELLS1, &
VEDGES1, &
EEDGES1, &
ALATC1, &
ALONC1, &
ALATE1, &
ALONE1, &
ALATV1, &
ALONV1, &
XCELL1, &
YCELL1, &
ZCELL1, &
XEDGE1, &
YEDGE1, &
ZEDGE1, &
XVRTX1, &
YVRTX1, &
ZVRTX1, &
EWGHTS1, &
DVEDGE1, &
DCEDGE1, &
EANGLE1, &
CAREAS1, &
VAREAS1, &
KAREAS1, &
MSHDEN1, &
ONSPHERE1, &
MESH_ID1, &
MESH_SPEC1, &
REARTH1 )
INTEGER, INTENT( IN ) :: NCELLS !! # of primary cells in the mesh
INTEGER, INTENT( IN ) :: NEDGES !! # of edges in the mesh
INTEGER, INTENT( IN ) :: NVRTXS !! # of vertices in the mesh
INTEGER, INTENT( IN ) :: NVLVLS !! # of vertical levels
INTEGER, INTENT( IN ) :: NVORDR !! max vertex-order: # of cells/edges per vertex
INTEGER, INTENT( IN ) :: NBNDYC !! max # of vertices/edges per cell
INTEGER, INTENT( IN ) :: NBNDY2 !! 2 * max # of vertices/edges per cell
INTEGER, INTENT( IN ) :: CELLID1( NCELLS ) !! (NCELLS): global cell-ID
INTEGER, INTENT( IN ) :: EDGEID1( NEDGES ) !! (NEDGES): global edge-ID
INTEGER, INTENT( IN ) :: VRTXID1( NVRTXS ) !! (NVRTXS): global vertex-ID
INTEGER, INTENT( IN ) :: NBNDYE1( NCELLS ) !! (NCELLS): # of edges per cell
INTEGER, INTENT( IN ) :: NEDGEE1( NEDGES ) !! (NEDGES): # of edges per edge 1(0 or 1)
INTEGER, INTENT( IN ) :: BYCELL1( NBNDYC, NCELLS ) !! (NBNDYC,NCELLS): boundary-cell indices for this cell
INTEGER, INTENT( IN ) :: BYEDGE1( NBNDYC, NCELLS ) !! (NBNDYC,NCELLS): boundary-edge indices for this cell
INTEGER, INTENT( IN ) :: BYVRTX1( NBNDYC, NCELLS ) !! (NBNDYC,NCELLS): boundary-vertex indices for this cell
INTEGER, INTENT( IN ) :: ECELLS1( 2,NEDGES ) !! (2,NEDGES): cell indices for this edge
INTEGER, INTENT( IN ) :: EVRTXS1( 2,NEDGES ) !! (2,NEDGES): vertex indices for this edge
INTEGER, INTENT( IN ) :: VCELLS1( NVORDR,NVRTXS ) !! (NVORDR,NVRTXS): Cell indices that radiate from a given vertex.
INTEGER, INTENT( IN ) :: VEDGES1( NVORDR,NVRTXS ) !! (NVORDR,NVRTXS): Edge indices that radiate from a given vertex
INTEGER, INTENT( IN ) :: EEDGES1( 2*NBNDYC,NEDGES ) !! (NBNDY2,NEDGES): Edge indices used to reconstruct tangential velocities.
REAL(8), INTENT( IN ) :: ALATC1( NCELLS ) !! (NCELLS): latitude-degrees for cell-centers
REAL(8), INTENT( IN ) :: ALONC1( NCELLS ) !! (NCELLS): longitude-degrees for cell-centers
REAL(8), INTENT( IN ) :: ALATE1( NEDGES ) !! (NEDGES): latitude-degrees for edge-centers
REAL(8), INTENT( IN ) :: ALONE1( NEDGES ) !! (NEDGES): longitude-degrees for edge-centers
REAL(8), INTENT( IN ) :: ALATV1( NVRTXS ) !! (NVRTXS): latitude-degrees for vertices
REAL(8), INTENT( IN ) :: ALONV1( NVRTXS ) !! (NVRTXS): longitude-degrees for vertices
REAL(8), INTENT( IN ) :: XCELL1( NCELLS ) !! (NCELLS): X-coordinate for cell-centers
REAL(8), INTENT( IN ) :: YCELL1( NCELLS ) !! (NCELLS): Y-coordinate for cell-centers
REAL(8), INTENT( IN ) :: ZCELL1( NCELLS ) !! (NCELLS): Z-coordinate for cell-centers
REAL(8), INTENT( IN ) :: XEDGE1( NEDGES ) !! (NCELLS): X-coordinate for edge-centers
REAL(8), INTENT( IN ) :: YEDGE1( NEDGES ) !! (NCELLS): Y-coordinate for edge-centers
REAL(8), INTENT( IN ) :: ZEDGE1( NEDGES ) !! (NCELLS): Z-coordinate for edge-centers
REAL(8), INTENT( IN ) :: XVRTX1( NVRTXS ) !! (NVRTXS): X-coordinate for vertices
REAL(8), INTENT( IN ) :: YVRTX1( NVRTXS ) !! (NVRTXS): Y-coordinate for vertices
REAL(8), INTENT( IN ) :: ZVRTX1( NVRTXS ) !! (NVRTXS): Z-coordinate for vertice
REAL(8), INTENT( IN ) :: EWGHTS1( 2*NBNDYC,NEDGES ) !! (NBNDY2,NEDGES): weights used to reconstruct tangential velocities.
REAL(8), INTENT( IN ) :: DVEDGE1( NEDGES ) !! (NEDGES): edge lengths 1(M)
REAL(8), INTENT( IN ) :: DCEDGE1( NEDGES ) !! (NEDGES): distance between the cell-centers that saddle a given edge 1(M)
REAL(8), INTENT( IN ) :: EANGLE1( NEDGES ) !! (NEDGES): angle from edge-normal vector to Easting
REAL(8), INTENT( IN ) :: CAREAS1( NCELLS ) !! (NCELLS): cell areas 1(M^2)
REAL(8), INTENT( IN ) :: VAREAS1( NVRTXS ) !! (NVRTXS): dual-mesh triangle areas 1(M^2)
REAL(8), INTENT( IN ) :: KAREAS1( NVORDR, NVRTXS ) !! (NBNDYC,NCELLS): kite-area: intersection of cell with dual-cell centered at vertex
REAL(8), INTENT( IN ) :: MSHDEN1( NCELLS ) !! (NCELLS): mesh density 1(none)
CHARACTER*64, OPTIONAL, INTENT( IN ) :: ONSPHERE1
CHARACTER*64, OPTIONAL, INTENT( IN ) :: MESH_ID1
CHARACTER*64, OPTIONAL, INTENT( IN ) :: MESH_SPEC1
REAL(8), OPTIONAL, INTENT( IN ) :: REARTH1 !! Earth-radius (M)
TheBack to ContentsINITREARTH()form is used to re-initialize the Earth-radius for use inSPHEREDIST(), for the case that function is used in stand-alone form.Preconditions for
INITMPGRID(), first form:The first
USE MODMPASFIO
setenv FNAME path for input initialization file- for the first form
setenv MPAS_VERBOSE T or F- to turn on verbose logging
setenv MPAS_CHKFILL T or F- to turn on input-data checking (to detect "gaps" in the input data)
INITMPGRID()form is the "usual" one; it initializes all of the internal data structures from an input MPAS-format fileFNAME, and must be called before one calls any of the other public routines exceptSPHEREDIST(). This is the usual form (and the form you should use): the MPAS grid specification is very complicated, and the MPAS-model-initialization software will get it correct for you; you should not normally go around re-inventing this wheel.The second
INITMPGRID()form is a for use with Fortran-90 keyword arguments, where a complete MPAS-grid specification is provided by arguments instead of from an MPAS-file header.
This is not recommended (after all, the MPAS-initialization program's authors have done this extremely tedious operation once—why re-invent that wheel?) Note that for this form, you must correctly specify all of these arguments, and should do so by name, e.g.FLAG = INITMPGRID( & NCELLS = <value> & ... REARTH1 = <value> & )
SUBROUTINE SHUTMPGRID()
SUBROUTINE SHUTMPGRID( )
Flush and close all files opened byBack to ContentsINITMPGRID(),OPENMPAS()orCREATEMPAS(), and deallocate and clear internal state variables and other data structures.Preconditions:
Should be called prior to program termination.
USE MODMPASFIO
OPENMPAS()
LOGICAL FUNCTION OPENMPAS( FNAME,FMODE )
CHARACTER*(*), INTENT(IN ) :: FNAME !! logical file name
INTEGER , INTENT(IN ) :: FMODE !! mode: FSREAD3 or FSRDWR3
Preconditions:Back to ContentsThis routine opens an existing MPAS-format netCDF file either for read-only or read-write,
USE MODMPASFIO
- Call
INITMPGRID()before callingOPENMPAS()
setenv FNAME path for input data file
FMODE = FSREAD3, orFMODE = FSRDWR3, respectively.
CREATEMPAS()
This is a generic LOGICAL FUNCTION with two forms:
LOGICAL FUNCTION CREATEMPAS( FNAME, &
NVARS, VNAMES, VTYPES, VNDIMS, &
VDNAME, VUNITS )
LOGICAL FUNCTION CREATEMPAS( FNAME, NDIMS, DNAMES, DSIZES, &
NVARS, VNAMES, VTYPES, VNDIMS, &
VDNAME, VUNITS )
CHARACTER*(*), INTENT(IN ) :: FNAME !! logical file name
INTEGER , INTENT(IN ) :: NDIMS !! number of dimensions used
CHARACTER*(*), INTENT(IN ) :: DNAMES( NDIMS ) !! dimension-names
INTEGER , INTENT(IN ) :: DSIZES( NDIMS ) !! dimension-values
INTEGER , INTENT(IN ) :: NVARS !! number of (extra) output variables
CHARACTER*(*), INTENT(IN ) :: VNAMES( NVARS ) !! variable-names
INTEGER , INTENT(IN ) :: VTYPES( NVARS ) !! variabale-type M3REAL, etc...)
INTEGER , INTENT(IN ) :: VNDIMS( NVARS ) !! rank (number of dimensions)
CHARACTER*(*), INTENT(IN ) :: VDNAME( NDIMS,NVARS ) !! names for dimensions used for the variables
CHARACTER*(*), INTENT(IN ), OPTIONAL :: VUNITS( NVARS ) !! variable-units
Preconditions:Back to ContentsThis routine opens/creates a new MPAS-format netCDF file with the given (dimension+variables) specification, sets up all the standard dimensions and variables as well as any extra dimensions and variables specified by the argument-list, and writes the standard MPAS variables to that file. If (recommended!)
USE MODMPASFIO
- Call
INITMPGRID()before callingCREATEMPAS()
setenv FNAME path for output data file
FNAMEmust not already exist: this is "create a new file."
OPTIONALargumentVUNITSis present, its values will be written to the file-header as additional metadata, thereby following "NetCDF Best Practices".The first form (without the dimension specifications) assumes that all the dimension-specifications are taken from the file used by
INITMPGRID(), i.e., matching the module's internal grid description variables described above.NOTE: See also MODULE MODNCFIO routine
CREATENC()for reading "raw" netCDF file descriptions.DESCNCVAR()performs an open-read-close sequence, and therefore should not be used after anINITMPGRID()orOPENMPAS()call for the same file.
DESCMPAS()
This is a generic LOGICAL FUNCTION with two forms (these
are necessary because while "NetCDF Best Practices" specifies
that the units attribute should be available for all variables,
the MPAS-netCDF specification does not, and many (most?) MPAS-netCDF
files do not provide it.;-( )
LOGICAL FUNCTION DESCMPAS( FNAME, NRECS, NVARS, VNAMES, VTYPES, &
VNDIMS, VDSIZE, VDNAME )
LOGICAL FUNCTION DESCMPAS( FNAME, NRECS, NVARS, VNAMES, VTYPES, VUNITS, &
VNDIMS, VDSIZE, VDNAME )
CHARACTER*(*), INTENT(IN ) :: FNAME !! logical file name
INTEGER , INTENT( OUT) :: NVARS !! number of (extra) output variables
INTEGER , INTENT( OUT) :: NRECS !! number of time steps
CHARACTER*(*), INTENT( OUT) :: VNAMES( MXVARS3 ) !! variable-names
INTEGER , INTENT( OUT) :: VTYPES( MXVARS3 ) !! variable-type M3REAL, etc...)
CHARACTER*(*), INTENT( OUT) :: VUNITS( NVARS ) !! variable-units
INTEGER , INTENT( OUT) :: VNDIMS( MXVARS3 ) !! rank (number of dimensions)
INTEGER , INTENT( OUT) :: VDSIZE( 7,MXVARS3 ) !! list of dimensions
CHARACTER*(*), INTENT( OUT) :: VDNAME( 7,MXVARS3 ) !! names for dimensions used for the variables
Preconditions:Back to Contents
setenv FNAME path for data file
USE MODMPASFIO
- Call
INITMPGRID()and eitherOPENMPAS(FNAME...)orCREATEMPAS(FNAME...)before callingDESCMPAS(FNAME...)
This routine returns a variables-description for the file with logical name
FNAME. If the version with argumentVUNITSis used, but the units-attribute is not available in the file, the valueCMISS3="????????????????"is returned in that argument.NOTE: See also MODULE MODNCFIO routine
DESCNCVAR()for reading "raw" netCDF file descriptions.DESCNCVAR()performs an open-read-close sequence, and therefore should not be used after anINITMPGRID()orOPENMPAS()call for the same file.
READMPSTEPS()
LOGICAL FUNCTION READMPSTEPS( FNAME, VNAME, NSTEPS, MXSTEP, JDATES, JTIMES )
CHARACTER*(*), INTENT(IN ) :: FNAME !! logical file name
CHARACTER*(*), INTENT(IN ) :: VNAME !! time-variable name
INTEGER , INTENT(IN ) :: MXSTEP !! array-dimension
INTEGER , INTENT( OUT) :: NSTEPS !! number of actual timesteps
INTEGER , INTENT( OUT) :: JDATES( MXSTEP ) !! dates YYYYDDD
INTEGER , INTENT( OUT) :: JTIMES( MXSTEP ) !! times HHMMSS
Preconditions:Back to Contents
setenv FNAME path for data file
USE MODMPASFIO
- Call
INITMPGRID()and eitherOPENMPAS(FNAME...)orCREATEMPAS(FNAME...)before callingREADMPSTEPS().
- Date-&-time variable
VNAMEin fileFNAMEmust beCHARACTER*64, formatted according to MPAS ASCIIYYYY-MM-DD_HH:MM[:SS]year-month-day-hour-minute-seconds conventions (seconds field optional).This routine reads the number of steps and the array of date-&-time strings for the time-variable
VNAME(oftenxTimeorxtime) from the fileFNAME, converts them from MPAS ASCII-string convention to I/O API conventionINTEGER JDATE=YYYYDDD,JTIME=HHMMSS, and returns the number of steps asNSTEPS, and the dates and times as arraysJDATES(NSTEPS)andJTIMES(NSTEPS).NOTE that there is no MPAS-model guarantee that these form a regular time step sequence, nor that particular variables do not have "gaps". In order to get timestep-number
ISTEPfor use withREADMPAS(FNAME,ISTEP,VNAME...)orWRITEMPAS(FNAME,ISTEP,VNAME...), you may find it useful to use I/O API 2-tuple lookup-routineFIND2()as follows (where note thatISTEP>0for the case thatJDATE:JTIMEis available, and is negative otherwise):... IF ( .NOT.READMPSTEPS( FNAME, NSTEPS, MXSTEP, JDATES, JTIMES ) ) THEN ... ISTEP = FIND2( JDATE, JTIME, NSTEPS, JDATES, JTIMES ) IF ( ISTEP .GT. 0 ) THEN IF ( .NOT.READMPAS( FNAME, ISTEP, VNAME, ... ) THEN ... ELSE !! file does not have JDATE:JTIME... ...See also
MPSTR2DT()andMPDT2STR(),below.WRITEMPSTEP(), below.READNCVAR()for reading "raw" netCDF.WRITENCVAR()for writing "raw" netCDF.
WRITEMPSTEP()
LOGICAL FUNCTION WRITEMPSTEP( FNAME, VNAME, ISTEP, JDATE, JTIME )
CHARACTER*(*), INTENT(IN ) :: FNAME !! logical file name
CHARACTER*(*), INTENT(IN ) :: VNAME !! time-variable name
INTEGER , INTENT(IN ) :: ISTEP !! timestep-number (1,2,...)
INTEGER , INTENT(IN ) :: JDATE !! date YYYYDDD
INTEGER , INTENT(IN ) :: JTIME !! time HHMMSS
Preconditions:Back to ContentsThis routine writes the date-&-time string for indicated
setenv FNAME path for data file
USE MODMPASFIO
- Call
INITMPGRID()and eitherOPENMPAS(FNAME, FSRDWR3,...)orCREATEMPAS(FNAME...)before callingREADMPSTEPS().
- Date-&-time variable
VNAMEin fileFNAMEmust beCHARACTER*64, formatted according to MPAS ASCIIYYYY-MM-DD_HH:MM[:SS]year-month-day-hour-minute-seconds conventions.
JDATE:JTIMEto the time-variableVNAMEin the fileFNAME, converting from I/O API conventionINTEGER JDATE=YYYYDDD,JTIME=HHMMSSto MPASYYYY-MM-DD_HH:MM:SSASCII-string convention.NOTE The MPAS data structures being as they are, it is your responsibility to calculate the timestep-number
ISTEP. For "regular" time stepping (with a fixed time step), you may find it useful to use I/O API timestep-number routineJSTEP3()as follows:!! SDATE = starting date (YYYYDDD) for the file !! STIME = starting time (HHMMSS) for the file !! TSTEP = time step (HHMMSS) for the file ... ISTEP = JSTEP3( SDATE, STIME, TSTEP, JDATE, JTIME ) IF ( .NOT.WRITEMPSTEP( FNAME, 'xTime', ISTEP, JDATE, JTIME ) ) THEN ...See also
MPSTR2DT()andMPDT2STR(),below.WRITEMPSTEP(), below.READNCVAR()for reading "raw" netCDF.WRITENCVAR()for writing "raw" netCDF.
READMPAS()This is a genericLOGICAL FUNCTIONfor reading either entire variables, or time steps of variables, from MPAS-formatted netCDF files that have been already opened byOPENMPAS()orCREATEMPAS().There are 50 implemented forms, for 0-D (scalar), 1-D, 2-D, 3-D, and 4-D arrays of types
REAL*8("double"),REAL,INTEGER,INTEGER*2("short"), orINTEGER*1("byte"), either without ("whole-variable") or with time stepping:
!! time independent forms:
LOGICAL FUNCTION READMPAS( FNAME, VNAME, DSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, RSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, ISCALAR )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, SSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, BSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, DARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, RARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, IARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, SARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, BARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, DARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, RARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, IARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, SARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, BARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, DARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, RARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, IARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, SARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, BARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, DARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, RARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, IARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, SARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, BARRAY4 )
!! time stepped forms:
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, DSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, RSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, ISCALAR )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, SSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, BSCALAR )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, DARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, RARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, IARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, SARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, BARRAY1 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, DARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, RARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, IARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, SARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, BARRAY2 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, DARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, RARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, IARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, SARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, BARRAY3 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, DARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, RARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, IARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, SARRAY4 )
LOGICAL FUNCTION READMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, BARRAY4 )
CHARACTER*(*), INTENT(IN ) :: FNAME !! logical file name
CHARACTER*(*), INTENT(IN ) :: VNAME !! logical file name
INTEGER , INTENT(IN ) :: ISTEP !! time step number (1, 2, ...)
INTEGER , INTENT(IN ) :: NDIM1, NDIM2, NDIM3, NDIM4
REAL*8 , INTENT( OUT) :: DSCALAR
REAL , INTENT( OUT) :: RSCALAR
INTEGER , INTENT( OUT) :: ISCALAR
INTEGER*2 , INTENT( OUT) :: SSCALAR
INTEGER*1 , INTENT( OUT) :: BSCALAR
REAL*8 , INTENT( OUT) :: DARRAY1( NDIM1 )
REAL , INTENT( OUT) :: RARRAY1( NDIM1 )
INTEGER , INTENT( OUT) :: IARRAY1( NDIM1 )
INTEGER*2 , INTENT( OUT) :: SARRAY1( NDIM1 )
INTEGER*1 , INTENT( OUT) :: BARRAY1( NDIM1 )
REAL*8 , INTENT( OUT) :: DARRAY2( NDIM1, NDIM2 )
REAL , INTENT( OUT) :: RARRAY2( NDIM1, NDIM2 )
INTEGER , INTENT( OUT) :: IARRAY2( NDIM1, NDIM2 )
INTEGER*2 , INTENT( OUT) :: SARRAY2( NDIM1, NDIM2 )
INTEGER*1 , INTENT( OUT) :: BARRAY2( NDIM1, NDIM2 )
REAL*8 , INTENT( OUT) :: DARRAY3( NDIM1, NDIM2, NDIM3 )
REAL , INTENT( OUT) :: RARRAY3( NDIM1, NDIM2, NDIM3 )
INTEGER , INTENT( OUT) :: IARRAY3( NDIM1, NDIM2, NDIM3 )
INTEGER*2 , INTENT( OUT) :: SARRAY3( NDIM1, NDIM2, NDIM3 )
INTEGER*1 , INTENT( OUT) :: BARRAY3( NDIM1, NDIM2, NDIM3 )
REAL*8 , INTENT( OUT) :: DARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
REAL , INTENT( OUT) :: RARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
INTEGER , INTENT( OUT) :: IARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
INTEGER*2 , INTENT( OUT) :: SARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
INTEGER*1 , INTENT( OUT) :: BARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
Preconditions:Back to ContentsNOTE 1 The MPAS data structures being as they are, it is your responsibility to calculate the timestep-number
setenv FNAME <path for data file>
USE MODMPASFIO
- Call
INITMPGRID()and eitherOPENMPAS(FNAME...)orCREATEMPAS(FNAME...)before callingREADMPAS().
setenv MPAS_CHKFILL Yto enable input-data "hole-detection"
ISTEP. You may find it useful to use READMPSTEPS() together with I/O API 2-tuple lookup-routineFIND2(), or I/O API timestep-number routineJSTEP3()(for the regular-timesteps-only case), as follows:!! SDATE = starting date (YYYYDDD) !! STIME = starting time (HHMMSS) !! TSTEP = time step (HHMMSS) ... !! general case: IF ( .NOT.READMPSTEPS( FNAME, NSTEPS, MXSTEP, JDATES, JTIMES ) ) THEN ... ISTEP = FIND2( JDATE, JTIME, NSTEPS, JDATES, JTIMES ) IF ( ISTEP .GT. 0 ) THEN IF ( .NOT.READMPAS( FNAME, ISTEP, VNAME, ... ) THEN ... ELSE !! file does not have JDATE:JTIME... ... !! regular-timestep-case only: ISTEP = JSTEP3( SDATE, STIME, TSTEP, JDATE, JTIME ) IF ( .NOT.READMPAS( FNAME, ISTEP, VNAME, ... ...NOTE 2 It is also your responsibility to detect "gaps" in the input data, since the internal MPAS data structures are inadequate for doing so.will turn on such checking from withinsetenv MPAS_CHKFILL YREADMPAS(), at some computational costSee also
READMPSTEPS(), aboveWRITEMPSTEP(), above.READNCVAR()for reading "raw" netCDF.WRITENCVAR()for writing "raw" netCDF.
WRITEMPAS()This is a genericLOGICAL FUNCTIONfor writing either entire variables, or time steps of variables, to MPAS-formatted netCDF files that have been already opened byCREATEMPAS(), or byOPENMPAS()in read-write (FMODE=FSRDWR3) mode.There are 50 implemented forms, for 0-D (scalar), 1-D, 2-D, 3-D, and 4-D arrays of types
REAL*8("double"),REAL,INTEGER,INTEGER*2("short"), orINTEGER*1("byte"), either without ("whole-variable") or with time stepping:
!! time independent forms:
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, DSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, RSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, ISCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, SSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, BSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, DARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, RARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, IARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, SARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, BARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, DARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, RARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, IARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, SARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, BARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, DARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, RARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, IARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, SARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, BARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, DARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, RARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, IARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, SARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, BARRAY4 )
!! time stepped forms:
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, DSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, RSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, ISCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, SSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, BSCALAR )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, DARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, RARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, IARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, SARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, BARRAY1 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, DARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, RARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, IARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, SARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, BARRAY2 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, DARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, RARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, IARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, SARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, BARRAY3 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, DARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, RARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, IARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, SARRAY4 )
LOGICAL FUNCTION WRITEMPAS( FNAME, ISTEP, VNAME, NDIM1, NDIM2, NDIM3, NDIM4, BARRAY4 )
CHARACTER*(*), INTENT(IN ) :: FNAME !! logical file name
CHARACTER*(*), INTENT(IN ) :: VNAME !! logical file name
INTEGER , INTENT(IN ) :: ISTEP !! time step number (1, 2, ...)
INTEGER , INTENT(IN ) :: NDIM1, NDIM2, NDIM3, NDIM4
REAL*8 , INTENT(IN ) :: DSCALAR
REAL , INTENT(IN ) :: RSCALAR
INTEGER , INTENT(IN ) :: ISCALAR
INTEGER*2 , INTENT(IN ) :: SSCALAR
INTEGER*1 , INTENT(IN ) :: BSCALAR
REAL*8 , INTENT(IN ) :: DARRAY1( NDIM1 )
REAL , INTENT(IN ) :: RARRAY1( NDIM1 )
INTEGER , INTENT(IN ) :: IARRAY1( NDIM1 )
INTEGER*2 , INTENT(IN ) :: SARRAY1( NDIM1 )
INTEGER*1 , INTENT(IN ) :: BARRAY1( NDIM1 )
REAL*8 , INTENT(IN ) :: DARRAY2( NDIM1, NDIM2 )
REAL , INTENT(IN ) :: RARRAY2( NDIM1, NDIM2 )
INTEGER , INTENT(IN ) :: IARRAY2( NDIM1, NDIM2 )
INTEGER*2 , INTENT(IN ) :: SARRAY2( NDIM1, NDIM2 )
INTEGER*1 , INTENT(IN ) :: BARRAY2( NDIM1, NDIM2 )
REAL*8 , INTENT(IN ) :: DARRAY3( NDIM1, NDIM2, NDIM3 )
REAL , INTENT(IN ) :: RARRAY3( NDIM1, NDIM2, NDIM3 )
INTEGER , INTENT(IN ) :: IARRAY3( NDIM1, NDIM2, NDIM3 )
INTEGER*2 , INTENT(IN ) :: SARRAY3( NDIM1, NDIM2, NDIM3 )
INTEGER*1 , INTENT(IN ) :: BARRAY3( NDIM1, NDIM2, NDIM3 )
REAL*8 , INTENT(IN ) :: DARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
REAL , INTENT(IN ) :: RARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
INTEGER , INTENT(IN ) :: IARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
INTEGER*2 , INTENT(IN ) :: SARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
INTEGER*1 , INTENT(IN ) :: BARRAY4( NDIM1, NDIM2, NDIM3, NDIM4 )
Preconditions:Back to ContentsNOTES
setenv FNAME <path for data file>
USE MODMPASFIO
- Call
INITMPGRID()and eitherOPENMPAS(FNAME,FSRDWR3)orCREATEMPAS(FNAME...)before callingREADMPAS().
See also
- The MPAS data structures being as they are, it is your responsibility to calculate (or keep track of) the timestep-number
ISTEP. If you are using "regular" time stepping (i.e., with a constant time step), you may find it useful to use I/O API timestep-number routineJSTEP3()as follows (where note thatISTEP<0whenJDATE:JTIMEis not on the regular time step sequence defined bySDATE:STIME:TSTEP).:!! SDATE = starting date (YYYYDDD) for the file !! STIME = starting time (HHMMSS) for the file !! TSTEP = time step (HHMMSS) for the file ... ISTEP = JSTEP3( SDATE, STIME, TSTEP, JDATE, JTIME ) IF ( ISTEP .GT. 0 ) THEN IF ( .NOT.WRITEMPAS( FNAME, ISTEP, ... ... ELSE ...
- For per-time-step
WRITEMPAS()calls, when all the calls for a particular time step are complete, you should useWRITEMPSTEP()to record the related MPAS-formatted-ASCII date and time to the time-variable in the file.
- For time dependent variables, it is most efficient if you will write all the variables for a given time step in consecutive order, before going on to the next time step.
It is your responsibility to ensure that there are no "gaps" in the timesteps for the variables you write: MPAS-format internal data structures are inadequate to deal with that situation, if you (accidentally?) introduce such gaps.
CREATEMPAS()writes the "standard" variables to this file, so that MPAS-mandated task is already automated for you. You only need to useWRITEMPAS()for the "extra" variables you have defined.
- See also MODULE MODNCFIO routines
DESCNCVAR()for reading netCDF file descriptions, andREADNCVAR()for reading variables, or time steps of variables, from "raw" netCDF files.
Note that these routines do a open-action-close sequence within each routine, and so are not compatible with the corresponding MPAS-file operations provided by this module (which presumes that once files are opened, they stay opened until you callSHUTMPGRID()); you should only call them prior to usingOPENMPAS()(or you may have obscure "double-free" errors).
READMPSTEPS(), aboveWRITEMPSTEP(), above.READNCVAR()for reading "raw" netCDF.WRITENCVAR()for writing "raw" netCDF.
FINDCELL() and FINDVRTX()
These are generic LOGICAL FUNCTIONs with two forms each:
INTEGER FUNCTION FINDCELL( ALATF, ALONF )
INTEGER FUNCTION FINDCELL( ALATD, ALOND )
INTEGER FUNCTION FINDVRTX( ALATF, ALONF )
INTEGER FUNCTION FINDVRTX( ALATD, ALOND )
REAL , INTENT( IN ) :: ALATF, ALONF
REAL(8), INTENT( IN ) :: ALATD, ALOND
Back to ContentsThese routines are OpenMP thread-safe.
FINDCELL():- Find the subscript for the (unstructured) MPAS-grid cell that contains the point
(ALAT,ALON).
FINDVRTX():- Find the subscript for the vertex at the center of the dual-mesh cell that contains the point
(ALAT,ALON)or equivalently find the subscript for the vertex closest to the point(ALAT,ALON)
Preconditions:
USE MODMPASFIO
- Call
INITMPGRID()before callingFINDCELL()orFINDVRTX()
ARC2MPAS()
This is a generic LOGICAL FUNCTION with four forms:
LOGICAL FUNCTION ARC2MPAS( LAT1, LON1, LAT2, LON2, NMAX, NSEGS, CELLS, WGT2D )
LOGICAL FUNCTION ARC2MPAS( LAT1, LON1, HGT1, LAT2, LON2, HGT2, &
NLAYS, NMAX, ZGRID, &
NSEGS, CELLS, WGT3D )
LOGICAL FUNCTION ARC2MPAS( LAT1, LON1, HGT1, LAT2, LON2, HGT2, &
NLAYS, NMAX, ZGRID, &
NSEGS, CELLS, WGT2D, ZBOTS, ZTOPS, WGT3D )
LOGICAL FUNCTION ARC2MPAS( LAT1, LON1, HGT1, LAT2, LON2, HGT2, &
NLAYS, NMAX, ZGRID, &
NSEGS, CELLS, WGT2D, &
LAYLO, LAYHI, ZBOTS, ZTOPS, WGT3D )
REAL , INTENT(IN ) :: LAT1, LON1, HGT1, LAT2, LON2, HGT2
INTEGER, INTENT(IN ) :: NLAYS, NMAX
REAL , INTENT(IN ) :: ZGRID( NLAYS+1,MPCELLS ) !! layer-surfaces (M a.s.l)
INTEGER, INTENT( OUT) :: NSEGS
INTEGER, INTENT( OUT) :: CELLS( NMAX ) !! horizontal-grid cells
INTEGER, INTENT( OUT) :: LAYLO( NMAX ) !! vertical-grid layer-ranges
INTEGER, INTENT( OUT) :: LAYHI( NMAX )
REAL , INTENT( OUT) :: ZBOTS( NMAX ) !! vertical-grid elevation-ranges (M a.s.l)
REAL , INTENT( OUT) :: ZTOPS( NMAX )
REAL , INTENT( OUT) :: WGT2D( NMAX ) !! horizontal-grid segment-fractions
REAL , INTENT( OUT) :: WGT3D( NLAYS, NMAX ) !! 3-D-grid segment-fractions
Preconditions:Back to Contents
USE MODMPASFIO
- Call
INITMPGRID()before callingARC2MPAS()
- The last two forms, involving 3-D
WGT3D, require that the must supply the MPAS-grid 3-D variableZGRIDcontaining the level-surfaces for the 3-D model.Computes weights for line-segment trajectories, e.g., for mobile source or aircraft emissions.
This routine intersects a segment from
(LAT1,LON1[,HGT1])to(LAT2,LON2[,HGT2])with either the horizontal (2-D) MPAS grid structure or 3-D (horizontal+vertical) grid structure, and returns the corresponding weights and optionally the layer- and elevation-ranges[LAYLO:LAYHI]and[ZBOTS:ZTOPS]for allocating the segment to the cells it intersects (e.g., for use in SMOKE).The last form is the "instrumented 3-D" variant, which provides additional diagnostic data.
MPCELLMATX()MPCELLMATX()is a generic routine for computing sparse "incidence" matrices for either a 1-D or a 2-D array of points.
LOGICAL FUNCTION MPCELLMATX( NP, LAT1F, LON1F, KCELL1 )
LOGICAL FUNCTION MPCELLMATX( NP, LAT1D, LON1D, KCELL1 )
LOGICAL FUNCTION MPCELLMATX( NC, NR, LAT2F, LON2F, KCELL2 )
LOGICAL FUNCTION MPCELLMATX( NC, NR, LAT2D, LON2D, KCELL2 )
INTEGER, INTENT(IN ) :: NP, NC, NR
REAL, INTENT(IN ) :: LAT1F( NP )
REAL, INTENT(IN ) :: LON1F( NP )
REAL(8), INTENT(IN ) :: LAT1D( NP )
REAL(8), INTENT(IN ) :: LON1D( NP )
REAL, INTENT(IN ) :: LAT2F( NC,NR )
REAL, INTENT(IN ) :: LON2F( NC,NR )
REAL(8), INTENT(IN ) :: LAT2D( NC,NR )
REAL(8), INTENT(IN ) :: LON2D( NC,NR )
INTEGER, INTENT( OUT) :: KCELL1( NP ) !! indices
INTEGER, INTENT( OUT) :: KCELL2( NC,NR )
See alsoBack to ContentsMPBARYMATX()aandMPINTERP()below.Returns
TRUEiff the interpolation can be done properly (i.e., if each point is contained in some MPAS-grid cell): the basic property is that the points(LAT*,LON*)are contained in the MPAS-grid cell with subscriptKCELL*.
MPBARYMATX()MPBARYMATX()is a generic routine for computing 3-band sparse matrices for barycentric-linear interpolation to either a 1-D or a 2-D array of points, with optional weighting by the ratio of cell areas.
LOGICAL FUNCTION MPBARYMATX( NP, LAT1F, LON1F, KCELL1, WCELL1F )
LOGICAL FUNCTION MPBARYMATX( NP, LAT1D, LON1D, KCELL1, WCELL1D )
LOGICAL FUNCTION MPBARYMATX( NP, LAT1D, LON1D, KCELL1, WCELL1F )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2F, LON2F, KCELL2, WCELL2F )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2D, LON2D, KCELL2, WCELL2D )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2D, LON2D, KCELL2, WCELL2F )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2F, LON2F, AF, KCELL2, WCELL2F )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2D, LON2D, AD, KCELL2, WCELL2D )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2D, LON2D, AD, KCELL2, WCELL2F )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2F, LON2F, AGRDF, KCELL2, WCELL2F )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2D, LON2D, AGRDD, KCELL2, WCELL2D )
LOGICAL FUNCTION MPBARYMATX( NC, NR, LAT2D, LON2D, AGRDD, KCELL2, WCELL2F )
INTEGER, INTENT(IN ) :: NP, NC, NR
REAL, INTENT(IN ) :: LAT1F( NP )
REAL, INTENT(IN ) :: LON1F( NP )
REAL(8), INTENT(IN ) :: LAT1D( NP )
REAL(8), INTENT(IN ) :: LON1D( NP )
REAL, INTENT(IN ) :: LAT2F( NC,NR )
REAL, INTENT(IN ) :: LON2F( NC,NR )
REAL(8), INTENT(IN ) :: LAT2D( NC,NR )
REAL(8), INTENT(IN ) :: LON2D( NC,NR )
REAL, INTENT(IN ) :: AF !! output-grid cell areas (M^2)
REAL(8), INTENT(IN ) :: AD
REAL, INTENT(IN ) :: AGRDF( NC,NR )
REAL(8), INTENT(IN ) :: AGRDD( NC,NR )
INTEGER, INTENT( OUT) :: KCELL1( 3,NP ) !! indices
INTEGER, INTENT( OUT) :: KCELL2( 3,NC,NR )
REAL, INTENT( OUT) :: WCELL1F( 3,NP ) !! coefficients
REAL(8), INTENT( OUT) :: WCELL1D( 3,NP )
REAL, INTENT( OUT) :: WCELL2F( 3,NC,NR )
REAL(8), INTENT( OUT) :: WCELL2D( 3,NC,NR )
See alsoMPCELLMATX()above aandMPINTERP()below.Returns
TRUEiff the interpolation can be done properly (i.e., if each point is surrounded by MPAS-grid cell-centers).Calls using output-grid cell-area arguments
A = AF,AD,AGRDF,AGRDDre-normalize the output result by the ratio of the relevant cell-areasA*/CAREAso that the emissions-per-grid-cell-per... units for emissions are retained.The interpolation uses barycentric-coordinate linear interpolation on the dual-cell (triangle whose corners are MPAS-cell centers) that contains the requested point(s). See https://codeplea.com/triangular-interpolation for a reference on the interpolation formulas.
OpenMP parallel. Code for the matrix:variable multiplication interpolation algorithm has the following pattern (shown for the 1-D non-layered REAL case):
REAL :: Z( MPCELLS ) !! MPAS variable to be interpolated
REAL :: V( NP ) !! result
INTEGER :: P
...
DO P = 1, NP
K1 = KCELL1( 1,P )
K2 = KCELL1( 2,P )
K3 = KCELL1( 3,P )
W1 = WCELL1F( 1,P )
W2 = WCELL1F( 2,P )
W3 = WCELL1F( 3,P )
V( P ) = W1 * Z( K1 ) + W2 * Z( K2 ) + W3 * Z( K3 )
END DO
Back to Contents
MPBARYMULT()This is a generic routine for multiplying 3-band barycentric-linear interpolation matrices fromMPBARYMATX()(above) with possibly-layered MPAS-grid variables.Note that layered results keep the MPAS-convention layer-subscript-first subscript ordering.
OpenMP parallel.
SUBROUTINE MPBARYMULT( NP, KCELL1, WCELL1F, Z1F, V1F )
SUBROUTINE MPBARYMULT( NP, KCELL1, WCELL1D, Z1D, V1D )
SUBROUTINE MPBARYMULT( NP, NL, KCELL1, WCELL1F, ZLF, V1LF )
SUBROUTINE MPBARYMULT( NP, NL, KCELL1, WCELL1D, ZLD, V1LD )
SUBROUTINE MPBARYMULT( NC, NR, KCELL2, WCELL2F, Z1F, V2F )
SUBROUTINE MPBARYMULT( NC, NR, KCELL2, WCELL2D, Z1D, V2D )
SUBROUTINE MPBARYMULT( NC, NR, NL, KCELL2, WCELL2F, ZLF, V2LF )
SUBROUTINE MPBARYMULT( NC, NR, NL, KCELL2, WCELL2D, ZLD, V2LD )
INTEGER, INTENT(IN ) :: NP, NC, NR, NL
INTEGER, INTENT(IN ) :: KCELL1( 3,NP )
REAL, INTENT(IN ) :: WCELL1F( 3,NP )
REAL(8), INTENT(IN ) :: WCELL1D( 3,NP )
INTEGER, INTENT(IN ) :: KCELL2( 3,NC,NR )
REAL, INTENT(IN ) :: WCELL2F( 3,NC,NR )
REAL(8), INTENT(IN ) :: WCELL2D( 3,NC,NR )
REAL, INTENT(IN ) :: Z1F( MPCELLS ) !! MPAS variable to be interpolated
REAL(8), INTENT(IN ) :: Z1D( MPCELLS )
REAL, INTENT(IN ) :: ZLF( NL,MPCELLS )
REAL(8), INTENT(IN ) :: ZLD( NL,MPCELLS )
REAL, INTENT( OUT) :: V1F( NP ) !! result
REAL(8), INTENT( OUT) :: V1D( NP )
REAL, INTENT( OUT) :: V1LF( NL,NP )
REAL(8), INTENT( OUT) :: V1LD( NL,NP )
REAL, INTENT( OUT) :: V2F( NC,NR )
REAL(8), INTENT( OUT) :: V2D( NC,NR )
REAL, INTENT( OUT) :: V2LF( NC,NR,NL )
REAL(8), INTENT( OUT) :: V2LD( NC,NR,NL )
Back to Contents
MPINTERP()This is a generic routine for interpolating possibly-layered MPAS-grid variables to Lat-Lon specified points, sets of points, or (2-D) grids of points, with optional area-re-weighting for emissions data with "per-cell" units, using barycentric-linear interpolation.
INTEGER forms
LOGICAL FUNCTION MPINTERP( LATR0, LONR0, IZ, IV0 )
LOGICAL FUNCTION MPINTERP( LATD0, LOND0, IZ, IV0 )
LOGICAL FUNCTION MPINTERP( LATR0, LONR0, NLAY, IZL, IVL0 )
LOGICAL FUNCTION MPINTERP( LATD0, LOND0, NLAY, IZL, IVL0 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATR2, LONR2, IZ, IV2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATD2, LOND2, IZ, IV2 )
LOGICAL FUNCTION MPINTERP( NPTS, LATR1, LONR1, IZ, IV1 )
LOGICAL FUNCTION MPINTERP( NPTS, LATD1, LOND1, IZ, IV1 )
LOGICAL FUNCTION MPINTERP( NPTS, LATR1, LONR1, NLAY, IZL, IVL1 )
LOGICAL FUNCTION MPINTERP( NPTS, LATD1, LOND1, NLAY, IZL, IVL1 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATR2, LONR2, NLAY, IZL, IVL2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATD2, LOND2, NLAY, IZL, IVL2 )
Normal REAL forms
LOGICAL FUNCTION MPINTERP( LATR0, LONR0, RZ, RV0 )
LOGICAL FUNCTION MPINTERP( LATD0, LOND0, RZ, RV0 )
LOGICAL FUNCTION MPINTERP( LATR0, LONR0, NLAY, RZL, RVL0 )
LOGICAL FUNCTION MPINTERP( LATD0, LOND0, NLAY, RZL, RVL0 )
LOGICAL FUNCTION MPINTERP( NPTS, LATR1, LONR1, RZ, RV1 )
LOGICAL FUNCTION MPINTERP( NPTS, LATD1, LOND1, RZ, RV1 )
LOGICAL FUNCTION MPINTERP( NPTS, LATR1, LONR1, NLAY, RZL, RVL1 )
LOGICAL FUNCTION MPINTERP( NPTS, LATD1, LOND1, NLAY, RZL, RVL1 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATR2, LONR2, RZ, RV2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATD2, LOND2, RZ, RV2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATR2, LONR2, NLAY, RZL, RVL2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, LATD2, LOND2, NLAY, RZL, RVL2 )
Cell-area re-weighting REAL forms
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AR, LATR2, LONR2, RZ, RV2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AD, LATD2, LOND2, RZ, RV2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AR, LATR2, LONR2, NLAY, RZL, RVL2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AD, LATD2, LOND2, NLAY, RZL, RVL2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AGR, LATR2, LONR2, RZ, RV2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AGD, LATD2, LOND2, RZ, RV2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AGR, LATR2, LONR2, NLAY, RZL, RVL2 )
LOGICAL FUNCTION MPINTERP( NCOL, NROW, AGD, LATD2, LOND2, NLAY, RZL, RVL2 )
INTEGER, INTENT(IN ) :: NPTS
INTEGER, INTENT(IN ) :: NCOL, NROW
INTEGER, INTENT(IN ) :: NLAY !! number of (MPAS-variable) layers
REAL, INTENT(IN ) :: AR !! area of output-grid cell
REAL*8, INTENT(IN ) :: AD
REAL, INTENT(IN ) :: AGR( NCOL,NROW ) !! ...for Lat-Lon or non-regular grids
REAL*8, INTENT(IN ) :: AGD( NCOL,NROW )
REAL, INTENT(IN ) :: LATR0, LONR0
REAL*8, INTENT(IN ) :: LATD0, LOND0
REAL, INTENT(IN ) :: LATR1( NPTS ), LONR1( NPTS )
REAL*8, INTENT(IN ) :: LATD1( NPTS ), LOND1( NPTS )
REAL, INTENT(IN ) :: LATR2( NCOL,NROW ), LONR2( NCOL,NROW )
REAL*8, INTENT(IN ) :: LATD2( NCOL,NROW ), LOND2( NCOL,NROW )
INTEGER, INTENT(IN ) :: IZ( MPCELLS ) !! MPAS-grid variable to be interpolated
REAL, INTENT(IN ) :: RZ( MPCELLS ) !! MPAS-grid variable to be interpolated
INTEGER, INTENT(IN ) :: IZL( NLAY,MPCELLS ) !! layered...
REAL, INTENT(IN ) :: RZL( NLAY,MPCELLS ) !! layered...
REAL, INTENT( OUT) :: RV0 !! result of interpolation
INTEGER, INTENT( OUT) :: IV0
REAL, INTENT( OUT) :: RV1( NPTS )
INTEGER, INTENT( OUT) :: IV1( NPTS )
REAL, INTENT( OUT) :: RV2( NCOL,NROW )
INTEGER, INTENT( OUT) :: IV2( NCOL,NROW )
REAL, INTENT( OUT) :: RVL0( NLAY )
INTEGER, INTENT( OUT) :: IVL0( NLAY )
REAL, INTENT( OUT) :: RVL1( NLAY,NPTS )
INTEGER, INTENT( OUT) :: IVL1( NLAY,NPTS )
REAL, INTENT( OUT) :: RVL2( NLAY,NCOL,NROW )
INTEGER, INTENT( OUT) :: IVL2( NLAY,NCOL,NROW )
ReturnsBack to ContentsTRUEiff the interpolation can be done properly (i.e., if each point is surrounded by MPAS-grid cell-centers). Note thatZvariables are subscripted by the numberMPCELLSof MPAS grid cells, and that layered results keep the MPAS-convention layer-subscript-first subscript ordering.The interpolation uses barycentric-coordinate linear interpolation on the dual-cell (triangle whose corners are MPAS-cell centers) that contains the requested point(s). See https://codeplea.com/triangular-interpolation for a reference on the interpolation formulas.
Calls using output-grid cell-area arguments
A = AR,AD,AGR,AGDre-normalize the output result by the ratio of the relevant cell-areasA*/CAREAso that the emissions-per-grid-cell-per... units for emissions are retained.Note that this routine combines the functionality of
MPBARYMATX()andMPBARYMULT()(above) but always computes the interpolation weights and indices in-line. As a result it simpler to use, but less efficient for repeated uses than those routines.For all but the single-point case,
MPINTERP()is OpenMP-parallel.See also m3tools program mpastom3.
SPHEREDIST()
This is a (potentially-stand-alone) generic FUNCTION
for distance between two Lat-Lon specified or MPAS-Cartesian
specified points on the Earth, with eight forms:
REAL FUNCTION SPHEREDIST( RRADIUS, RLAT1, RLON1, RLAT2, RLON2 )
REAL(8) FUNCTION SPHEREDIST( DRADIUS, DLAT1, DLON1, DLAT2, DLON2 )
REAL FUNCTION SPHEREDIST( RLAT1, RLON1, RLAT2, RLON2 )
REAL(8) FUNCTION SPHEREDIST( DLAT1, DLON1, DLAT2, DLON2 )
REAL FUNCTION SPHEREDIST( RLAT1, RLON1, RLAT2, RLON2 )
REAL(8) FUNCTION SPHEREDIST( DLAT1, DLON1, DLAT2, DLON2 )
REAL FUNCTION SPHEREDIST( RRADIUS, RX1, RY1, RZ1, RX2, RY2, RZ2 )
REAL(8) FUNCTION SPHEREDIST( DRADIUS, DX1, DY1, DZ1, DX2, DY2, DZ2 )
REAL FUNCTION SPHEREDIST( RX1, RY1, RZ1, RX2, RY2, RZ2 )
REAL(8) FUNCTION SPHEREDIST( DX1, DY1, DZ1, DX2, DY2, DZ2 )
REAL , INTENT( IN ) :: RRADIUS !! Earth-radius (Meters)
REAL , INTENT( IN ) :: RLAT1, RLON1, RLAT2, RLON2
REAL(8), INTENT( IN ) :: DRADIUS !! Earth-radius (Meters)
REAL(8), INTENT( IN ) :: RX1, RY1, RZ1, RX2, RY2, RZ2
REAL , INTENT( IN ) :: DX1, DY1, DZ1, DX2, DY2, DZ2
Precondition:Back to ContentsUSE MODMPASFIOThe forms with a
*RADIUSarguments use a caller specified value for the radius of the Earth; the last two use the value fromMODULE—variableREARTH, which either was set byINITMPGRID()or was the module—default value (6370.0d3 meters). The forms with*LAT—*LONarguments expect the point to be specified by Lat—Lon coordinates; the forms with*X*—*Y*—*Z*arguments expect the point to be specified in the MPASX—Y—Z3—D Cartesian coordinate system.Internally, the routines use the (numerically stable) Haversine formula for distance on a sphere; see https://en.wikipedia.org/wiki/Great-circle_distance
See also
INITREARTH()above, if you want to change the Earth-radius value.
MPSTR2DT() and MPDT2STR()
SUBROUTINE MPSTR2DT( CBUF, JDATE, JTIME )
CHARACTER*(*), INTENT(IN ) :: CBUF
INTEGER , INTENT( OUT) :: JDATE, JTIME
SUBROUTINE MPDT2STR( IDATE, ITIME, CSTR )
INTEGER , INTENT(IN ) :: IDATE, ITIME
CHARACTER*(*), INTENT( OUT) :: CSTR
Precondition:Back to ContentsUSE MODMPASFIODates and times on MPAS formatted files are stored in
CHARACTER(LEN=64)strings, in the formatYYYY-MM-DD_HH:MM,YYYY-MM-DD_HH:MM:SSorYYYY-MM-DD_HH:MM:SS.SSSusing the (4-digit) year, and the 2-digit month, the 2-digit day-of-month, the 2-digit hour-of-day, the 2-digit minute-of-hour, and optionally the 2-digit seconds or 6-digit seconds.milliseconds. These subroutines convert back and forth between this MPAS format and the I/O APIINTEGERcoded Julian-date&timeYYYYDDD:HHMMSS.See also
READMPSTEPS()andWRITEMPSTEP(), above.
Permission is granted to copy, distribute and/or modify this document under the terms of the GNU Free Documentation License, Version 1.3 or any later version published by the Free Software Foundation; with no Invariant Sections, no Front-Cover Texts, and no Back-Cover Texts.To: Models-3/EDSS I/O API: The Help Pages
Send comments to
Carlie J. Coats, Jr.
carlie@jyarborough.com