mom_diabatic_aux module reference

Provides functions for some diabatic processes such as fraxil, brine rejection, tendency due to surface flux divergence.

More…

Data Types

diabatic_aux_cs

Control structure for diabatic_aux.

Functions/Subroutines

make_frazil()

Frazil formation keeps the temperature above the freezing point.

differential_diffuse_t_s()

This subroutine applies double diffusion to T & S, assuming no diapycnal mass fluxes, using a simple tridiagonal solver.

adjust_salt()

This subroutine keeps salinity from falling below a small but positive threshold.

tridiagts()

This is a simple tri-diagonal solver for T and S.

tridiagts_eulerian()

This is a simple tri-diagonal solver for T and S, with mixing across interfaces but no net transfer of mass.

find_uv_at_h()

This subroutine calculates u_h and v_h (velocities at thickness points), optionally using the entrainment amounts passed in as arguments.

set_pen_shortwave()

Estimate the optical properties of the water column and determine the penetrating shortwave radiation by band, extracting the relevant information from the fluxes type and storing it in the optics type for later application.

diagnosemldbydensitydifference()

Diagnose a mixed layer depth (MLD) determined by a given density difference with the surface.

diagnosemldbyenergy()

Diagnose a mixed layer depth (MLD) determined by the depth a given energy value would mix.

applyboundaryfluxesinout()

Update the thickness, temperature, and salinity due to thermodynamic boundary forcing (contained in fluxes type) applied to h, tvT and tvS, and calculate the TKE implications of this heating.

diabatic_aux_init()

This subroutine initializes the parameters and control structure of the diabatic_aux module.

diabatic_aux_end()

This subroutine initializes the control structure and any related memory for the diabatic_aux module.

Detailed Description

This module contains subroutines that apply various diabatic processes. Usually these subroutines are called from the MOM_diabatic module. All of these routines use appropriate limiters or logic to work properly with arbitrary layer thicknesses (including massless layers) and an arbitrarily large timestep.

The subroutine make_frazil facilitates the formation of frazil ice when the ocean water drops below the in situ freezing point by heating the water to the freezing point and accumulating the required heat for exchange with the sea-ice module.

The subroutine adjust_salt adds salt as necessary to keep the salinity above a specified minimum value, and keeps track of the cumulative additions. If the minimum salinity is the natural value of 0, this routine should never do anything.

The subroutine differential_diffuse_T_S solves a pair of tridiagonal equations for the diffusion of temperatures and salinities with differing diffusivities.

The subroutine triDiagTS solves a tridiagonal equations for the evolution of temperatures and salinities due to net entrainment by layers and a diffusion with the same diffusivity.

The subroutine triDiagTS_Eulerian solves a tridiagonal equations for the evolution of temperatures and salinities due to diffusion with the same diffusivity, but no net entrainment.

The subroutine find_uv_at_h interpolates velocities to thickness points, optionally also using tridiagonal equations to solve for the impacts of net entrainment or mixing of momentum between layers.

The subroutine set_pen_shortwave determines the optical properties of the water column and the net shortwave fluxes, and stores them in the optics type, working via calls to set_opacity.

The subroutine diagnoseMLDbyDensityDifference diagnoses a mixed layer depth based on a density difference criterion, and may also estimate the stratification of the water below this diagnosed mixed layer.

The subroutine diagnoseMLDbyEnergy diagnoses a mixed layer depth based on a mixing-energy criterion, as described by Reichl et al., 2022, JGR: Oceans, doi:10.1029/2021JC018140.

The subroutine applyBoundaryFluxesInOut updates the layer thicknesses, temperatures and salinities due to the application of the surface forcing. It may also calculate the implied turbulent kinetic energy requirements for this forcing to be mixed over the model’s finite vertical resolution in the surface layers.

Type Documentation

type mom_diabatic_aux/diabatic_aux_cs

Control structure for diabatic_aux.

Type fields:
  • % do_rivermix [logical] :: Provide additional TKE to mix river runoff at the river mouths to a depth of “rivermix_depth”.

  • % rivermix_depth [real] :: The depth to which rivers are mixed if do_rivermix = T [Z ~> m].

  • % dsalt_frac_max [real] :: An upper limit on the fraction of the salt in a layer that can be lost to the net surface salt fluxes within a timestep [nondim].

  • % reclaim_frazil [logical] :: If true, try to use any frazil heat deficit to to cool the topmost layer down to the freezing point. The default is true.

  • % pressure_dependent_frazil [logical] :: If true, use a pressure dependent freezing temperature when making frazil. The default is false, which will be faster but is inappropriate with ice-shelf cavities.

  • % ignore_fluxes_over_land [logical] :: If true, the model does not check if fluxes are applied over land points. This flag must be used when the ocean is coupled with sea ice and ice shelves and use_ePBL = true.

  • % use_river_heat_content [logical] :: If true, assumes that ice-ocean boundary has provided a river heat content. Otherwise, runoff is added with a temperature of the local SST.

  • % use_calving_heat_content [logical] :: If true, assumes that ice-ocean boundary has provided a calving heat content. Otherwise, calving is added with a temperature of the local SST.

  • % var_pen_sw [logical] :: If true, use one of the CHL_A schemes to determine the e-folding depth of incoming shortwave radiation.

  • % sbc_chl [type(external_field)] :: A handle used in time interpolation of chlorophyll read from a file.

  • % chl_from_file [logical] :: If true, chl_a is read from a file.

  • % do_brine_plume [logical] :: If true, insert salt flux below the surface according to a parameterization by

  • % brine_plume_n [integer] :: The exponent in the brine plume parameterization.

  • % plume_strength [real] :: Fraction of the available brine to take to the bottom of the mixed layer [nondim].

  • % time [type(time_type),pointer] :: A pointer to the ocean model’s clock.

  • % diag [type( diag_ctrl ),pointer] :: Structure used to regulate timing of diagnostic output.

  • % id_createdh [integer] :: Diagnostic ID of mass added to avoid grounding.

  • % id_brine_lay [integer] :: Diagnostic ID of which layer receives the brine.

  • % id_pensw_diag [integer] :: Diagnostic ID of Penetrative shortwave heating (flux convergence)

  • % id_penswflux_diag [integer] :: Diagnostic ID of Penetrative shortwave flux.

  • % id_nonpensw_diag [integer] :: Diagnostic ID of Non-penetrative shortwave heating.

  • % id_chl [integer] :: Diagnostic ID of chlorophyll-A handles for opacity.

  • % createdh [real(:,:),allocatable] :: The amount of volume added in order to avoid grounding [H T-1 ~> m s-1].

  • % pensw_diag [real(:,:,:),allocatable] :: Heating in a layer from convergence of penetrative SW [Q R Z T-1 ~> W m-2].

  • % penswflux_diag [real(:,:,:),allocatable] :: Penetrative SW flux at base of grid layer [Q R Z T-1 ~> W m-2].

  • % nonpensw_diag [real(:,:),allocatable] :: Non-downwelling SW radiation at ocean surface [Q R Z T-1 ~> W m-2].

Function/Subroutine Documentation

subroutine mom_diabatic_aux/make_frazil(h, tv, G, GV, US, CS, p_surf, halo)

Frazil formation keeps the temperature above the freezing point. This subroutine warms any water that is colder than the (currently surface) freezing point up to the freezing point and accumulates the required heat (in [Q R Z ~> J m-2]) in tvfrazil.

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

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

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

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

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

  • cs :: [in] The control structure returned by a previous call to diabatic_aux_init.

  • p_surf :: [in] The pressure at the ocean surface [R L2 T-2 ~> Pa].

  • halo :: [in] Halo width over which to calculate frazil

Call to:

id_clock_frazil

subroutine mom_diabatic_aux/differential_diffuse_t_s(h, T, S, Kd_T, Kd_S, tv, dt, G, GV)

This subroutine applies double diffusion to T & S, assuming no diapycnal mass fluxes, using a simple tridiagonal solver.

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

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

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

  • t :: [inout] Potential temperature [C ~> degC].

  • s :: [inout] Salinity [PSU] or [gSalt/kg], generically [S ~> ppt].

  • kd_t :: [in] The extra diffusivity of temperature due to

  • kd_s :: [in] The extra diffusivity of salinity due to

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

  • dt :: [in] Time increment [T ~> s].

subroutine mom_diabatic_aux/adjust_salt(h, tv, G, GV, CS)

This subroutine keeps salinity from falling below a small but positive threshold. This usually occurs when the ice model attempts to extract more salt then is actually available to it from the ocean.

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

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

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

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

  • cs :: [in] The control structure returned by a previous call to diabatic_aux_init.

subroutine mom_diabatic_aux/tridiagts(G, GV, is, ie, js, je, hold, ea, eb, T, S)

This is a simple tri-diagonal solver for T and S. “Simple” means it only uses arrays hold, ea and eb.

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

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

  • is :: [in] The start i-index to work on.

  • ie :: [in] The end i-index to work on.

  • js :: [in] The start j-index to work on.

  • je :: [in] The end j-index to work on.

  • hold :: [in] The layer thicknesses before entrainment, [H ~> m or kg m-2].

  • ea :: [in] The amount of fluid entrained from the layer above within this time step [H ~> m or kg m-2]

  • eb :: [in] The amount of fluid entrained from the layer below within this time step [H ~> m or kg m-2]

  • t :: [inout] Layer potential temperatures [C ~> degC].

  • s :: [inout] Layer salinities [S ~> ppt].

subroutine mom_diabatic_aux/tridiagts_eulerian(G, GV, is, ie, js, je, hold, ent, T, S)

This is a simple tri-diagonal solver for T and S, with mixing across interfaces but no net transfer of mass.

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

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

  • is :: [in] The start i-index to work on.

  • ie :: [in] The end i-index to work on.

  • js :: [in] The start j-index to work on.

  • je :: [in] The end j-index to work on.

  • hold :: [in] The layer thicknesses before entrainment, [H ~> m or kg m-2].

  • ent :: [in] The amount of fluid mixed across an interface within this time step [H ~> m or kg m-2]

  • t :: [inout] Layer potential temperatures [C ~> degC].

  • s :: [inout] Layer salinities [S ~> ppt].

subroutine mom_diabatic_aux/find_uv_at_h(u, v, h, u_h, v_h, G, GV, US, ea, eb, zero_mix)

This subroutine calculates u_h and v_h (velocities at thickness points), optionally using the entrainment amounts passed in as arguments.

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

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

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

  • u :: [in] The zonal velocity [L T-1 ~> m s-1]

  • v :: [in] The meridional velocity [L T-1 ~> m s-1]

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

  • u_h :: [out] Zonal velocity interpolated to h points [L T-1 ~> m s-1].

  • v_h :: [out] Meridional velocity interpolated to h points [L T-1 ~> m s-1].

  • ea :: [in] The amount of fluid entrained from the layer

  • eb :: [in] The amount of fluid entrained from the layer

  • zero_mix :: [in] If true, do the calculation of u_h and v_h as though ea and eb were being supplied with uniformly zero values.

Call to:

id_clock_uv_at_h

subroutine mom_diabatic_aux/set_pen_shortwave(optics, fluxes, G, GV, US, CS, opacity, tracer_flow_CSp)

Estimate the optical properties of the water column and determine the penetrating shortwave radiation by band, extracting the relevant information from the fluxes type and storing it in the optics type for later application. This routine is effectively a wrapper for set_opacity with added error handling and diagnostics.

Parameters:
  • optics :: An optics structure that has will contain information about shortwave fluxes and absorption.

  • fluxes :: [inout] points to forcing fields unused fields have NULL pointers

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

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

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

  • cs :: Control structure for diabatic_aux

  • opacity :: The control structure for the opacity module.

  • tracer_flow_csp :: A pointer to the control structure organizing the tracer modules.

Call to:

mom_tracer_flow_control::get_chl_from_model

subroutine mom_diabatic_aux/diagnosemldbydensitydifference(id_MLD, h, tv, densityDiff, G, GV, US, diagPtr, id_N2subML, id_MLDsq, dz_subML)

Diagnose a mixed layer depth (MLD) determined by a given density difference with the surface. This routine is appropriate in MOM_diabatic_aux due to its position within the time stepping.

Parameters:
  • g :: [in] Grid type

  • gv :: [in] ocean vertical grid structure

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

  • id_mld :: [in] Handle (ID) of MLD diagnostic

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

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

  • densitydiff :: [in] Density difference to determine MLD [R ~> kg m-3]

  • diagptr :: Diagnostics structure

  • id_n2subml :: [in] Optional handle (ID) of subML stratification

  • id_mldsq :: [in] Optional handle (ID) of squared MLD

  • dz_subml :: [in] The distance over which to calculate N2subML or 50 m if missing [Z ~> m]

Called from:

mom_diabatic_driver::diabatic

subroutine mom_diabatic_aux/diagnosemldbyenergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr)

Diagnose a mixed layer depth (MLD) determined by the depth a given energy value would mix. This routine is appropriate in MOM_diabatic_aux due to its position within the time stepping.

Parameters:
  • id_mld :: [in] Energy output diagnostic IDs

  • g :: [in] Grid type

  • gv :: [in] ocean vertical grid structure

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

  • mixing_energy :: [in] Energy values for up to 3 MLDs [R Z3 T-2 ~> J m-2]

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

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

  • diagptr :: Diagnostics structure

Called from:

mom_diabatic_driver::diabatic

subroutine mom_diabatic_aux/applyboundaryfluxesinout(CS, G, GV, US, dt, fluxes, optics, nsw, h, tv, aggregate_FW_forcing, evap_CFL_limit, minimum_forcing_depth, cTKE, dSV_dT, dSV_dS, SkinBuoyFlux, MLD)

Update the thickness, temperature, and salinity due to thermodynamic boundary forcing (contained in fluxes type) applied to h, tvT and tvS, and calculate the TKE implications of this heating.

Parameters:
  • cs :: Control structure for diabatic_aux

  • g :: [in] Grid structure

  • gv :: [in] ocean vertical grid structure

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

  • dt :: [in] Time-step over which forcing is applied [T ~> s]

  • fluxes :: [inout] Surface fluxes container

  • optics :: Optical properties container

  • nsw :: [in] The number of frequency bands of penetrating shortwave radiation

  • h :: [inout] Layer thickness [H ~> m or kg m-2]

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

  • aggregate_fw_forcing :: [in] If False, treat in/out fluxes separately.

  • evap_cfl_limit :: [in] The largest fraction of a layer that can be evaporated in one time-step [nondim].

  • minimum_forcing_depth :: [in] The smallest depth over which heat and freshwater fluxes is applied [H ~> m or kg m-2].

  • ctke :: [out] Turbulent kinetic energy requirement to mix

  • dsv_dt :: [out] Partial derivative of specific volume with

  • dsv_ds :: [out] Partial derivative of specific volume with

  • skinbuoyflux :: [out] Buoyancy flux at surface [Z2 T-3 ~> m2 s-3].

  • mld :: Mixed layer depth for brine plumes [Z ~> m]

Call to:

mom_opacity::absorbremainingsw mom_forcing_type::extractfluxes1d mom_forcing_type::forcing_singlepointprint

subroutine mom_diabatic_aux/diabatic_aux_init(Time, G, GV, US, param_file, diag, CS, useALEalgorithm, use_ePBL)

This subroutine initializes the parameters and control structure of the diabatic_aux module.

Parameters:
  • time :: [in] The current model time.

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

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

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

  • param_file :: [in] A structure to parse for run-time parameters

  • diag :: [inout] A structure used to regulate diagnostic output

  • cs :: A pointer to the control structure for the diabatic_aux module, which is initialized here.

  • usealealgorithm :: [in] If true, use the ALE algorithm rather than layered mode.

  • use_epbl :: [in] If true, use the implicit energetics planetary boundary layer scheme to determine the diffusivity in the surface boundary layer.

Call to:

id_clock_frazil id_clock_uv_at_h

subroutine mom_diabatic_aux/diabatic_aux_end(CS)

This subroutine initializes the control structure and any related memory for the diabatic_aux module.

Parameters:

cs :: The control structure returned by a previous call to diabatic_aux_init; it is deallocated here.