mom_lateral_boundary_diffusion module reference

Calculates and applies diffusive fluxes as a parameterization of lateral mixing (non-neutral) by mesoscale eddies near the top and bottom (to be implemented) boundary layers of the ocean.

More…

Data Types

lateral_boundary_diffusion_cs

Sets parameters for lateral boundary mixing module.

Functions/Subroutines

lateral_boundary_diffusion_init()

Initialization routine that reads runtime parameters and sets up pointers to other control structures that might be needed for lateral boundary diffusion.

lateral_boundary_diffusion()

Driver routine for calculating lateral diffusive fluxes near the top and bottom boundaries.

bulk_average()

harmonic_mean()

Calculate the harmonic mean of two quantities See Harmonic Mean .

boundary_k_range()

Find the k-index range corresponding to the layers that are within the boundary-layer region.

fluxes_layer_method()

Calculate the lateral boundary diffusive fluxes using the layer by layer method.

fluxes_bulk_method()

Apply the lateral boundary diffusive fluxes calculated from a ‘bulk model’ See Bulk layer approach (Method #2) .

near_boundary_unit_tests()

Unit tests for near-boundary horizontal mixing.

test_layer_fluxes()

Returns true if output of near-boundary unit tests does not match correct computed values and conditionally writes results to stream.

test_boundary_k_range()

Return true if output of unit tests for boundary_k_range does not match answers.

Detailed Description

The Lateral Boundary Diffusion (LBD) framework

The LBD framework accounts for the effects of diabatic mesoscale fluxes within surface and bottom boundary layers. Unlike the equivalent adiabatic fluxes, which is applied along neutral density surfaces, LBD is purely horizontal.

The bottom boundary layer fluxes remain to be implemented, although most of the steps needed to do so have already been added and tested.

Boundary lateral diffusion can be applied using one of the three methods:

A brief summary of these methods is provided below.

Along layer approach (Method #1)

This is the recommended and more straight forward method where diffusion is applied layer by layer using only information from neighboring cells.

Step #1: compute vertical indices containing boundary layer (boundary_k_range). For the TOP boundary layer, these are:

k_top, k_bot, zeta_top, zeta_bot

Step #2: calculate the diffusive flux at each layer:

\[F_{k} = -KHTR \times h_{eff}(k) \times (\phi_R(k) - \phi_L(k)),\]

where h_eff is the harmonic mean of the layer thickness in the left and right columns. This method does not require a limiter since KHTR is already limted based on a diffusive CFL condition prior to the call of this module.

Step #3: option to linearly decay the flux from k_bot_min to k_bot_max:

If LBD_LINEAR_TRANSITION = True and k_bot_diff > 1, the diffusive flux will decay linearly between the top interface of the layer containing the minimum boundary layer depth (k_bot_min) and the lower interface of the layer containing the maximum layer depth (k_bot_max).

Bulk layer approach (Method #2)

Apply the lateral boundary diffusive fluxes calculated from a ‘bulk model’.This is a lower order representation (Kraus-Turner like approach) which assumes that eddies are acting along well mixed layers (i.e., eddies do not know care about vertical tracer gradients within the boundary layer).

Step #1: compute vertical indices containing boundary layer (boundary_k_range). For the TOP boundary layer, these are:

k_top, k_bot, zeta_top, zeta_bot

Step #2: compute bulk averages (thickness weighted) tracer averages (phi_L and phi_R), then calculate the bulk diffusive flux (F_{bulk}):

\[F_{bulk} = -KHTR \times h_{eff} \times (\phi_R - \phi_L),\]

where h_eff is the harmonic mean of the boundary layer depth in the left and right columns (

\[HBL_L\]

and

\[HBL_R\]

, respectively).

Step #3: decompose F_bulk onto individual layers:

\[F_{layer}(k) = F_{bulk} \times h_{frac}(k) ,\]

where h_{frac} is

\[h_{frac}(k) = h_u(k) \times \frac{1}{\sum(h_u)}.\]

h_u is the harmonic mean of thicknesses at each layer. Special care (layer reconstruction) must be taken at k_min = min(k_botL, k_bot_R).

Step #4: option to linearly decay the flux from k_bot_min to k_bot_max:

If LBD_LINEAR_TRANSITION = True and k_bot_diff > 1, the diffusive flux will decay linearly between the top interface of the layer containing the minimum boundary layer depth (k_bot_min) and the lower interface of the layer containing the maximum layer depth (k_bot_max).

Step #5: limit the tracer flux so that 1) only down-gradient fluxes are applied, and 2) the flux cannot be larger than F_max, which is defined using the tracer gradient:

\[F_{max} = -0.2 \times [(V_R(k) \times \phi_R(k)) - (V_L(k) \times \phi_L(k))],\]

where V is the cell volume. Why 0.2? t=0 t=inf 0 .2 0 1 0 .2.2.2 0 .2

Harmonic Mean

The harmonic mean (HM) betwen h1 and h2 is defined as:

\[HM = \frac{2 \times h1 \times h2}{h1 + h2}\]

Type Documentation

type mom_lateral_boundary_diffusion/lateral_boundary_diffusion_cs

Sets parameters for lateral boundary mixing module.

Type fields
  • % method [integer] :: Determine which of the three methods calculate and apply near boundary layer fluxes.

  • % deg [integer] :: Degree of polynomial reconstruction.

  • % surface_boundary_scheme [integer] :: Which boundary layer scheme to use.

  • % limiter [logical] :: Controls wether a flux limiter is applied. Only valid when method = 2.

  • % linear [logical] :: If True, apply a linear transition at the base/top of the boundary. The flux will be fully applied at k=k_min and zero at k=k_max.

  • % remap_cs [type(remapping_cs)] :: Control structure to hold remapping configuration.

  • % kpp_csp [type(kpp_cs),pointer] :: KPP control structure needed to get BLD.

  • % energetic_pbl_csp [type(energetic_pbl_cs),pointer] :: ePBL control structure needed to get BLD

  • % diag [type(diag_ctrl),pointer] :: A structure that is used to regulate the timing of diagnostic output.

Function/Subroutine Documentation

function mom_lateral_boundary_diffusion/lateral_boundary_diffusion_init(Time, G, param_file, diag, diabatic_CSp, CS) [logical]

Initialization routine that reads runtime parameters and sets up pointers to other control structures that might be needed for lateral boundary diffusion.

Parameters
  • time :: [in] Time structure

  • g :: [in] Grid structure

  • param_file :: [in] Parameter file structure

  • diag :: [inout] Diagnostics control structure

  • diabatic_csp :: KPP control structure needed to get BLD

  • cs :: Lateral boundary mixing control structure

Call to

mom_diabatic_driver::extract_diabatic_member mdl mom_error_handler::mom_error mom_remapping::remappingdefaultscheme mom_remapping::remappingschemesdoc

subroutine mom_lateral_boundary_diffusion/lateral_boundary_diffusion(G, GV, US, h, Coef_x, Coef_y, dt, Reg, CS)

Driver routine for calculating lateral diffusive fluxes near the top and bottom boundaries. Two different methods are available: Method 1: lower order representation, calculate fluxes from bulk layer integrated quantities. Method 2: more straight forward, diffusion is applied layer by layer using only information from neighboring cells.

Parameters
  • g :: [inout] Grid type

  • gv :: [in] ocean vertical grid structure

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

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

  • coef_x :: [in] dt * Kh * dy / dx at u-points [L2 ~> m2]

  • coef_y :: [in] dt * Kh * dx / dy at v-points [L2 ~> m2]

  • dt :: [in] Tracer time step * I_numitts (I_numitts in tracer_hordiff)

  • reg :: Tracer registry

  • cs :: [in] Control structure for this module

Call to

mom_energetic_pbl::energetic_pbl_get_mld fluxes_bulk_method fluxes_layer_method mom_cvmix_kpp::kpp_get_bld surface

Called from

mom_tracer_hor_diff::tracer_hordiff

function mom_lateral_boundary_diffusion/bulk_average(boundary, nk, deg, h, hBLT, phi, ppoly0_E, ppoly0_coefs, method, k_top, zeta_top, k_bot, zeta_bot) [real]
Parameters
  • boundary :: SURFACE or BOTTOM [nondim]

  • nk :: Number of layers [nondim]

  • deg :: Degree of polynomial [nondim]

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

  • hblt :: Depth of the boundary layer [H ~> m or kg m-2]

  • phi :: Scalar quantity

  • ppoly0_e :: Edge value of polynomial

  • ppoly0_coefs :: Coefficients of polynomial

  • method :: Remapping scheme to use

  • k_top :: Index of the first layer within the boundary

  • zeta_top :: Fraction of the layer encompassed by the bottom boundary layer (0 if none, 1. if all). For the surface, this is always 0. because integration starts at the surface [nondim]

  • k_bot :: Index of the last layer within the boundary

  • zeta_bot :: Fraction of the layer encompassed by the surface boundary layer (0 if none, 1. if all). For the bottom boundary layer, this is always 1. because integration starts at the bottom [nondim]

Call to

mom_remapping::average_value_ppoly bottom mom_error_handler::mom_error surface

Called from

fluxes_bulk_method

function mom_lateral_boundary_diffusion/harmonic_mean(h1, h2) [real]

Calculate the harmonic mean of two quantities See Harmonic Mean .

Parameters
  • h1 :: Scalar quantity

  • h2 :: Scalar quantity

Called from

fluxes_bulk_method fluxes_layer_method

subroutine mom_lateral_boundary_diffusion/boundary_k_range(boundary, nk, h, hbl, k_top, zeta_top, k_bot, zeta_bot)

Find the k-index range corresponding to the layers that are within the boundary-layer region.

Parameters
  • boundary :: [in] SURFACE or BOTTOM [nondim]

  • nk :: [in] Number of layers [nondim]

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

  • hbl :: [in] Thickness of the boundary layer [H ~> m or kg m-2] If surface, with respect to zbl_ref = 0. If bottom, with respect to zbl_ref = SUM(h)

  • k_top :: [out] Index of the first layer within the boundary

  • zeta_top :: [out] Distance from the top of a layer to the intersection of the top extent of the boundary layer (0 at top, 1 at bottom) [nondim]

  • k_bot :: [out] Index of the last layer within the boundary

  • zeta_bot :: [out] Distance of the lower layer to the boundary layer depth (0 at top, 1 at bottom) [nondim]

Call to

bottom mom_error_handler::mom_error surface

Called from

fluxes_bulk_method fluxes_layer_method near_boundary_unit_tests mom_neutral_diffusion::neutral_diffusion_calc_coeffs

subroutine mom_lateral_boundary_diffusion/fluxes_layer_method(boundary, nk, deg, h_L, h_R, hbl_L, hbl_R, area_L, area_R, phi_L, phi_R, ppoly0_coefs_L, ppoly0_coefs_R, ppoly0_E_L, ppoly0_E_R, method, khtr_u, F_layer, linear_decay)

Calculate the lateral boundary diffusive fluxes using the layer by layer method. See Along layer approach (Method #1) .

Parameters
  • boundary :: [in] Which boundary layer SURFACE or BOTTOM [nondim]

  • nk :: [in] Number of layers [nondim]

  • deg :: [in] order of the polynomial reconstruction [nondim]

  • h_l :: [in] Layer thickness (left) [H ~> m or kg m-2]

  • h_r :: [in] Layer thickness (right) [H ~> m or kg m-2]

  • hbl_l :: [in] Thickness of the boundary boundary layer (left) [H ~> m or kg m-2]

  • hbl_r :: [in] Thickness of the boundary boundary layer (right) [H ~> m or kg m-2]

  • area_l :: [in] Area of the horizontal grid (left) [L2 ~> m2]

  • area_r :: [in] Area of the horizontal grid (right) [L2 ~> m2]

  • phi_l :: [in] Tracer values (left) [conc]

  • phi_r :: [in] Tracer values (right) [conc]

  • ppoly0_coefs_l :: [in] Tracer reconstruction (left) [conc]

  • ppoly0_coefs_r :: [in] Tracer reconstruction (right) [conc]

  • ppoly0_e_l :: [in] Polynomial edge values (left) [nondim]

  • ppoly0_e_r :: [in] Polynomial edge values (right) [nondim]

  • method :: [in] Method of polynomial integration [nondim]

  • khtr_u :: [in] Horizontal diffusivities times delta t at a velocity point [L2 ~> m2]

  • f_layer :: [out] Layerwise diffusive flux at U- or V-point [H L2 conc ~> m3 conc]

  • linear_decay :: [in] If True, apply a linear transition at the base of the boundary layer

Call to

mom_remapping::average_value_ppoly bottom boundary_k_range harmonic_mean surface

Called from

lateral_boundary_diffusion near_boundary_unit_tests

subroutine mom_lateral_boundary_diffusion/fluxes_bulk_method(boundary, nk, deg, h_L, h_R, hbl_L, hbl_R, area_L, area_R, phi_L, phi_R, ppoly0_coefs_L, ppoly0_coefs_R, ppoly0_E_L, ppoly0_E_R, method, khtr_u, F_bulk, F_layer, F_limit, linear_decay)

Apply the lateral boundary diffusive fluxes calculated from a ‘bulk model’ See Bulk layer approach (Method #2) .

Parameters
  • boundary :: [in] Which boundary layer SURFACE or BOTTOM [nondim]

  • nk :: [in] Number of layers [nondim]

  • deg :: [in] order of the polynomial reconstruction [nondim]

  • h_l :: [in] Layer thickness (left) [H ~> m or kg m-2]

  • h_r :: [in] Layer thickness (right) [H ~> m or kg m-2]

  • hbl_l :: [in] Thickness of the boundary boundary layer (left) [H ~> m or kg m-2]

  • hbl_r :: [in] Thickness of the boundary boundary layer (left) [H ~> m or kg m-2]

  • area_l :: [in] Area of the horizontal grid (left) [L2 ~> m2]

  • area_r :: [in] Area of the horizontal grid (right) [L2 ~> m2]

  • phi_l :: [in] Tracer values (left) [conc]

  • phi_r :: [in] Tracer values (right) [conc]

  • ppoly0_coefs_l :: [in] Tracer reconstruction (left) [conc]

  • ppoly0_coefs_r :: [in] Tracer reconstruction (right) [conc]

  • ppoly0_e_l :: [in] Polynomial edge values (left) [nondim]

  • ppoly0_e_r :: [in] Polynomial edge values (right) [nondim]

  • method :: [in] Method of polynomial integration [nondim]

  • khtr_u :: [in] Horizontal diffusivities times delta t at a velocity point [L2 ~> m2]

  • f_bulk :: [out] The bulk mixed layer lateral flux [H L2 conc ~> m3 conc]

  • f_layer :: [out] Layerwise diffusive flux at U- or V-point [H L2 conc ~> m3 conc]

  • f_limit :: [in] If True, apply a limiter

  • linear_decay :: [in] If True, apply a linear transition at the base of the boundary layer

Call to

bottom boundary_k_range bulk_average harmonic_mean surface

Called from

lateral_boundary_diffusion near_boundary_unit_tests

function mom_lateral_boundary_diffusion/near_boundary_unit_tests(verbose) [logical]

Unit tests for near-boundary horizontal mixing.

Parameters

verbose :: [in] If true, output additional information for debugging unit tests

Call to

bottom boundary_k_range fluxes_bulk_method fluxes_layer_method surface test_boundary_k_range test_layer_fluxes

Called from

mom_unit_tests::unit_tests

function mom_lateral_boundary_diffusion/test_layer_fluxes(verbose, nk, test_name, F_calc, F_ans) [logical]

Returns true if output of near-boundary unit tests does not match correct computed values and conditionally writes results to stream.

Parameters
  • verbose :: [in] If true, write results to stdout

  • test_name :: [in] Brief description of the unit test

  • nk :: [in] Number of layers

  • f_calc :: [in] Fluxes of the unitless tracer from the algorithm [s^-1]

  • f_ans :: [in] Fluxes of the unitless tracer calculated by hand [s^-1]

Called from

near_boundary_unit_tests

function mom_lateral_boundary_diffusion/test_boundary_k_range(k_top, zeta_top, k_bot, zeta_bot, k_top_ans, zeta_top_ans, k_bot_ans, zeta_bot_ans, test_name, verbose) [logical]

Return true if output of unit tests for boundary_k_range does not match answers.

Parameters
  • k_top :: Index of cell containing top of boundary

  • zeta_top :: Nondimension position

  • k_bot :: Index of cell containing bottom of boundary

  • zeta_bot :: Nondimension position

  • k_top_ans :: Index of cell containing top of boundary

  • zeta_top_ans :: Nondimension position

  • k_bot_ans :: Index of cell containing bottom of boundary

  • zeta_bot_ans :: Nondimension position

  • test_name :: Name of the unit test

  • verbose :: If true always print output

Called from

near_boundary_unit_tests