mom_verticalgrid module reference
Provides a transparent vertical ocean grid type and supporting routines.
Data Types
Describes the vertical ocean grid, including unit conversion factors. |
Functions/Subroutines
Allocates and initializes the ocean model vertical grid structure. |
|
Returns the model's thickness units, usually m or kg/m^2. |
|
Returns the model's thickness flux units, usually m^3/s or kg/s. |
|
Returns the model's tracer flux units. |
|
This sets the coordinate data for the "layer mode" of the isopycnal model. |
|
Deallocates the model's vertical grid structure. |
Detailed Description
Provides a transparent vertical ocean grid type and supporting routines.
Type Documentation
- type mom_verticalgrid/verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
- Type fields:
% ke ::
integerThe number of layers/levels in the vertical.% max_depth ::
realThe maximum depth of the ocean [Z ~> m].% g_earth ::
realThe gravitational acceleration [L2 Z-1 T-2 ~> m s-2].% g_earth_z_t2 ::
realThe gravitational acceleration in alternatively rescaled units [Z T-2 ~> m s-2].% rho0 ::
realThe density used in the Boussinesq approximation or nominal density used to convert depths into mass units [R ~> kg m-3].% zaxisunits ::
character(len=40)The units that vertical coordinates are written in.% zaxislongname ::
character(len=40)Coordinate name to appear in files, e.g. “Target Potential Density” or “Height”.% slayer ::
real, dimension(:), allocatableCoordinate values of layer centers, in unscaled units that depend on the vertical coordinate, such as [kg m-3] for an isopycnal or some hybrid coordinates, [m] for a Z* coordinate, or [nondim] for a sigma coordinate.% sinterface ::
real, dimension(:), allocatableCoordinate values on interfaces, in the same unscale units as sLayer [various].% direction ::
integerDirection defaults to 1, positive up.% boussinesq ::
logicalIf true, make the Boussinesq approximation.% semi_boussinesq ::
logicalIf true, do non-Boussinesq pressure force calculations and use mass-based “thicknesses, but use Rho0 to convert layer thicknesses into certain height changes. This only applies if BOUSSINESQ is false.% angstrom_h ::
realA one-Angstrom thickness in the model thickness units [H ~> m or kg m-2].% angstrom_z ::
realA one-Angstrom thickness in the model depth units [Z ~> m].% angstrom_m ::
realA one-Angstrom thickness [m].% h_subroundoff ::
realA thickness that is so small that it can be added to a thickness of Angstrom or larger without changing it at the bit level [H ~> m or kg m-2]. If Angstrom is 0 or exceedingly small, this is negligible compared to 1e-17 m.% dz_subroundoff ::
realA thickness in height units that is so small that it can be added to a vertical distance of Angstrom_Z or 1e-17 m without changing it at the bit level [Z ~> m]. This is the height equivalent of H_subroundoff.% g_prime ::
real, dimension(:), allocatableThe reduced gravity at each interface [L2 Z-1 T-2 ~> m s-2].% rlay ::
real, dimension(:), allocatableThe target coordinate value (potential density) in each layer [R ~> kg m-3].% nkml ::
integerThe number of layers at the top that should be treated as parts of a homogeneous region.% nk_rho_varies ::
integerThe number of layers at the top where the density does not track any target density.% h_to_kg_m2 ::
realA constant that translates thicknesses from the units of thickness to kg m-2 [kg m-2 H-1 ~> kg m-3 or 1].% kg_m2_to_h ::
realA constant that translates thicknesses from kg m-2 to the units of thickness [H m2 kg-1 ~> m3 kg-1 or 1].% m_to_h ::
realA constant that translates distances in m to the units of thickness [H m-1 ~> 1 or kg m-3].% h_to_m ::
realA constant that translates distances in the units of thickness to m [m H-1 ~> 1 or m3 kg-1].% h_to_pa ::
realA constant that translates the units of thickness to pressure [Pa H-1 ~> kg m-2 s-2 or m s-2].% h_to_z ::
realA constant that translates thickness units to the units of depth [Z H-1 ~> 1 or m3 kg-1].% z_to_h ::
realA constant that translates depth units to thickness units depth [H Z-1 ~> 1 or kg m-3].% h_to_rz ::
realA constant that translates thickness units to the units of mass per unit area [R Z H-1 ~> kg m-3 or 1].% rz_to_h ::
realA constant that translates mass per unit area units to thickness units [H R-1 Z-1 ~> m3 kg-2 or 1].% h_to_mks ::
realA constant that translates thickness units to its MKS unit (m or kg m-2) based on GVBoussinesq [m H-1 ~> 1] or [kg m-2 H-1 ~> 1].% m2_s_to_hz_t ::
realThe combination of conversion factors that converts kinematic viscosities in m2 s-1 to the internal units of the kinematic (in Boussinesq mode) or dynamic viscosity [H Z s T-1 m-2 ~> 1 or kg m-3].% hz_t_to_m2_s ::
realThe combination of conversion factors that converts the viscosities from their internal representation into a kinematic viscosity in m2 s-1 [T m2 H-1 Z-1 s-1 ~> 1 or m3 kg-1].% hz_t_to_mks ::
realThe combination of conversion factors that converts the viscosities from their internal representation into their unnscaled MKS units (m2 s-1 or Pa s), depending on whether the model is Boussinesq [T m2 H-1 Z-1 s-1 ~> 1] or [T Pa s H-1 Z-1 ~> 1].
Function/Subroutine Documentation
- subroutine mom_verticalgrid/verticalgridinit(param_file, GV, US)
Allocates and initializes the ocean model vertical grid structure.
- Parameters:
param_file ::
param_file[in] Parameter file handle/typegv :: The container for vertical grid data
us :: [in] A dimensional unit scaling type
- Call to:
- Called from:
- function mom_verticalgrid/get_thickness_units(GV)
Returns the model’s thickness units, usually m or kg/m^2.
- Return:
undefined :: The vertical thickness units
- Parameters:
gv :: [in] The ocean’s vertical grid structure
- Called from:
mom_ale::ale_register_diagsmom_diabatic_driver::diabatic_driver_initmom_geothermal::geothermal_initmom_oda_incupd::init_oda_incupd_diagsmom_diagnostics::mom_diagnostics_initmom_oda_incupd::output_oda_incupd_incmom::register_diagsmom_offline_main::register_diags_offline_transportmom_diagnostics::register_transport_diagsmom::set_restart_fieldsmom_set_visc::set_visc_register_restarts
- function mom_verticalgrid/get_flux_units(GV)
Returns the model’s thickness flux units, usually m^3/s or kg/s.
- Return:
undefined :: The thickness flux units
- Parameters:
gv :: [in] The ocean’s vertical grid structure
- Called from:
mom_dynamics_split_rk2::initialize_dyn_split_rk2mom_dynamics_split_rk2b::initialize_dyn_split_rk2bmom_dynamics_unsplit::initialize_dyn_unsplitmom_dynamics_unsplit_rk2::initialize_dyn_unsplit_rk2mom_diagnostics::mom_diagnostics_initmom_dynamics_split_rk2::register_restarts_dyn_split_rk2mom_dynamics_split_rk2b::register_restarts_dyn_split_rk2bmom_dynamics_unsplit::register_restarts_dyn_unsplitmom_dynamics_unsplit_rk2::register_restarts_dyn_unsplit_rk2mom::set_restart_fields
- function mom_verticalgrid/get_tr_flux_units(GV, tr_units, tr_vol_conc_units, tr_mass_conc_units)
Returns the model’s tracer flux units.
- Return:
undefined :: The model’s flux units for a tracer.
- Parameters:
gv :: [in] The ocean’s vertical grid structure.
tr_units ::
tr_units[in] Units for a tracer, for example Celsius or PSU.tr_vol_conc_units ::
tr_vol_conc_units[in] The concentration units per unit volume, for example if the units are umol m-3, tr_vol_conc_units would be umol.tr_mass_conc_units ::
tr_mass_conc_units[in] The concentration units per unit mass of sea water, for example if the units are mol kg-1, tr_vol_conc_units would be mol.
- Call to:
- Called from:
- subroutine mom_verticalgrid/setverticalgridaxes(Rlay, GV, scale)
This sets the coordinate data for the “layer mode” of the isopycnal model.
- Parameters:
gv :: [inout] The container for vertical grid data
rlay :: [in] The layer target density [R ~> kg m-3]
scale ::
scale[in] A unit scaling factor for Rlay to convert it into the units of sInterface, usually [kg m-3 R-1 ~> 1] when used in layer mode.
- Called from:
- subroutine mom_verticalgrid/verticalgridend(GV)
Deallocates the model’s vertical grid structure.
- Parameters:
gv :: The ocean’s vertical grid structure