namespace mom_diabatic_aux

Overview

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

namespace mom_diabatic_aux {

// global variables

integer id_clock_uv_at_h;

// global functions

subroutine, public make_frazil(h h, tv tv, G G, GV GV, CS CS, p_surf p_surf, halo halo);
subroutine, public differential_diffuse_t_s(h h, tv tv, visc visc, dt dt, G G, GV GV);
subroutine, public adjust_salt(h h, tv tv, G G, GV GV, CS CS, halo halo);

subroutine, public insert_brine(
    h h,
    tv tv,
    G G,
    GV GV,
    fluxes fluxes,
    nkmb nkmb,
    CS CS,
    dt dt,
    id_brine_lay id_brine_lay
    );

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

subroutine, public find_uv_at_h(
    u u,
    v v,
    h h,
    u_h u_h,
    v_h v_h,
    G G,
    GV GV,
    US US,
    ea ea,
    eb eb
    );

subroutine, public set_pen_shortwave(
    optics optics,
    fluxes fluxes,
    G G,
    GV GV,
    CS CS,
    opacity_CSp opacity_CSp,
    tracer_flow_CSp tracer_flow_CSp
    );

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

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

subroutine, public diabatic_aux_init(
    Time Time,
    G G,
    GV GV,
    US US,
    param_file param_file,
    diag diag,
    CS CS,
    useALEalgorithm useALEalgorithm,
    use_ePBL use_ePBL
    );

subroutine, public diabatic_aux_end(CS CS);

} // namespace mom_diabatic_aux

Detailed Documentation

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

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.

Global Variables

integer id_clock_uv_at_h

CPU time clock IDs.

Global Functions

subroutine, public make_frazil(
    h h,
    tv tv,
    G G,
    GV GV,
    CS CS,
    p_surf p_surf,
    halo 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 J m-2) in tvfrazil.

Parameters:

g

The ocean’s grid structure

gv

The ocean’s vertical grid structure

h

Layer thicknesses [H ~> m or kg m-2]

tv

Structure containing pointers to any available thermodynamic fields.

cs

The control structure returned by a previous call to diabatic_aux_init.

p_surf

The pressure at the ocean surface [Pa].

halo

Halo width over which to calculate frazil

subroutine, public differential_diffuse_t_s(
    h h,
    tv tv,
    visc visc,
    dt dt,
    G G,
    GV GV
    )

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

Parameters:

g

The ocean’s grid structure

gv

The ocean’s vertical grid structure

h

Layer thicknesses [H ~> m or kg m-2]

tv

Structure containing pointers to any available thermodynamic fields.

visc

Structure containing vertical viscosities, bottom boundary layer properies, and related fields.

dt

Time increment [T ~> s].

subroutine, public adjust_salt(h h, tv tv, G G, GV GV, CS CS, halo 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

The ocean’s grid structure

gv

The ocean’s vertical grid structure

h

Layer thicknesses [H ~> m or kg m-2]

tv

Structure containing pointers to any available thermodynamic fields.

cs

The control structure returned by a previous call to diabatic_aux_init.

halo

Halo width over which to work

subroutine, public insert_brine(
    h h,
    tv tv,
    G G,
    GV GV,
    fluxes fluxes,
    nkmb nkmb,
    CS CS,
    dt dt,
    id_brine_lay id_brine_lay
    )

Insert salt from brine rejection into the first layer below the mixed layer which both contains mass and in which the change in layer density remains stable after the addition of salt via brine rejection.

Parameters:

g

The ocean’s grid structure

gv

The ocean’s vertical grid structure

h

Layer thicknesses [H ~> m or kg m-2]

tv

Structure containing pointers to any available thermodynamic fields

fluxes

A structure of thermodynamic surface fluxes

nkmb

The number of layers in the mixed and buffer layers

cs

The control structure returned by a previous call to diabatic_aux_init

dt

The thermodynamic time step [s].

id_brine_lay

The handle for a diagnostic of which layer receivees the brine.

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

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

Parameters:

g

The ocean’s grid structure

gv

The ocean’s vertical grid structure

is

The start i-index to work on.

ie

The end i-index to work on.

js

The start j-index to work on.

je

The end j-index to work on.

hold

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

ea

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

eb

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

t

Layer potential temperatures [degC].

s

Layer salinities [ppt].

subroutine, public find_uv_at_h(
    u u,
    v v,
    h h,
    u_h u_h,
    v_h v_h,
    G G,
    GV GV,
    US US,
    ea ea,
    eb eb
    )

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

Parameters:

g

The ocean’s grid structure

gv

The ocean’s vertical grid structure

us

A dimensional unit scaling type

u

The zonal velocity [m s-1]

v

The meridional velocity [m s-1]

h

Layer thicknesses [H ~> m or kg m-2]

u_h

Zonal velocity interpolated to h points [m s-1].

v_h

Meridional velocity interpolated to h points [m s-1].

ea

The amount of fluid entrained from the layer

eb

The amount of fluid entrained from the layer

subroutine, public set_pen_shortwave(
    optics optics,
    fluxes fluxes,
    G G,
    GV GV,
    CS CS,
    opacity_CSp opacity_CSp,
    tracer_flow_CSp tracer_flow_CSp
    )

Parameters:

optics

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

fluxes

points to forcing fields unused fields have NULL ptrs

g

The ocean’s grid structure.

gv

The ocean’s vertical grid structure.

cs

Control structure for diabatic_aux

opacity_csp

The control structure for the opacity module.

tracer_flow_csp

A pointer to the control structure of the tracer modules.

subroutine, public diagnosemldbydensitydifference(
    id_MLD id_MLD,
    h h,
    tv tv,
    densityDiff densityDiff,
    G G,
    GV GV,
    US US,
    diagPtr diagPtr,
    id_N2subML id_N2subML,
    id_MLDsq id_MLDsq,
    dz_subML 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

Grid type

gv

ocean vertical grid structure

us

A dimensional unit scaling type

id_mld

Handle (ID) of MLD diagnostic

h

Layer thickness [H ~> m or kg m-2]

tv

Structure containing pointers to any available thermodynamic fields.

densitydiff

Density difference to determine MLD [kg m-3]

diagptr

Diagnostics structure

id_n2subml

Optional handle (ID) of subML stratification

id_mldsq

Optional handle (ID) of squared MLD

dz_subml

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

subroutine, public applyboundaryfluxesinout(
    CS CS,
    G G,
    GV GV,
    US US,
    dt dt,
    fluxes fluxes,
    optics optics,
    nsw nsw,
    h h,
    tv tv,
    aggregate_FW_forcing aggregate_FW_forcing,
    evap_CFL_limit evap_CFL_limit,
    minimum_forcing_depth minimum_forcing_depth,
    cTKE cTKE,
    dSV_dT dSV_dT,
    dSV_dS dSV_dS,
    SkinBuoyFlux 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

Grid structure

gv

ocean vertical grid structure

us

A dimensional unit scaling type

dt

Time-step over which forcing is applied [s]

fluxes

Surface fluxes container

optics

Optical properties container

nsw

The number of frequency bands of penetrating shortwave radiation

h

Layer thickness [H ~> m or kg m-2]

tv

Structure containing pointers to any available thermodynamic fields.

aggregate_fw_forcing

If False, treat in/out fluxes separately.

evap_cfl_limit

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

minimum_forcing_depth

The smallest depth over which heat and freshwater fluxes is applied [m].

ctke

Turbulent kinetic energy requirement to mix

dsv_dt

Partial derivative of specific volume with

dsv_ds

Partial derivative of specific volume with

skinbuoyflux

Buoyancy flux at surface [Z2 T-3 ~> m2 s-3].

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

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

Parameters:

time

The current model time.

g

The ocean’s grid structure

gv

The ocean’s vertical grid structure

us

A dimensional unit scaling type

param_file

A structure to parse for run-time parameters

diag

A structure used to regulate diagnostic output

cs

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

usealealgorithm

If true, use the ALE algorithm rather than layered mode.

use_epbl

If true, use the implicit energetics planetary boundary layer scheme to determine the diffusivity in the surface boundary layer.

subroutine, public diabatic_aux_end(CS 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.