mom_grid module reference
Provides the ocean grid type.
Data Types
Ocean grid type. |
Functions/Subroutines
MOM_grid_init initializes the ocean grid array sizes and grid memory. |
|
set_derived_metrics calculates metric terms that are derived from other metrics. |
|
Adcroft_reciprocal(x) = 1/x for |x|>0 or 0 for x=0. |
|
Returns true if the coordinates (x,y) are within the h-cell (i,j) |
|
Store an integer indicating which direction to work on first. |
|
Return global shape of horizontal grid. |
|
Allocate memory used by the |
|
Release memory used by the |
Detailed Description
Grid metrics and their inverses are labelled according to their staggered location on a Arakawa C (or B) grid.
Metrics centered on h- or T-points are labelled T, e.g. dxT is the distance across the cell in the x-direction.
Metrics centered on u-points are labelled Cu (C-grid u location). e.g. dyCu is the y-distance between two corners of a T-cell.
Metrics centered on v-points are labelled Cv (C-grid v location). e.g. dyCv is the y-distance between two -points.
Metrics centered on q-points are labelled Bu (B-grid u,v location). e.g. areaBu is the area centered on a q-point.
The labelling of distances (grid metrics) at various staggered
location on an T-cell and around a q-point.”
Areas centered at T-, u-, v- and q- points are areaT, areaCu, areaCv and areaBu respectively.
The reciprocal of metrics are pre-calculated and also stored in the ocean_grid_type() with a I prepended to the name. For example, with a I prepended to the name. For example, 1./areaT is called IareaT, and 1./dyCv is IdyCv.
Geographic latitude and longitude (or model coordinates if not on a sphere) are stored in geoLatT, geoLonT for T-points. u-, v- and q- point coordinates are follow same pattern of replacing T with Cu, Cv and Bu respectively.
Each location also has a 2D mask indicating whether the entire column is land or ocean. mask2dT is 1 if the column is wet or 0 if the T-cell is land. mask2dCu is 1 if both neighboring columns are ocean, and 0 if either is land. OBCmasku is 1 if both neighboring columns are ocean, and 0 if either is land of if this is OBC point.
Type Documentation
- type mom_grid/ocean_grid_type
Ocean grid type. See
mom_grid()for details. for details.- Type fields:
% domain ::
type(mom_domain_type), pointerOcean model domain.% domain_aux ::
type(mom_domain_type), pointerA non-symmetric auxiliary domain type.% hi ::
type(hor_index_type)Horizontal index ranges.% hid2 ::
type(hor_index_type)Horizontal index ranges for level-2-downsampling.% isc ::
integerThe start i-index of cell centers within the computational domain.% iec ::
integerThe end i-index of cell centers within the computational domain.% jsc ::
integerThe start j-index of cell centers within the computational domain.% jec ::
integerThe end j-index of cell centers within the computational domain.% isd ::
integerThe start i-index of cell centers within the data domain.% ied ::
integerThe end i-index of cell centers within the data domain.% jsd ::
integerThe start j-index of cell centers within the data domain.% jed ::
integerThe end j-index of cell centers within the data domain.% isg ::
integerThe start i-index of cell centers within the global domain.% ieg ::
integerThe end i-index of cell centers within the global domain.% jsg ::
integerThe start j-index of cell centers within the global domain.% jeg ::
integerThe end j-index of cell centers within the global domain.% iscb ::
integerThe start i-index of cell vertices within the computational domain.% iecb ::
integerThe end i-index of cell vertices within the computational domain.% jscb ::
integerThe start j-index of cell vertices within the computational domain.% jecb ::
integerThe end j-index of cell vertices within the computational domain.% isdb ::
integerThe start i-index of cell vertices within the data domain.% iedb ::
integerThe end i-index of cell vertices within the data domain.% jsdb ::
integerThe start j-index of cell vertices within the data domain.% jedb ::
integerThe end j-index of cell vertices within the data domain.% isgb ::
integerThe start i-index of cell vertices within the global domain.% iegb ::
integerThe end i-index of cell vertices within the global domain.% jsgb ::
integerThe start j-index of cell vertices within the global domain.% jegb ::
integerThe end j-index of cell vertices within the global domain.% isd_global ::
integerThe value of isd in the global index space (decomposition invariant).% jsd_global ::
integerThe value of isd in the global index space (decomposition invariant).% idg_offset ::
integerThe offset between the corresponding global and local i-indices.% jdg_offset ::
integerThe offset between the corresponding global and local j-indices.% ke ::
integerThe number of layers in the vertical.% symmetric ::
logicalTrue if symmetric memory is used.% nonblocking_updates ::
logicalIf true, non-blocking halo updates are allowed. The default is .false. (for now).% first_direction ::
integerAn integer that indicates which direction is to be updated first in directionally split parts of the calculation. This can be altered during the course of the run via calls to set_first_direction.% mask2dt ::
real, dimension(:, :), allocatable0 for land points and 1 for ocean points on the h-grid [nondim].% geolatt ::
real, dimension(:, :), allocatableThe geographic latitude at tracer (h) points [degrees_N] or [km] or [m].% geolont ::
real, dimension(:, :), allocatableThe geographic longitude at tracer (h) points [degrees_E] or [km] or [m].% dxt ::
real, dimension(:, :), allocatabledxT is delta x at h points [L ~> m].% idxt ::
real, dimension(:, :), allocatable1/dxT [L-1 ~> m-1].% dyt ::
real, dimension(:, :), allocatabledyT is delta y at h points [L ~> m].% idyt ::
real, dimension(:, :), allocatableIdyT is 1/dyT [L-1 ~> m-1].% areat ::
real, dimension(:, :), allocatableThe area of an h-cell [L2 ~> m2].% iareat ::
real, dimension(:, :), allocatable1/areaT [L-2 ~> m-2].% sin_rot ::
real, dimension(:, :), allocatableThe sine of the angular rotation between the local model grid’s northward.% cos_rot ::
real, dimension(:, :), allocatableThe cosine of the angular rotation between the local model grid’s northward.% mask2dcu ::
real, dimension(:, :), allocatable0 for boundary points and 1 for ocean points on the u grid [nondim].% obcmaskcu ::
real, dimension(:, :), allocatable0 for boundary or OBC points and 1 for ocean points on the u grid [nondim].% geolatcu ::
real, dimension(:, :), allocatableThe geographic latitude at u points [degrees_N] or [km] or [m].% geoloncu ::
real, dimension(:, :), allocatableThe geographic longitude at u points [degrees_E] or [km] or [m].% dxcu ::
real, dimension(:, :), allocatabledxCu is delta x at u points [L ~> m].% idxcu ::
real, dimension(:, :), allocatable1/dxCu [L-1 ~> m-1].% idxcu_obcmask ::
real, dimension(:, :), allocatable1/dxCu or 0 at boundary or OBC points [L-1 ~> m-1].% dycu ::
real, dimension(:, :), allocatabledyCu is delta y at u points [L ~> m].% idycu ::
real, dimension(:, :), allocatable1/dyCu [L-1 ~> m-1].% dy_cu ::
real, dimension(:, :), allocatableThe unblocked lengths of the u-faces of the h-cell [L ~> m].% iareacu ::
real, dimension(:, :), allocatableThe masked inverse areas of u-grid cells [L-2 ~> m-2].% areacu ::
real, dimension(:, :), allocatableThe areas of the u-grid cells [L2 ~> m2].% mask2dcv ::
real, dimension(:, :), allocatable0 for boundary points and 1 for ocean points on the v grid [nondim].% obcmaskcv ::
real, dimension(:, :), allocatable0 for boundary or OBC points and 1 for ocean points on the v grid [nondim].% geolatcv ::
real, dimension(:, :), allocatableThe geographic latitude at v points [degrees_N] or [km] or [m].% geoloncv ::
real, dimension(:, :), allocatableThe geographic longitude at v points [degrees_E] or [km] or [m].% dxcv ::
real, dimension(:, :), allocatabledxCv is delta x at v points [L ~> m].% idxcv ::
real, dimension(:, :), allocatable1/dxCv [L-1 ~> m-1].% dycv ::
real, dimension(:, :), allocatabledyCv is delta y at v points [L ~> m].% idycv ::
real, dimension(:, :), allocatable1/dyCv [L-1 ~> m-1].% idycv_obcmask ::
real, dimension(:, :), allocatable1/dxCv or 0 at boundary or OBC points [L-1 ~> m-1].% dx_cv ::
real, dimension(:, :), allocatableThe unblocked lengths of the v-faces of the h-cell [L ~> m].% iareacv ::
real, dimension(:, :), allocatableThe masked inverse areas of v-grid cells [L-2 ~> m-2].% areacv ::
real, dimension(:, :), allocatableThe areas of the v-grid cells [L2 ~> m2].% porous_dminu ::
real, dimension(:, :), allocatableminimum topographic height (deepest) of U-face [Z ~> m]% porous_dmaxu ::
real, dimension(:, :), allocatablemaximum topographic height (shallowest) of U-face [Z ~> m]% porous_davgu ::
real, dimension(:, :), allocatableaverage topographic height of U-face [Z ~> m]% porous_dminv ::
real, dimension(:, :), allocatableminimum topographic height (deepest) of V-face [Z ~> m]% porous_dmaxv ::
real, dimension(:, :), allocatablemaximum topographic height (shallowest) of V-face [Z ~> m]% porous_davgv ::
real, dimension(:, :), allocatableaverage topographic height of V-face [Z ~> m]% mask2dbu ::
real, dimension(:, :), allocatable0 for boundary points and 1 for ocean points on the q grid [nondim].% geolatbu ::
real, dimension(:, :), allocatableThe geographic latitude at q points [degrees_N] or [km] or [m].% geolonbu ::
real, dimension(:, :), allocatableThe geographic longitude at q points [degrees_E] or [km] or [m].% dxbu ::
real, dimension(:, :), allocatabledxBu is delta x at q points [L ~> m].% idxbu ::
real, dimension(:, :), allocatable1/dxBu [L-1 ~> m-1].% dybu ::
real, dimension(:, :), allocatabledyBu is delta y at q points [L ~> m].% idybu ::
real, dimension(:, :), allocatable1/dyBu [L-1 ~> m-1].% areabu ::
real, dimension(:, :), allocatableareaBu is the area of a q-cell [L2 ~> m2]% iareabu ::
real, dimension(:, :), allocatableIareaBu = 1/areaBu [L-2 ~> m-2].% gridlatt ::
real, dimension(:), pointerThe latitude of T points for the purpose of labeling the output axes,.% gridlatb ::
real, dimension(:), pointerThe latitude of B points for the purpose of labeling the output axes,.% gridlont ::
real, dimension(:), pointerThe longitude of T points for the purpose of labeling the output axes,.% gridlonb ::
real, dimension(:), pointerThe longitude of B points for the purpose of labeling the output axes,.% x_axis_units ::
character(len=40)The units that are used in labeling the x coordinate axes.% y_axis_units ::
character(len=40)The units that are used in labeling the y coordinate axes.% x_ax_unit_short ::
character(len=40)A short description of the x-axis units for documenting parameter units.% y_ax_unit_short ::
character(len=40)A short description of the y-axis units for documenting parameter units.% bathyt ::
real, dimension(:, :), allocatableOcean bottom depth, referenced to Z_ref at tracer points. bathyT is in.% meansl ::
real, dimension(:, :), allocatableSpatially varying time mean sea level, referenced to Z_ref at tracer points.% z_ref ::
realA reference value for all geometric height fields, such as bathyT [Z ~> m].% bathymetry_at_vel ::
logicalIf true, there are separate values for the basin depths at velocity points. Otherwise the effects of of topography are entirely determined from thickness points.% dblock_u ::
real, dimension(:, :), allocatableTopographic depths at u-points at which the flow is blocked [Z ~> m].% dopen_u ::
real, dimension(:, :), allocatableTopographic depths at u-points at which the flow is open at width dy_Cu [Z ~> m].% dblock_v ::
real, dimension(:, :), allocatableTopographic depths at v-points at which the flow is blocked [Z ~> m].% dopen_v ::
real, dimension(:, :), allocatableTopographic depths at v-points at which the flow is open at width dx_Cv [Z ~> m].% coriolisbu ::
real, dimension(:, :), allocatableThe Coriolis parameter at corner points [T-1 ~> s-1].% coriolis2bu ::
real, dimension(:, :), allocatableThe square of the Coriolis parameter at corner points [T-2 ~> s-2].% df_dx ::
real, dimension(:, :), allocatableDerivative d/dx f (Coriolis parameter) at h-points [T-1 L-1 ~> s-1 m-1].% df_dy ::
real, dimension(:, :), allocatableDerivative d/dy f (Coriolis parameter) at h-points [T-1 L-1 ~> s-1 m-1].% areat_global ::
realGlobal sum of h-cell area [L2 ~> m2].% iareat_global ::
realGlobal sum of inverse h-cell area (1/areaT_global) [L-2 ~> m-2].% us ::
type(unit_scale_type), pointerA dimensional unit scaling type.% nblocks ::
integerThe number of sub-PE blocks on this PE.% block ::
type(hor_index_type), dimension(:), pointerIndex ranges for each block.% grid_unit_to_l ::
realA factor that converts a the geoLat and geoLon variables and related variables like len_lat and len_lon into rescaled horizontal distance units on a Cartesian grid, in [L km ~> 1000] or [L m-1 ~> 1] or is 0 for a non-Cartesian grid.% south_lat ::
realThe latitude (or y-coordinate) of the first v-line [degrees_N] or [km] or [m].% west_lon ::
realThe longitude (or x-coordinate) of the first u-line [degrees_E] or [km] or [m].% len_lat ::
realThe latitudinal (or y-coord) extent of physical domain [degrees_N] or [km] or [m].% len_lon ::
realThe longitudinal (or x-coord) extent of physical domain [degrees_E] or [km] or [m].% rad_earth_l ::
realThe radius of the planet in rescaled units [L ~> m].% max_depth ::
realThe maximum depth of the ocean in depth units [Z ~> m].
Function/Subroutine Documentation
- subroutine mom_grid/mom_grid_init(G, param_file, US, HI, global_indexing, bathymetry_at_vel)
MOM_grid_init initializes the ocean grid array sizes and grid memory.
- Parameters:
g :: [inout] The horizontal grid type
param_file ::
param_file[in] Parameter file handleus :: A dimensional unit scaling type
hi :: [in] A hor_index_type for array extents
global_indexing ::
global_indexing[in] If true use global index values instead of having the data domain on each processor start at 1.bathymetry_at_vel ::
bathymetry_at_vel[in] If true, there are separate values for the ocean bottom depths at velocity points. Otherwise the effects of topography are entirely determined from thickness points.
- Call to:
allocate_metricsmom_hor_index::hor_index_initmom_error_handler::mom_errormom_error_handler::mom_mesg- Called from:
mom_oda_driver_mod::init_odamom_ice_shelf::initialize_ice_shelf
- subroutine mom_grid/set_derived_metrics(G, US)
set_derived_metrics calculates metric terms that are derived from other metrics.
- Parameters:
g :: [inout] The horizontal grid structure
us :: [in] A dimensional unit scaling type
- Call to:
- Called from:
- function mom_grid/adcroft_reciprocal(val)
Adcroft_reciprocal(x) = 1/x for |x|>0 or 0 for x=0.
- Parameters:
val ::
val[in] The value being inverted [A].- Return:
undefined :: The Adcroft reciprocal of val [A-1].
- Called from:
- function mom_grid/ispointincell(G, i, j, x, y)
Returns true if the coordinates (x,y) are within the h-cell (i,j)
- Parameters:
g :: [in] Grid type
i ::
i[in] i index of cell to testj ::
j[in] j index of cell to testx ::
x[in] x coordinate of point [degrees_E]y ::
y[in] y coordinate of point [degrees_N]
- subroutine mom_grid/set_first_direction(G, y_first)
Store an integer indicating which direction to work on first.
- Parameters:
g :: [inout] The ocean’s grid structure
y_first ::
y_first[in] The first direction to store
- Called from:
- subroutine mom_grid/get_global_grid_size(G, niglobal, njglobal)
Return global shape of horizontal grid.
- Parameters:
g :: [inout] The horizontal grid type
niglobal ::
niglobal[out] i-index global size of gridnjglobal ::
njglobal[out] j-index global size of grid
- subroutine mom_grid/allocate_metrics(G)
Allocate memory used by the
ocean_grid_type()and related structures. and related structures.- Parameters:
g :: [inout] The horizontal grid type
- Called from:
- subroutine mom_grid/mom_grid_end(G)
Release memory used by the
ocean_grid_type()and related structures. and related structures.- Parameters:
g :: [inout] The horizontal grid type