mom_hor_bnd_diffusion module reference

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

More…

Data Types

hbd_cs

Sets parameters for horizontal boundary mixing module.

Functions/Subroutines

hor_bnd_diffusion_init()

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

hor_bnd_diffusion()

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

hbd_grid()

Build the HBD grid where tracers will be remapped to.

harmonic_mean()

Calculate the harmonic mean of two quantities See Harmonic Mean.

find_minimum()

Returns the location of the minimum value in a 1D array between indices s and e.

swap()

Swaps the values of its two formal arguments.

sort()

Receives a 1D array x and sorts it into ascending order.

unique()

Returns the unique values in a 1D array.

merge_interfaces()

Given layer thicknesses (and corresponding interfaces) and BLDs in two adjacent columns, return a set of 1-d layer thicknesses whose interfaces cover all interfaces in the left and right columns plus the two BLDs.

flux_limiter()

Calculates the maximum flux that can leave a cell and uses that to apply a limiter to F_layer.

boundary_k_range()

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

fluxes_layer_method()

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

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.

hbd_grid_test()

Same as hbd_grid, but only used in the unit tests.

hor_bnd_diffusion_end()

Deallocates hor_bnd_diffusion control structure.

Detailed Description

The Horizontal Boundary Diffusion (HBD) framework

The HBD 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, HBD is purely horizontal. To assure that diffusive fluxes are strictly horizontal regardless of the vertical coordinate system, this method relies on regridding/remapping techniques.

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

Horizontal boundary diffusion is applied as follows:

  1. remap tracer to a z* grid (HBD grid) 2) calculate diffusive tracer fluxes (F) in the HBD grid using a layer by layer approach (Along layer approach) 3) remap fluxes to the native grid 4) update tracer by adding the divergence of F

Along layer approach

Here diffusion is applied layer by layer using only information from neighboring cells.

Step #1: define vertical grid using interfaces and surface boundary layers from left and right columns (see merge_interfaces).

Step #2: 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.

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

If HBD_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 #4: remap the fluxes back to the native grid. This is done at velocity points, whose vertical grid is determined using harmonic mean. To assure monotonicity, tracer fluxes are limited 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_hor_bnd_diffusion/hbd_cs

Sets parameters for horizontal boundary mixing module.

Type fields:
  • % debug [logical] :: If true, write verbose checksums for debugging.

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

  • % hbd_nk [integer] :: Maximum number of levels in the HBD grid [nondim].

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

  • % limiter [logical] :: Controls whether a flux limiter is applied in the native grid (default is true).

  • % limiter_remap [logical] :: Controls whether a flux limiter is applied in the remapped grid (default is false).

  • % 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.

  • % h_subroundoff [real] :: A thickness that is so small that it can be added to a thickness of Angstrom or larger without changing it at the bit level [H ~> m or kg m-2]. If Angstrom is 0 or exceedingly small, this is negligible compared to 1e-17 m.

  • % hbd_grd_u [real(:,:,:),allocatable] :: HBD thicknesses at t-points adjacent to u-points [H ~> m or kg m-2].

  • % hbd_grd_v [real(:,:,:),allocatable] :: HBD thicknesses at t-points adjacent to v-points (left and right) [H ~> m or kg m-2].

  • % hbd_u_kmax [integer(:,:),allocatable] :: Maximum vertical index in hbd_grd_u [nondim].

  • % hbd_v_kmax [integer(:,:),allocatable] :: Maximum vertical index in hbd_grd_v [nondim].

  • % 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_hor_bnd_diffusion/hor_bnd_diffusion_init(Time, G, GV, US, 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 horizontal boundary diffusion.

Parameters:
  • time :: [in] Time structure

  • g :: [in] Grid structure

  • gv :: [in] ocean vertical grid structure

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

  • param_file :: [in] Parameter file structure

  • diag :: [inout] Diagnostics control structure

  • diabatic_csp :: KPP control structure needed to get BLD

  • cs :: Horizontal boundary mixing control structure

Call to:

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

subroutine mom_hor_bnd_diffusion/hor_bnd_diffusion(G, GV, US, h, Coef_x, Coef_y, dt, Reg, CS)

Driver routine for calculating horizontal diffusive fluxes near the top and bottom boundaries. Diffusion is applied using only information from neighboring cells, as follows: 1) remap tracer to a z* grid (HBD grid) 2) calculate diffusive tracer fluxes (F) in the HBD grid using a layer by layer approach 3) remap fluxes to the native grid 4) update tracer by adding the divergence of F.

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) [T ~> s]

  • reg :: Tracer registry

  • cs :: Control structure for this module

Call to:

mom_energetic_pbl::energetic_pbl_get_mld fluxes_layer_method mom_spatial_means::global_mass_integral hbd_grid id_clock_hbd mom_cvmix_kpp::kpp_get_bld mom_error_handler::mom_mesg surface

Called from:

mom_tracer_hor_diff::tracer_hordiff

subroutine mom_hor_bnd_diffusion/hbd_grid(boundary, G, GV, hbl, h, CS)

Build the HBD grid where tracers will be remapped to.

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

  • g :: [inout] Grid type

  • gv :: [in] ocean vertical grid structure

  • hbl :: Boundary layer depth [H ~> m or kg m-2]

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

  • cs :: Horizontal diffusion control structure

Call to:

merge_interfaces mom_error_handler::mom_error

Called from:

hor_bnd_diffusion

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

Calculate the harmonic mean of two quantities See Harmonic Mean.

Parameters:
  • h1 :: Scalar quantity [arbitrary]

  • h2 :: Scalar quantity [arbitrary]

Called from:

fluxes_layer_method

function mom_hor_bnd_diffusion/find_minimum(x, s, e) [integer]

Returns the location of the minimum value in a 1D array between indices s and e.

Parameters:
  • s :: [in] start index

  • e :: [in] end index

  • x :: [in] 1D array to be checked [arbitrary]

Called from:

sort

subroutine mom_hor_bnd_diffusion/swap(a, b)

Swaps the values of its two formal arguments.

Parameters:
  • a :: [inout] First value to be swapped [arbitrary]

  • b :: [inout] Second value to be swapped [arbitrary]

Called from:

sort

subroutine mom_hor_bnd_diffusion/sort(x, n)

Receives a 1D array x and sorts it into ascending order.

Parameters:
  • n :: [in] Number of points in the array

  • x :: [inout] 1D array to be sorted [arbitrary]

Call to:

find_minimum swap

Called from:

merge_interfaces near_boundary_unit_tests

subroutine mom_hor_bnd_diffusion/unique(val, n, val_unique, val_max)

Returns the unique values in a 1D array.

Parameters:
  • n :: [in] Number of points in the array.

  • val :: [in] 1D array to be checked [arbitrary]

  • val_unique :: [inout] Returned 1D array with unique values [arbitrary]

  • val_max :: [in] sets the maximum value in val_unique to this value [arbitrary]

Call to:

mom_error_handler::mom_error

Called from:

merge_interfaces near_boundary_unit_tests

subroutine mom_hor_bnd_diffusion/merge_interfaces(nk, h_L, h_R, hbl_L, hbl_R, H_subroundoff, h)

Given layer thicknesses (and corresponding interfaces) and BLDs in two adjacent columns, return a set of 1-d layer thicknesses whose interfaces cover all interfaces in the left and right columns plus the two BLDs. This can be used to accurately remap tracer tendencies in both columns.

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

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

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

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

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

  • h_subroundoff :: [in] GVH_subroundoff [H ~> m or kg m-2]

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

Call to:

sort unique

Called from:

hbd_grid hbd_grid_test near_boundary_unit_tests

subroutine mom_hor_bnd_diffusion/flux_limiter(F_layer, area_L, area_R, phi_L, phi_R, h_L, h_R)

Calculates the maximum flux that can leave a cell and uses that to apply a limiter to F_layer.

Parameters:
  • f_layer :: [inout] Tracer flux to be checked [H L2 conc ~> m3 conc]

  • area_l :: [in] Area of left cell [L2 ~> m2]

  • area_r :: [in] Area of right cell [L2 ~> m2]

  • h_l :: [in] Thickness of left cell [H ~> m or kg m-2]

  • h_r :: [in] Thickness of right cell [H ~> m or kg m-2]

  • phi_l :: [in] Tracer concentration in the left cell [conc]

  • phi_r :: [in] Tracer concentration in the right cell [conc]

Called from:

fluxes_layer_method

subroutine mom_hor_bnd_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:

mom_neutral_diffusion::compute_tapering_coeffs fluxes_layer_method near_boundary_unit_tests mom_neutral_diffusion::neutral_diffusion_calc_coeffs

subroutine mom_hor_bnd_diffusion/fluxes_layer_method(boundary, ke, hbl_L, hbl_R, h_L, h_R, phi_L, phi_R, khtr_u, F_layer, area_L, area_R, nk, dz_top, CS)

Calculate the horizontal boundary diffusive fluxes using the layer by layer method. See Along layer approach.

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

  • ke :: [in] Number of layers in the native grid [nondim]

  • 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]

  • h_l :: [in] Thicknesses in the native grid (left) [H ~> m or kg m-2]

  • h_r :: [in] Thicknesses in the native grid (right) [H ~> m or kg m-2]

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

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

  • khtr_u :: [in] Horizontal diffusivities times the time step at a velocity point and vertical interfaces [L2 ~> m2]

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

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

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

  • nk :: [in] Number of layers in the HBD grid [nondim]

  • dz_top :: [in] The HBD z grid [H ~> m or kg m-2]

  • cs :: Horizontal diffusion control structure

Call to:

boundary_k_range flux_limiter harmonic_mean surface

Called from:

hor_bnd_diffusion near_boundary_unit_tests

function mom_hor_bnd_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_layer_method hbd_grid_test merge_interfaces sort mom_io::stdout surface test_boundary_k_range test_layer_fluxes unique

Called from:

mom_unit_tests::unit_tests

function mom_hor_bnd_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 or other quantity from the algorithm [arbitrary]

  • f_ans :: [in] Expected value calculated by hand [arbitrary]

Call to:

mom_io::stdout

Called from:

near_boundary_unit_tests

function mom_hor_bnd_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 :: Fractional position in the cell of the top boundary [nondim]

  • k_bot :: Index of cell containing bottom of boundary

  • zeta_bot :: Fractional position in the cell of the bottom boundary [nondim]

  • k_top_ans :: Expected index of cell containing top of boundary

  • zeta_top_ans :: Expected fractional position of the top boundary [nondim]

  • k_bot_ans :: Expected index of cell containing bottom of boundary

  • zeta_bot_ans :: Expected fractional position of the bottom boundary [nondim]

  • test_name :: Name of the unit test

  • verbose :: If true always print output

Call to:

mom_io::stdout

Called from:

near_boundary_unit_tests

subroutine mom_hor_bnd_diffusion/hbd_grid_test(boundary, hbl_L, hbl_R, h_L, h_R, CS)

Same as hbd_grid, but only used in the unit tests.

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

  • hbl_l :: [in] Boundary layer depth, left [H ~> m or kg m-2]

  • hbl_r :: [in] Boundary layer depth, right [H ~> m or kg m-2]

  • h_l :: [in] Layer thickness in the native grid, left [H ~> m or kg m-2]

  • h_r :: [in] Layer thickness in the native grid, right [H ~> m or kg m-2]

  • cs :: Horizontal diffusion control structure

Call to:

merge_interfaces mom_error_handler::mom_error

Called from:

near_boundary_unit_tests

subroutine mom_hor_bnd_diffusion/hor_bnd_diffusion_end(CS)

Deallocates hor_bnd_diffusion control structure.

Parameters:

cs :: Horizontal boundary diffusion control structure

Called from:

mom_tracer_hor_diff::tracer_hor_diff_end