mom_variables module reference

Provides transparent structures with groups of MOM6 variables and supporting routines.

More…

Data Types

p3d

A structure for creating arrays of pointers to 3D arrays.

p2d

A structure for creating arrays of pointers to 2D arrays.

surface

Pointers to various fields which may be used describe the surface state of MOM, and which will be returned to the calling program.

thermo_var_ptrs

Pointers to an assortment of thermodynamic fields that may be available, including potential temperature, salinity, heat capacity, and the equation of state control structure.

ocean_internal_state

Pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM.F90.

accel_diag_ptrs

Pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances.

cont_diag_ptrs

Pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances.

vertvisc_type

Vertical viscosities, drag coefficients, and related fields.

bt_cont_type

Container for information about the summed layer transports and how they will vary as the barotropic velocity is changed.

porous_barrier_type

Container for grids modifying cell metric at porous barriers.

Functions/Subroutines

allocate_surface_state()

Allocates the fields for the surface (return) properties of the ocean model.

deallocate_surface_state()

Deallocates the elements of a surface state type.

rotate_surface_state()

Rotate the surface state fields from the input to the model indices.

alloc_bt_cont_type()

Allocates the arrays contained within a BT_cont_type and initializes them to 0.

dealloc_bt_cont_type()

Deallocates the arrays contained within a BT_cont_type.

mom_thermovar_chksum()

Diagnostic checksums on various elements of a thermo_var_ptrs() type for debugging.

Detailed Description

Provides transparent structures with groups of MOM6 variables and supporting routines.

Type Documentation

type mom_variables/p3d

A structure for creating arrays of pointers to 3D arrays.

Type fields:
  • % p [real(:,:,:),pointer] :: A pointer to a 3D array.

type mom_variables/p2d

A structure for creating arrays of pointers to 2D arrays.

Type fields:
  • % p [real(:,:),pointer] :: A pointer to a 2D array.

type mom_variables/surface

Pointers to various fields which may be used describe the surface state of MOM, and which will be returned to the calling program.

Type fields:
  • % sst [real(:,:),allocatable] :: The sea surface temperature [C ~> degC].

  • % sss [real(:,:),allocatable] :: The sea surface salinity [S ~> psu or gSalt/kg].

  • % sfc_density [real(:,:),allocatable] :: The mixed layer density [R ~> kg m-3].

  • % hml [real(:,:),allocatable] :: The mixed layer depth [Z ~> m].

  • % u [real(:,:),allocatable] :: The mixed layer zonal velocity [L T-1 ~> m s-1].

  • % v [real(:,:),allocatable] :: The mixed layer meridional velocity [L T-1 ~> m s-1].

  • % sea_lev [real(:,:),allocatable] :: The sea level [Z ~> m]. If a reduced surface gravity is.

  • % frazil [real(:,:),allocatable] :: The energy needed to heat the ocean column to the freezing point during.

  • % melt_potential [real(:,:),allocatable] :: Instantaneous amount of heat that can be used to melt sea ice [Q R Z ~> J m-2].

  • % ocean_mass [real(:,:),allocatable] :: The total mass of the ocean [R Z ~> kg m-2].

  • % ocean_heat [real(:,:),allocatable] :: The total heat content of the ocean in [C R Z ~> degC kg m-2].

  • % ocean_salt [real(:,:),allocatable] :: The total salt content of the ocean in [1e-3 S R Z ~> kgSalt m-2].

  • % taux_shelf [real(:,:),allocatable] :: The zonal stresses on the ocean under shelves [R L Z T-2 ~> Pa].

  • % tauy_shelf [real(:,:),allocatable] :: The meridional stresses on the ocean under shelves [R L Z T-2 ~> Pa].

  • % t_is_cont [logical] :: If true, the temperature variable SST is actually the conservative temperature in [C ~> degC].

  • % s_is_abss [logical] :: If true, the salinity variable SSS is actually the absolute salinity in [S ~> gSalt kg-1].

  • % tr_fields [type(coupler_2d_bc_type)] :: A structure that may contain an array of named fields describing tracer-related quantities.

  • % arrays_allocated [logical] :: A flag that indicates whether the surface type has had its memory allocated.

type mom_variables/thermo_var_ptrs

Pointers to an assortment of thermodynamic fields that may be available, including potential temperature, salinity, heat capacity, and the equation of state control structure.

Type fields:
  • % t [real(:,:,:),pointer] :: Potential temperature [C ~> degC].

  • % s [real(:,:,:),pointer] :: Salinity [PSU] or [gSalt/kg], generically [S ~> ppt].

  • % p_surf [real(:,:),pointer] :: Ocean surface pressure used in equation of state calculations [R L2 T-2 ~> Pa].

  • % eqn_of_state [type( eos_type ),pointer] :: Type that indicates the equation of state to use.

  • % p_ref [real] :: The coordinate-density reference pressure [R L2 T-2 ~> Pa]. This is the pressure used to calculate Rml from T and S when eqn_of_state is associated.

  • % c_p [real] :: The heat capacity of seawater [Q C-1 ~> J degC-1 kg-1]. When conservative temperature is used, this is constant and exactly 3991.86795711963 J degC-1 kg-1.

  • % t_is_cont [logical] :: If true, the temperature variable tvT is actually the conservative temperature [degC].

  • % s_is_abss [logical] :: If true, the salinity variable tvS is actually the absolute salinity in units of [gSalt kg-1].

  • % min_salinity [real] :: The minimum value of salinity when BOUND_SALINITY=True [S ~> ppt].

  • % spv_avg [real(:,:,:),allocatable] :: The layer averaged in situ specific volume [R-1 ~> m3 kg-1].

  • % valid_spv_halo [integer] :: If positive, the valid halo size for SpV_avg, or if negative SpV_avg is not currently set.

  • % frazil [real(:,:),pointer] :: The energy needed to heat the ocean column to the freezing point since calculate_surface_state was2 last called [Q Z R ~> J m-2].

  • % salt_deficit [real(:,:),pointer] :: The salt needed to maintain the ocean column at a minimum salinity of MIN_SALINITY since the last time that calculate_surface_state was called, [S R Z ~> gSalt m-2].

  • % tempxpme [real(:,:),pointer] :: The net inflow of water into the ocean times the temperature at which this inflow occurs since the last call to calculate_surface_state [C R Z ~> degC kg m-2]. This should be prescribed in the forcing fields, but as it often is not, this is a useful heat budget diagnostic.

  • % internal_heat [real(:,:),pointer] :: Any internal or geothermal heat sources that have been applied to the ocean since the last call to calculate_surface_state [C R Z ~> degC kg m-2].

  • % vart [real(:,:,:),pointer] :: SGS variance of potential temperature [C2 ~> degC2].

  • % vars [real(:,:,:),pointer] :: SGS variance of salinity [S2 ~> ppt2].

  • % covarts [real(:,:,:),pointer] :: SGS covariance of salinity and potential temperature [C S ~> degC ppt].

  • % tr_t [type( tracer_type ),pointer] :: pointer to temp in tracer registry

  • % tr_s [type( tracer_type ),pointer] :: pointer to salinty in tracer registry

type mom_variables/ocean_internal_state

Pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM.F90.

Type fields:
  • % t [real(:,:,:),pointer] :: Pointer to the temperature state variable [C ~> degC].

  • % s [real(:,:,:),pointer] :: Pointer to the salinity state variable [S ~> ppt] (i.e., PSU or g/kg)

  • % u [real(:,:,:),pointer] :: Pointer to the zonal velocity [L T-1 ~> m s-1].

  • % v [real(:,:,:),pointer] :: Pointer to the meridional velocity [L T-1 ~> m s-1].

  • % h [real(:,:,:),pointer] :: Pointer to the layer thicknesses [H ~> m or kg m-2].

  • % uh [real(:,:,:),pointer] :: Pointer to zonal transports [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % vh [real(:,:,:),pointer] :: Pointer to meridional transports [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % cau [real(:,:,:),pointer] :: Pointer to the zonal Coriolis and Advective acceleration [L T-2 ~> m s-2].

  • % cav [real(:,:,:),pointer] :: Pointer to the meridional Coriolis and Advective acceleration [L T-2 ~> m s-2].

  • % pfu [real(:,:,:),pointer] :: Pointer to the zonal Pressure force acceleration [L T-2 ~> m s-2].

  • % pfv [real(:,:,:),pointer] :: Pointer to the meridional Pressure force acceleration [L T-2 ~> m s-2].

  • % diffu [real(:,:,:),pointer] :: Pointer to the zonal acceleration due to lateral viscosity [L T-2 ~> m s-2].

  • % diffv [real(:,:,:),pointer] :: Pointer to the meridional acceleration due to lateral viscosity [L T-2 ~> m s-2].

  • % pbce [real(:,:,:),pointer] :: Pointer to the baroclinic pressure force dependency on free surface movement.

  • % u_accel_bt [real(:,:,:),pointer] :: Pointer to the zonal barotropic-solver acceleration [L T-2 ~> m s-2].

  • % v_accel_bt [real(:,:,:),pointer] :: Pointer to the meridional barotropic-solver acceleration [L T-2 ~> m s-2].

  • % u_av [real(:,:,:),pointer] :: Pointer to zonal velocity averaged over the timestep [L T-1 ~> m s-1].

  • % v_av [real(:,:,:),pointer] :: Pointer to meridional velocity averaged over the timestep [L T-1 ~> m s-1].

  • % u_prev [real(:,:,:),pointer] :: Pointer to zonal velocity at the end of the last timestep [L T-1 ~> m s-1].

  • % v_prev [real(:,:,:),pointer] :: Pointer to meridional velocity at the end of the last timestep [L T-1 ~> m s-1].

type mom_variables/accel_diag_ptrs

Pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances.

Type fields:
  • % diffu [real(:,:,:),pointer] :: Zonal acceleration due to along isopycnal viscosity [L T-2 ~> m s-2].

  • % diffv [real(:,:,:),pointer] :: Meridional acceleration due to along isopycnal viscosity [L T-2 ~> m s-2].

  • % cau [real(:,:,:),pointer] :: Zonal Coriolis and momentum advection accelerations [L T-2 ~> m s-2].

  • % cav [real(:,:,:),pointer] :: Meridional Coriolis and momentum advection accelerations [L T-2 ~> m s-2].

  • % pfu [real(:,:,:),pointer] :: Zonal acceleration due to pressure forces [L T-2 ~> m s-2].

  • % pfv [real(:,:,:),pointer] :: Meridional acceleration due to pressure forces [L T-2 ~> m s-2].

  • % du_dt_visc [real(:,:,:),pointer] :: Zonal acceleration due to vertical viscosity [L T-2 ~> m s-2].

  • % dv_dt_visc [real(:,:,:),pointer] :: Meridional acceleration due to vertical viscosity [L T-2 ~> m s-2].

  • % du_dt_visc_gl90 [real(:,:,:),pointer] :: Zonal acceleration due to GL90 vertical viscosity.

  • % dv_dt_visc_gl90 [real(:,:,:),pointer] :: Meridional acceleration due to GL90 vertical viscosity.

  • % du_dt_str [real(:,:,:),pointer] :: Zonal acceleration due to the surface stress (included.

  • % dv_dt_str [real(:,:,:),pointer] :: Meridional acceleration due to the surface stress (included.

  • % du_dt_dia [real(:,:,:),pointer] :: Zonal acceleration due to diapycnal mixing [L T-2 ~> m s-2].

  • % dv_dt_dia [real(:,:,:),pointer] :: Meridional acceleration due to diapycnal mixing [L T-2 ~> m s-2].

  • % u_accel_bt [real(:,:,:),pointer] :: Pointer to the zonal barotropic-solver acceleration [L T-2 ~> m s-2].

  • % v_accel_bt [real(:,:,:),pointer] :: Pointer to the meridional barotropic-solver acceleration [L T-2 ~> m s-2].

  • % du_other [real(:,:,:),pointer] :: Zonal velocity changes due to any other processes that are not due to any explicit accelerations [L T-1 ~> m s-1].

  • % dv_other [real(:,:,:),pointer] :: Meridional velocity changes due to any other processes that are not due to any explicit accelerations [L T-1 ~> m s-1].

  • % gradkeu [real(:,:,:),pointer] :: gradKEu = - d/dx(u2) [L T-2 ~> m s-2]

  • % gradkev [real(:,:,:),pointer] :: gradKEv = - d/dy(u2) [L T-2 ~> m s-2]

  • % rv_x_v [real(:,:,:),pointer] :: rv_x_v = rv * v at u [L T-2 ~> m s-2]

  • % rv_x_u [real(:,:,:),pointer] :: rv_x_u = rv * u at v [L T-2 ~> m s-2]

  • % diag_hfrac_u [real(:,:,:),pointer] :: Fractional layer thickness at u points [nondim].

  • % diag_hfrac_v [real(:,:,:),pointer] :: Fractional layer thickness at v points [nondim].

  • % diag_hu [real(:,:,:),pointer] :: layer thickness at u points, modulated by the viscous remnant and fractional open areas [H ~> m or kg m-2]

  • % diag_hv [real(:,:,:),pointer] :: layer thickness at v points, modulated by the viscous remnant and fractional open areas [H ~> m or kg m-2]

  • % visc_rem_u [real(:,:,:),pointer] :: viscous remnant at u points [nondim]

  • % visc_rem_v [real(:,:,:),pointer] :: viscous remnant at v points [nondim]

type mom_variables/cont_diag_ptrs

Pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances.

Type fields:
  • % uh [real(:,:,:),pointer] :: Resolved zonal layer thickness fluxes, [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % vh [real(:,:,:),pointer] :: Resolved meridional layer thickness fluxes, [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % uh_smooth [real(:,:,:),pointer] :: Interface height smoothing induced zonal volume fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % vh_smooth [real(:,:,:),pointer] :: Interface height smoothing induced meridional volume fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % uhgm [real(:,:,:),pointer] :: Isopycnal height diffusion induced zonal volume fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % vhgm [real(:,:,:),pointer] :: Isopycnal height diffusion induced meridional volume fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % diapyc_vel [real(:,:,:),pointer] :: The net diapycnal velocity [H T-1 ~> m s-1 or kg m-2 s-1].

type mom_variables/vertvisc_type

Vertical viscosities, drag coefficients, and related fields.

Type fields:
  • % prandtl_turb [real] :: The Prandtl number for the turbulent diffusion that is captured in Kd_shear [nondim].

  • % bbl_thick_u [real(:,:),allocatable] :: The bottom boundary layer thickness at the u-points [Z ~> m].

  • % bbl_thick_v [real(:,:),allocatable] :: The bottom boundary layer thickness at the v-points [Z ~> m].

  • % kv_bbl_u [real(:,:),allocatable] :: The bottom boundary layer viscosity at the u-points [H Z T-1 ~> m2 s-1 or Pa s].

  • % kv_bbl_v [real(:,:),allocatable] :: The bottom boundary layer viscosity at the v-points [H Z T-1 ~> m2 s-1 or Pa s].

  • % ustar_bbl [real(:,:),allocatable] :: The turbulence velocity in the bottom boundary layer at.

  • % tke_bbl [real(:,:),allocatable] :: A term related to the bottom boundary layer source of turbulent kinetic.

  • % taux_shelf [real(:,:),allocatable] :: The zonal stresses on the ocean under shelves [R Z L T-2 ~> Pa].

  • % tauy_shelf [real(:,:),allocatable] :: The meridional stresses on the ocean under shelves [R Z L T-2 ~> Pa].

  • % tbl_thick_shelf_u [real(:,:),allocatable] :: Thickness of the viscous top boundary layer under ice shelves at u-points [Z ~> m].

  • % tbl_thick_shelf_v [real(:,:),allocatable] :: Thickness of the viscous top boundary layer under ice shelves at v-points [Z ~> m].

  • % kv_tbl_shelf_u [real(:,:),allocatable] :: Viscosity in the viscous top boundary layer under ice shelves at u-points [H Z T-1 ~> m2 s-1 or Pa s].

  • % kv_tbl_shelf_v [real(:,:),allocatable] :: Viscosity in the viscous top boundary layer under ice shelves at v-points [H Z T-1 ~> m2 s-1 or Pa s].

  • % nkml_visc_u [real(:,:),allocatable] :: The number of layers in the viscous surface mixed layer at u-points [nondim]. This is not an integer because there may be fractional layers, and it is stored in terms of layers, not depth, to facilitate the movement of the viscous boundary layer with the flow.

  • % nkml_visc_v [real(:,:),allocatable] :: The number of layers in the viscous surface mixed layer at v-points [nondim].

  • % ray_u [real(:,:,:),allocatable] :: The Rayleigh drag velocity to be applied to each layer at u-points [H T-1 ~> m s-1 or Pa s m-1].

  • % ray_v [real(:,:,:),allocatable] :: The Rayleigh drag velocity to be applied to each layer at v-points [H T-1 ~> m s-1 or Pa s m-1].

  • % mld [real(:,:),pointer] :: Instantaneous active mixing layer depth [Z ~> m].

  • % sfc_buoy_flx [real(:,:),pointer] :: Surface buoyancy flux (derived) [Z2 T-3 ~> m2 s-3].

  • % kd_shear [real(:,:,:),pointer] :: The shear-driven turbulent diapycnal diffusivity at the interfaces between layers in tracer columns [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • % kv_shear [real(:,:,:),pointer] :: The shear-driven turbulent vertical viscosity at the interfaces between layers in tracer columns [H Z T-1 ~> m2 s-1 or Pa s].

  • % kv_shear_bu [real(:,:,:),pointer] :: The shear-driven turbulent vertical viscosity at the interfaces between layers in corner columns [H Z T-1 ~> m2 s-1 or Pa s].

  • % kv_slow [real(:,:,:),pointer] :: The turbulent vertical viscosity component due to “slow” processes (e.g., tidal, background, convection etc) [H Z T-1 ~> m2 s-1 or Pa s].

  • % tke_turb [real(:,:,:),pointer] :: The turbulent kinetic energy per unit mass at the interfaces [Z2 T-2 ~> m2 s-2]. This may be at the tracer or corner points.

type mom_variables/bt_cont_type

Container for information about the summed layer transports and how they will vary as the barotropic velocity is changed.

Type fields:
  • % fa_u_ee [real(:,:),allocatable] :: The effective open face area for zonal barotropic transport drawing from locations far to the east [H L ~> m2 or kg m-1].

  • % fa_u_e0 [real(:,:),allocatable] :: The effective open face area for zonal barotropic transport drawing from nearby to the east [H L ~> m2 or kg m-1].

  • % fa_u_w0 [real(:,:),allocatable] :: The effective open face area for zonal barotropic transport drawing from nearby to the west [H L ~> m2 or kg m-1].

  • % fa_u_ww [real(:,:),allocatable] :: The effective open face area for zonal barotropic transport drawing from locations far to the west [H L ~> m2 or kg m-1].

  • % ubt_ww [real(:,:),allocatable] :: uBT_WW is the barotropic velocity [L T-1 ~> m s-1], beyond which the marginal open face area is FA_u_WW. uBT_WW must be non-negative.

  • % ubt_ee [real(:,:),allocatable] :: uBT_EE is a barotropic velocity [L T-1 ~> m s-1], beyond which the marginal open face area is FA_u_EE. uBT_EE must be non-positive.

  • % fa_v_nn [real(:,:),allocatable] :: The effective open face area for meridional barotropic transport drawing from locations far to the north [H L ~> m2 or kg m-1].

  • % fa_v_n0 [real(:,:),allocatable] :: The effective open face area for meridional barotropic transport drawing from nearby to the north [H L ~> m2 or kg m-1].

  • % fa_v_s0 [real(:,:),allocatable] :: The effective open face area for meridional barotropic transport drawing from nearby to the south [H L ~> m2 or kg m-1].

  • % fa_v_ss [real(:,:),allocatable] :: The effective open face area for meridional barotropic transport drawing from locations far to the south [H L ~> m2 or kg m-1].

  • % vbt_ss [real(:,:),allocatable] :: vBT_SS is the barotropic velocity, [L T-1 ~> m s-1], beyond which the marginal open face area is FA_v_SS. vBT_SS must be non-negative.

  • % vbt_nn [real(:,:),allocatable] :: vBT_NN is the barotropic velocity, [L T-1 ~> m s-1], beyond which the marginal open face area is FA_v_NN. vBT_NN must be non-positive.

  • % h_u [real(:,:,:),allocatable] :: An effective thickness at zonal faces, taking into account the effects of vertical viscosity and fractional open areas [H ~> m or kg m-2]. This is primarily used as a non-normalized weight in determining the depth averaged accelerations for the barotropic solver.

  • % h_v [real(:,:,:),allocatable] :: An effective thickness at meridional faces, taking into account the effects of vertical viscosity and fractional open areas [H ~> m or kg m-2]. This is primarily used as a non-normalized weight in determining the depth averaged accelerations for the barotropic solver.

  • % pass_polarity_bt [type(group_pass_type)] :: Structure for polarity group halo updates.

  • % pass_fa_uv [type(group_pass_type)] :: Structure for face area group halo updates.

type mom_variables/porous_barrier_type

Container for grids modifying cell metric at porous barriers.

Type fields:
  • % por_face_areau [real(:,:,:),allocatable] :: fractional open area of U-faces [nondim]

  • % por_face_areav [real(:,:,:),allocatable] :: fractional open area of V-faces [nondim]

  • % por_layer_widthu [real(:,:,:),allocatable] :: fractional open width of U-faces [nondim]

  • % por_layer_widthv [real(:,:,:),allocatable] :: fractional open width of V-faces [nondim]

Function/Subroutine Documentation

subroutine mom_variables/allocate_surface_state(sfc_state, G, use_temperature, do_integrals, gas_fields_ocn, use_meltpot, use_iceshelves, omit_frazil)

Allocates the fields for the surface (return) properties of the ocean model. Unused fields are unallocated.

Parameters:
  • g :: [in] ocean grid structure

  • sfc_state :: [inout] ocean surface state type to be allocated.

  • use_temperature :: [in] If true, allocate the space for thermodynamic variables.

  • do_integrals :: [in] If true, allocate the space for vertically integrated fields.

  • gas_fields_ocn :: [in] If present, this type describes the ocean

  • use_meltpot :: [in] If true, allocate the space for melt potential

  • use_iceshelves :: [in] If true, allocate the space for the stresses under ice shelves.

  • omit_frazil :: [in] If present and false, do not allocate the space to pass frazil fluxes to the coupler

Called from:

rotate_surface_state

subroutine mom_variables/deallocate_surface_state(sfc_state)

Deallocates the elements of a surface state type.

Parameters:

sfc_state :: [inout] ocean surface state type to be deallocated here.

subroutine mom_variables/rotate_surface_state(sfc_state_in, sfc_state, G, turns)

Rotate the surface state fields from the input to the model indices.

Call to:

allocate_surface_state

Called from:

mom::extract_surface_state mom_ice_shelf::initialize_ice_shelf mom_ice_shelf::shelf_calc_flux mom::step_mom

subroutine mom_variables/alloc_bt_cont_type(BT_cont, G, GV, alloc_faces)

Allocates the arrays contained within a BT_cont_type and initializes them to 0.

Parameters:
  • bt_cont :: The BT_cont_type whose elements will be allocated

  • g :: [in] The ocean’s grid structure

  • gv :: [in] The ocean’s vertical grid structure.

  • alloc_faces :: [in] If present and true, allocate memory for effective face thicknesses.

Call to:

mom_error_handler::mom_error

subroutine mom_variables/dealloc_bt_cont_type(BT_cont)

Deallocates the arrays contained within a BT_cont_type.

Parameters:

bt_cont :: The BT_cont_type whose elements will be deallocated.

subroutine mom_variables/mom_thermovar_chksum(mesg, tv, G, US)

Diagnostic checksums on various elements of a thermo_var_ptrs() type for debugging. type for debugging.

Parameters:
  • mesg :: [in] A message that appears in the checksum lines

  • tv :: [in] A structure pointing to various thermodynamic variables

  • g :: [in] The ocean’s grid structure

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

Called from:

mom_diabatic_driver::diabatic_ale mom_diabatic_driver::diabatic_ale_legacy mom_diabatic_driver::layered_diabatic