mom_diabatic_aux module reference¶
Provides functions for some diabatic processes such as fraxil, brine rejection, tendency due to surface flux divergence.
Data Types¶
Control structure for diabatic_aux. |
Functions/Subroutines¶
Frazil formation keeps the temperature above the freezing point. |
|
This subroutine applies double diffusion to T & S, assuming no diapycal mass fluxes, using a simple triadiagonal solver. |
|
This subroutine keeps salinity from falling below a small but positive threshold. |
|
This is a simple tri-diagonal solver for T and S. |
|
This subroutine calculates u_h and v_h (velocities at thickness points), optionally using the entrainment amounts passed in as arguments. |
|
Diagnose a mixed layer depth (MLD) determined by a given density difference with the surface. |
|
Diagnose a mixed layer depth (MLD) determined by the depth a given energy value would mix. |
|
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. |
|
This subroutine initializes the parameters and control structure of the diabatic_aux module. |
|
This subroutine initializes the control structure and any related memory for the diabatic_aux module. |
Detailed Description¶
This module contains the subroutines that, along with the subroutines that it calls, implements diapycnal mass and momentum fluxes and a bulk mixed layer. The diapycnal diffusion can be used without the bulk mixed layer.
diabatic first determines the (diffusive) diapycnal mass fluxes based on the convergence of the buoyancy fluxes within each layer. The dual-stream entrainment scheme of MacDougall and Dewar (JPO, 1997) is used for combined diapycnal advection and diffusion, calculated implicitly and potentially with the Richardson number dependent mixing, as described by Hallberg (MWR, 2000). Diapycnal advection is fundamentally the residual of diapycnal diffusion, so the fully implicit upwind differencing scheme that is used is entirely appropriate. The downward buoyancy flux in each layer is determined from an implicit calculation based on the previously calculated flux of the layer above and an estimated flux in the layer below. This flux is subject to the following conditions: (1) the flux in the top and bottom layers are set by the boundary conditions, and (2) no layer may be driven below an Angstrom thick- ness. If there is a bulk mixed layer, the buffer layer is treat- ed as a fixed density layer with vanishingly small diffusivity.
diabatic takes 5 arguments: the two velocities (u and v), the thicknesses (h), a structure containing the forcing fields, and the length of time over which to act (dt). The velocities and thickness are taken as inputs and modified within the subroutine. There is no limit on the time step.
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].%
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 false.%
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
[integer] :: An integer handle used in time interpolation of chlorophyll read from a file.%
chl_from_file
[logical] :: If true, chl_a is read from a file.%
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, tv, visc, dt, G, GV)¶ This subroutine applies double diffusion to T & S, assuming no diapycal mass fluxes, using a simple triadiagonal 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]
tv :: [inout] Structure containing pointers to any available thermodynamic fields.
visc :: [in] Structure containing vertical viscosities, bottom boundary layer properies, and related fields.
dt :: [in] Time increment [T ~> s].
-
subroutine
mom_diabatic_aux/
adjust_salt
(h, tv, G, GV, CS, halo)¶ 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.
halo :: [in] Halo width over which to work
-
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 [degC].
s :: [inout] Layer salinities [ppt].
-
subroutine
mom_diabatic_aux/
find_uv_at_h
(u, v, h, u_h, v_h, G, GV, US, ea, eb)¶ 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
- Call to
id_clock_uv_at_h
-
subroutine
mom_diabatic_aux/
set_pen_shortwave
(optics, fluxes, G, GV, US, CS, opacity_CSp, tracer_flow_CSp)¶ - 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 ptrs
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_csp :: The control structure for the opacity module.
tracer_flow_csp :: A pointer to the control structure organizing the tracer modules.
- Call to
- Called from
mom_diabatic_driver::diabatic_ale
mom_diabatic_driver::diabatic_ale_legacy
mom_diabatic_driver::layered_diabatic
-
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_driver 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]
-
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_driver due to its position within the time stepping.
- Parameters
id_mld :: [in] Energy output diag 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
h :: [in] Layer thickness [H ~> m or kg m-2]
tv :: [in] Structure containing pointers to any available thermodynamic fields.
diagptr :: Diagnostics structure
-
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)¶ 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].
- 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.