mom_verticalgrid module reference

Provides a transparent vertical ocean grid type and supporting routines.

More…

Data Types

verticalgrid_type

Describes the vertical ocean grid, including unit conversion factors.

Functions/Subroutines

verticalgridinit()

Allocates and initializes the ocean model vertical grid structure.

get_thickness_units()

Returns the model’s thickness units, usually m or kg/m^2.

get_flux_units()

Returns the model’s thickness flux units, usually m^3/s or kg/s.

get_tr_flux_units()

Returns the model’s tracer flux units.

setverticalgridaxes()

This sets the coordinate data for the “layer mode” of the isopycnal model.

verticalgridend()

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 [integer] :: The number of layers/levels in the vertical.

  • % max_depth [real] :: The maximum depth of the ocean [Z ~> m].

  • % mks_g_earth [real] :: The gravitational acceleration in unscaled MKS units [m s-2].

  • % g_earth [real] :: The gravitational acceleration [L2 Z-1 T-2 ~> m s-2].

  • % rho0 [real] :: The 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(:),allocatable] :: Coordinate 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(:),allocatable] :: Coordinate values on interfaces, in the same unscale units as sLayer [various].

  • % direction [integer] :: Direction defaults to 1, positive up.

  • % boussinesq [logical] :: If true, make the Boussinesq approximation.

  • % semi_boussinesq [logical] :: If 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 [real] :: A one-Angstrom thickness in the model thickness units [H ~> m or kg m-2].

  • % angstrom_z [real] :: A one-Angstrom thickness in the model depth units [Z ~> m].

  • % angstrom_m [real] :: A one-Angstrom thickness [m].

  • % h_subroundoff [real] :: A 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 [real] :: A 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(:),allocatable] :: The reduced gravity at each interface [L2 Z-1 T-2 ~> m s-2].

  • % rlay [real(:),allocatable] :: The target coordinate value (potential density) in each layer [R ~> kg m-3].

  • % nkml [integer] :: The number of layers at the top that should be treated as parts of a homogeneous region.

  • % nk_rho_varies [integer] :: The number of layers at the top where the density does not track any target density.

  • % h_to_kg_m2 [real] :: A 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 [real] :: A 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 [real] :: A constant that translates distances in m to the units of thickness [H m-1 ~> 1 or kg m-3].

  • % h_to_m [real] :: A constant that translates distances in the units of thickness to m [m H-1 ~> 1 or m3 kg-1].

  • % h_to_pa [real] :: A constant that translates the units of thickness to pressure [Pa H-1 = kg m-1 s-2 H-1 ~> kg m-2 s-2 or m s-2].

  • % h_to_z [real] :: A constant that translates thickness units to the units of depth [Z H-1 ~> 1 or m3 kg-1].

  • % z_to_h [real] :: A constant that translates depth units to thickness units depth [H Z-1 ~> 1 or kg m-3].

  • % h_to_rz [real] :: A 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 [real] :: A constant that translates mass per unit area units to thickness units [H R-1 Z-1 ~> m3 kg-2 or 1].

  • % h_to_mks [real] :: A 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 [real] :: The 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 [real] :: The 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 [real] :: The 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 :: [in] Parameter file handle/type

  • gv :: The container for vertical grid data

  • us :: [in] A dimensional unit scaling type

Call to:

mom_error_handler::mom_error

Called from:

mom_oda_driver_mod::init_oda

function mom_verticalgrid/get_thickness_units(GV) [character(len=48)]

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_diags mom_diabatic_driver::diabatic_driver_init mom_geothermal::geothermal_init mom_oda_incupd::init_oda_incupd_diags mom_diagnostics::mom_diagnostics_init mom_oda_incupd::output_oda_incupd_inc mom::register_diags mom_offline_main::register_diags_offline_transport mom_diagnostics::register_transport_diags mom::set_restart_fields

function mom_verticalgrid/get_flux_units(GV) [character(len=48)]

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_rk2 mom_dynamics_split_rk2b::initialize_dyn_split_rk2b mom_dynamics_unsplit::initialize_dyn_unsplit mom_dynamics_unsplit_rk2::initialize_dyn_unsplit_rk2 mom_diagnostics::mom_diagnostics_init mom_dynamics_split_rk2::register_restarts_dyn_split_rk2 mom_dynamics_split_rk2b::register_restarts_dyn_split_rk2b mom_dynamics_unsplit::register_restarts_dyn_unsplit mom_dynamics_unsplit_rk2::register_restarts_dyn_unsplit_rk2 mom::set_restart_fields

function mom_verticalgrid/get_tr_flux_units(GV, tr_units, tr_vol_conc_units, tr_mass_conc_units) [character(len=48)]

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 :: [in] Units for a tracer, for example Celsius or PSU.

  • 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 :: [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:

mom_error_handler::mom_error

Called from:

mom::initialize_mom

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 :: [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:

mom_coord_initialization::mom_initialize_coord

subroutine mom_verticalgrid/verticalgridend(GV)

Deallocates the model’s vertical grid structure.

Parameters:

gv :: The ocean’s vertical grid structure