mom_interface_heights module reference

Functions for calculating interface heights, including free surface height.

More…

Functions/Subroutines

find_eta_3d()

Calculates the heights of all interfaces between layers, using the appropriate form for consistency with the calculation of the pressure gradient forces.

find_eta_2d()

Calculates the free surface height, using the appropriate form for consistency with the calculation of the pressure gradient forces.

calc_derived_thermo()

Calculate derived thermodynamic quantities for re-use later.

find_col_avg_spv()

Determine the column average specific volumes.

find_rho_bottom()

Determine the in situ density averaged over a specified distance from the bottom, calculating it as the inverse of the mass-weighted average specific volume.

dz_to_thickness_tv()

Converts thickness from geometric height units to thickness units, perhaps via an inversion of the integral of the density in pressure using variables stored in the thermo_var_ptrs type when in non-Boussinesq mode.

dz_to_thickness_eos()

Converts thickness from geometric height units to thickness units, working via an inversion of the integral of the density in pressure when in non-Boussinesq mode.

dz_to_thickness_simple()

Converts thickness from geometric height units to thickness units, perhaps using a simple conversion factor that may be problematic in non-Boussinesq mode.

thickness_to_dz_3d()

Converts layer thicknesses in thickness units to the vertical distance between edges in height units, perhaps by multiplication by the precomputed layer-mean specific volume stored in an array in the thermo_var_ptrs type when in non-Boussinesq mode.

thickness_to_dz_jslice()

Converts a vertical i- / k- slice of layer thicknesses in thickness units to the vertical distance between edges in height units, perhaps by multiplication by the precomputed layer-mean specific volume stored in an array in the thermo_var_ptrs type when in non-Boussinesq mode.

Detailed Description

Functions for calculating interface heights, including free surface height.

Function/Subroutine Documentation

subroutine mom_interface_heights/find_eta_3d(h, tv, G, GV, US, eta, eta_bt, halo_size, dZref)

Calculates the heights of all interfaces between layers, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, these height may be dilated for consistency with the corresponding time-average quantity from the barotropic calculation.

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

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

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

  • h :: [in] Layer thicknesses [H ~> m or kg m-2]

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

  • eta :: [out] layer interface heights [Z ~> m]

  • eta_bt :: [in] optional barotropic variable that gives the “correct” free surface height (Boussinesq) or total water column mass per unit area (non-Boussinesq). This is used to dilate the layer thicknesses when calculating interface heights [H ~> m or kg m-2]. In Boussinesq mode, eta_bt and GbathyT use the same reference height.

  • halo_size :: [in] width of halo points on which to calculate eta.

  • dzref :: [in] The difference in the reference height between GbathyT and eta [Z ~> m]. The default is 0.

subroutine mom_interface_heights/find_eta_2d(h, tv, G, GV, US, eta, eta_bt, halo_size, dZref)

Calculates the free surface height, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, the sea surface height may be adjusted for consistency with the corresponding time-average quantity from the barotropic calculation.

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

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

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

  • h :: [in] Layer thicknesses [H ~> m or kg m-2]

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

  • eta :: [out] free surface height relative to mean sea level (z=0) often [Z ~> m].

  • eta_bt :: [in] optional barotropic variable that gives the “correct” free surface height (Boussinesq) or total water column mass per unit area (non-Boussinesq) [H ~> m or kg m-2]. In Boussinesq mode, eta_bt and GbathyT use the same reference height.

  • halo_size :: [in] width of halo points on which to calculate eta.

  • dzref :: [in] The difference in the reference height between GbathyT and eta [Z ~> m]. The default is 0.

subroutine mom_interface_heights/calc_derived_thermo(tv, h, G, GV, US, halo, debug)

Calculate derived thermodynamic quantities for re-use later.

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

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

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

  • tv :: [inout] A structure pointing to various thermodynamic variables, some of which will be set here.

  • h :: [in] Layer thicknesses [H ~> m or kg m-2].

  • halo :: [in] Width of halo within which to calculate thicknesses

  • debug :: [in] If present and true, write debugging checksums

Call to:

mom_density_integrals::avg_specific_vol

Called from:

mom_ale::ale_regrid_accelerated mom_hybgen_unmix::hybgen_unmix mom::initialize_mom mom_diabatic_driver::layered_diabatic mom_state_initialization::mom_initialize_state mom::step_mom mom::step_mom_thermo mom::step_mom_tracer_dyn mom::step_offline mom_offline_main::update_offline_fields

subroutine mom_interface_heights/find_col_avg_spv(h, SpV_avg, tv, G, GV, US, halo_size)

Determine the column average specific volumes.

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

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

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

  • h :: [in] Layer thicknesses [H ~> m or kg m-2]

  • spv_avg :: [inout] Column average specific volume [R-1 ~> m3 kg-1]

  • tv :: [in] Structure containing pointers to any available thermodynamic fields.

  • halo_size :: [in] width of halo points on which to work

Call to:

mom_error_handler::mom_error

Called from:

mom_dynamics_split_rk2::step_mom_dyn_split_rk2 mom_dynamics_split_rk2b::step_mom_dyn_split_rk2b

subroutine mom_interface_heights/find_rho_bottom(h, dz, pres_int, dz_avg, tv, j, G, GV, US, Rho_bot)

Determine the in situ density averaged over a specified distance from the bottom, calculating it as the inverse of the mass-weighted average specific volume.

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

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

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

  • h :: [in] Layer thicknesses [H ~> m or kg m-2]

  • dz :: [in] Height change across layers [Z ~> m]

  • pres_int :: [in] Pressure at each interface [R L2 T-2 ~> Pa]

  • dz_avg :: [in] The vertical distance over which to average [Z ~> m]

  • tv :: [in] Structure containing pointers to any available thermodynamic fields.

  • j :: [in] j-index of row to work on

  • rho_bot :: [out] Near-bottom density [R ~> kg m-3].

Call to:

mom_eos::average_specific_vol mom_eos::eos_domain mom_error_handler::mom_error

Called from:

mom_set_diffusivity::find_n2 mom_int_tide_input::find_n2_bottom

subroutine mom_interface_heights/dz_to_thickness_tv(dz, tv, h, G, GV, US, halo_size)

Converts thickness from geometric height units to thickness units, perhaps via an inversion of the integral of the density in pressure using variables stored in the thermo_var_ptrs type when in non-Boussinesq mode.

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

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

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

  • dz :: [in] Geometric layer thicknesses in height units [Z ~> m]

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

  • h :: [inout] Output thicknesses in thickness units [H ~> m or kg m-2].

  • halo_size :: [in] Width of halo within which to calculate thicknesses

subroutine mom_interface_heights/dz_to_thickness_eos(dz, Temp, Saln, EoS, h, G, GV, US, halo_size, p_surf)

Converts thickness from geometric height units to thickness units, working via an inversion of the integral of the density in pressure when in non-Boussinesq mode.

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

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

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

  • dz :: [in] Geometric layer thicknesses in height units [Z ~> m]

  • temp :: [in] Input layer temperatures [C ~> degC]

  • saln :: [in] Input layer salinities [S ~> ppt]

  • eos :: [in] Equation of state structure

  • h :: [inout] Output thicknesses in thickness units [H ~> m or kg m-2].

  • halo_size :: [in] Width of halo within which to calculate thicknesses

  • p_surf :: [in] Surface pressures [R L2 T-2 ~> Pa]

Called from:

mom_interface_heights::dz_to_thickness::dz_to_thickness_tv

subroutine mom_interface_heights/dz_to_thickness_simple(dz, h, G, GV, US, halo_size, layer_mode)

Converts thickness from geometric height units to thickness units, perhaps using a simple conversion factor that may be problematic in non-Boussinesq mode.

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

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

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

  • dz :: [in] Geometric layer thicknesses in height units [Z ~> m]

  • h :: [inout] Output thicknesses in thickness units [H ~> m or kg m-2].

  • halo_size :: [in] Width of halo within which to calculate thicknesses

  • layer_mode :: [in] If present and true, do the conversion that is appropriate in pure isopycnal layer mode with no state variables or equation of state. Otherwise use a simple constant rescaling factor and avoid the use of GVRlay.

Called from:

mom_tracer_initialization_from_z::mom_initialize_tracer_from_z

subroutine mom_interface_heights/thickness_to_dz_3d(h, tv, dz, G, GV, US, halo_size)

Converts layer thicknesses in thickness units to the vertical distance between edges in height units, perhaps by multiplication by the precomputed layer-mean specific volume stored in an array in the thermo_var_ptrs type when in non-Boussinesq mode.

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

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

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

  • h :: [in] Input thicknesses in thickness units [H ~> m or kg m-2].

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

  • dz :: [inout] Geometric layer thicknesses in height units [Z ~> m]

  • halo_size :: [in] Width of halo within which to calculate thicknesses

subroutine mom_interface_heights/thickness_to_dz_jslice(h, tv, dz, j, G, GV, halo_size)

Converts a vertical i- / k- slice of layer thicknesses in thickness units to the vertical distance between edges in height units, perhaps by multiplication by the precomputed layer-mean specific volume stored in an array in the thermo_var_ptrs type when in non-Boussinesq mode.

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

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

  • h :: [in] Input thicknesses in thickness units [H ~> m or kg m-2].

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

  • dz :: [inout] Geometric layer thicknesses in height units [Z ~> m]

  • j :: [in] The second (j-) index of the input thicknesses to work with

  • halo_size :: [in] Width of halo within which to calculate thicknesses