namespace mom_bulk_mixed_layer

Overview

Build mixed layer parameterization. More…

namespace mom_bulk_mixed_layer {

// global variables

integer id_clock_detrain =0;

// global functions

subroutine, public bulkmixedlayer(
    h_3d h_3d,
    u_3d u_3d,
    v_3d v_3d,
    tv tv,
    fluxes fluxes,
    dt_in_T dt_in_T,
    ea ea,
    eb eb,
    G G,
    GV GV,
    US US,
    CS CS,
    optics optics,
    Hml Hml,
    aggregate_FW_forcing aggregate_FW_forcing,
    dt_diag dt_diag,
    last_call last_call
    );

subroutine, public bulkmixedlayer_init(
    Time Time,
    G G,
    GV GV,
    US US,
    param_file param_file,
    diag diag,
    CS CS
    );

} // namespace mom_bulk_mixed_layer

Detailed Documentation

Build mixed layer parameterization.

By Robert Hallberg, 1997 - 2005.

This file contains the subroutine (bulkmixedlayer) that implements a Kraus-Turner-like bulk mixed layer, based on the work of various people, as described in the review paper by Niiler and Kraus (1979), with particular attention to the form proposed by Oberhuber (JPO, 1993, 808-829), with an extension to a refied bulk mixed layer as described in Hallberg (Aha Huliko’a, 2003). The physical processes portrayed in this subroutine include convective adjustment and mixed layer entrainment and detrainment. Penetrating shortwave radiation and an exponential decay of TKE fluxes are also supported by this subroutine. Several constants can alternately be set to give a traditional Kraus-Turner mixed layer scheme, although that is not the preferred option. The physical processes and arguments are described in detail below.

Global Variables

integer id_clock_detrain =0

CPU clock IDs.

Global Functions

subroutine, public bulkmixedlayer(
    h_3d h_3d,
    u_3d u_3d,
    v_3d v_3d,
    tv tv,
    fluxes fluxes,
    dt_in_T dt_in_T,
    ea ea,
    eb eb,
    G G,
    GV GV,
    US US,
    CS CS,
    optics optics,
    Hml Hml,
    aggregate_FW_forcing aggregate_FW_forcing,
    dt_diag dt_diag,
    last_call last_call
    )

This subroutine partially steps the bulk mixed layer model. The following processes are executed, in the order listed.

  1. Undergo convective adjustment into mixed layer.

  2. Apply surface heating and cooling.

  3. Starting from the top, entrain whatever fluid the TKE budget permits. Penetrating shortwave radiation is also applied at this point.

  4. If there is any unentrained fluid that was formerly in the mixed layer, detrain this fluid into the buffer layer. This is equivalent to the mixed layer detraining to the Monin- Obukhov depth.

  5. Divide the fluid in the mixed layer evenly into CSnkml pieces.

  6. Split the buffer layer if appropriate. Layers 1 to nkml are the mixed layer, nkml+1 to nkml+nkbl are the buffer layers. The results of this subroutine are mathematically identical if there are multiple pieces of the mixed layer with the same density or if there is just a single layer. There is no stability limit on the time step.

The key parameters for the mixed layer are found in the control structure. These include mstar, nstar, nstar2, pen_SW_frac, pen_SW_scale, and TKE_decay. For the Oberhuber (1993) mixed layer, the values of these are: pen_SW_frac = 0.42, pen_SW_scale = 15.0 m, mstar = 1.25, nstar = 1, TKE_decay = 2.5, conv_decay = 0.5 TKE_decay is 1/kappa in eq. 28 of Oberhuber (1993), while conv_decay is 1/mu. Conv_decay has been eliminated in favor of the well-calibrated form for the efficiency of penetrating convection from Wang (2003). For a traditional Kraus-Turner mixed layer, the values are: pen_SW_frac = 0.0, pen_SW_scale = 0.0 m, mstar = 1.25, nstar = 0.4, TKE_decay = 0.0, conv_decay = 0.0

Parameters:

g

The ocean’s grid structure.

gv

The ocean’s vertical grid structure.

us

A dimensional unit scaling type

h_3d

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

u_3d

Zonal velocities interpolated to h points

v_3d

Zonal velocities interpolated to h points

tv

A structure containing pointers to any available thermodynamic fields. Absent fields have NULL ptrs.

fluxes

A structure containing pointers to any possible forcing fields. Unused fields have NULL ptrs.

dt_in_t

Time increment [T ~> s].

ea

The amount of fluid moved downward into a

eb

The amount of fluid moved upward into a

cs

The control structure returned by a previous call to mixedlayer_init.

optics

The structure containing the inverse of the vertical absorption decay scale for penetrating shortwave radiation [m-1].

hml

Active mixed layer depth [m].

aggregate_fw_forcing

If true, the net incoming and outgoing surface freshwater fluxes are combined before being applied, instead of being applied separately.

dt_diag

The diagnostic time step, which may be less than dt if there are two callse to mixedlayer [T ~> s].

last_call

if true, this is the last call to mixedlayer in the current time step, so diagnostics will be written. The default is .true.

subroutine, public bulkmixedlayer_init(
    Time Time,
    G G,
    GV GV,
    US US,
    param_file param_file,
    diag diag,
    CS CS
    )

This subroutine initializes the MOM bulk mixed layer module.

Parameters:

time

The model’s clock with the current 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 that is used to regulate diagnostic output.

cs

A pointer that is set to point to the control structure for this module.