|
MOM6
|
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.
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.
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).
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
The harmonic mean (HM) betwen h1 and h2 is defined as:
\[ HM = \frac{2 \times h1 \times h2}{h1 + h2} \]
Data Types | |
| type | lateral_boundary_diffusion_cs |
| Sets parameters for lateral boundary mixing module. More... | |
Functions/Subroutines | |
| logical function, public | lateral_boundary_diffusion_init (Time, G, param_file, diag, diabatic_CSp, CS) |
| Initialization routine that reads runtime parameters and sets up pointers to other control structures that might be needed for lateral boundary diffusion. More... | |
| subroutine, public | 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. More... | |
| real function | bulk_average (boundary, nk, deg, h, hBLT, phi, ppoly0_E, ppoly0_coefs, method, k_top, zeta_top, k_bot, zeta_bot) |
| real function | harmonic_mean (h1, h2) |
| Calculate the harmonic mean of two quantities See Harmonic Mean. More... | |
| subroutine, public | 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. More... | |
| subroutine | 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). More... | |
| subroutine | 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). More... | |
| logical function, public | near_boundary_unit_tests (verbose) |
| Unit tests for near-boundary horizontal mixing. More... | |
| logical function | test_layer_fluxes (verbose, nk, test_name, F_calc, F_ans) |
| Returns true if output of near-boundary unit tests does not match correct computed values and conditionally writes results to stream. More... | |
| logical function | 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) |
| Return true if output of unit tests for boundary_k_range does not match answers. More... | |
Variables | |
| integer, parameter, public | surface = -1 |
| Set a value that corresponds to the surface bopundary. | |
| integer, parameter, public | bottom = 1 |
| Set a value that corresponds to the bottom boundary. | |
| character(len=40) | mdl = "MOM_lateral_boundary_diffusion" |
| Name of this module. | |
| subroutine, public mom_lateral_boundary_diffusion::boundary_k_range | ( | integer, intent(in) | boundary, |
| integer, intent(in) | nk, | ||
| real, dimension(nk), intent(in) | h, | ||
| real, intent(in) | hbl, | ||
| integer, intent(out) | k_top, | ||
| real, intent(out) | zeta_top, | ||
| integer, intent(out) | k_bot, | ||
| real, intent(out) | zeta_bot | ||
| ) |
Find the k-index range corresponding to the layers that are within the boundary-layer region.
| [in] | boundary | SURFACE or BOTTOM [nondim] |
| [in] | nk | Number of layers [nondim] |
| [in] | h | Layer thicknesses of the column [H ~> m or kg m-2] |
| [in] | hbl | 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) |
| [out] | k_top | Index of the first layer within the boundary |
| [out] | zeta_top | 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] |
| [out] | k_bot | Index of the last layer within the boundary |
| [out] | zeta_bot | Distance of the lower layer to the boundary layer depth (0 at top, 1 at bottom) [nondim] |
Definition at line 374 of file MOM_lateral_boundary_diffusion.F90.
|
private |
| 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] |
Definition at line 310 of file MOM_lateral_boundary_diffusion.F90.
|
private |
Apply the lateral boundary diffusive fluxes calculated from a 'bulk model' See Bulk layer approach (Method #2).
| [in] | boundary | Which boundary layer SURFACE or BOTTOM [nondim] |
| [in] | nk | Number of layers [nondim] |
| [in] | deg | order of the polynomial reconstruction [nondim] |
| [in] | h_l | Layer thickness (left) [H ~> m or kg m-2] |
| [in] | h_r | Layer thickness (right) [H ~> m or kg m-2] |
| [in] | hbl_l | Thickness of the boundary boundary layer (left) [H ~> m or kg m-2] |
| [in] | hbl_r | Thickness of the boundary boundary layer (left) [H ~> m or kg m-2] |
| [in] | area_l | Area of the horizontal grid (left) [L2 ~> m2] |
| [in] | area_r | Area of the horizontal grid (right) [L2 ~> m2] |
| [in] | phi_l | Tracer values (left) [conc] |
| [in] | phi_r | Tracer values (right) [conc] |
| [in] | ppoly0_coefs_l | Tracer reconstruction (left) [conc] |
| [in] | ppoly0_coefs_r | Tracer reconstruction (right) [conc] |
| [in] | ppoly0_e_l | Polynomial edge values (left) [nondim] |
| [in] | ppoly0_e_r | Polynomial edge values (right) [nondim] |
| [in] | method | Method of polynomial integration [nondim] |
| [in] | khtr_u | Horizontal diffusivities times delta t at a velocity point [L2 ~> m2] |
| [out] | f_bulk | The bulk mixed layer lateral flux [H L2 conc ~> m3 conc] |
| [out] | f_layer | Layerwise diffusive flux at U- or V-point [H L2 conc ~> m3 conc] |
| [in] | f_limit | If True, apply a limiter |
| [in] | linear_decay | If True, apply a linear transition at the base of the boundary layer |
Definition at line 595 of file MOM_lateral_boundary_diffusion.F90.
|
private |
Calculate the lateral boundary diffusive fluxes using the layer by layer method. See Along layer approach (Method #1).
| [in] | boundary | Which boundary layer SURFACE or BOTTOM [nondim] |
| [in] | nk | Number of layers [nondim] |
| [in] | deg | order of the polynomial reconstruction [nondim] |
| [in] | h_l | Layer thickness (left) [H ~> m or kg m-2] |
| [in] | h_r | Layer thickness (right) [H ~> m or kg m-2] |
| [in] | hbl_l | Thickness of the boundary boundary layer (left) [H ~> m or kg m-2] |
| [in] | hbl_r | Thickness of the boundary boundary layer (right) [H ~> m or kg m-2] |
| [in] | area_l | Area of the horizontal grid (left) [L2 ~> m2] |
| [in] | area_r | Area of the horizontal grid (right) [L2 ~> m2] |
| [in] | phi_l | Tracer values (left) [conc] |
| [in] | phi_r | Tracer values (right) [conc] |
| [in] | ppoly0_coefs_l | Tracer reconstruction (left) [conc] |
| [in] | ppoly0_coefs_r | Tracer reconstruction (right) [conc] |
| [in] | ppoly0_e_l | Polynomial edge values (left) [nondim] |
| [in] | ppoly0_e_r | Polynomial edge values (right) [nondim] |
| [in] | method | Method of polynomial integration [nondim] |
| [in] | khtr_u | Horizontal diffusivities times delta t at a velocity point [L2 ~> m2] |
| [out] | f_layer | Layerwise diffusive flux at U- or V-point [H L2 conc ~> m3 conc] |
| [in] | linear_decay | If True, apply a linear transition at the base of the boundary layer |
Definition at line 441 of file MOM_lateral_boundary_diffusion.F90.
|
private |
Calculate the harmonic mean of two quantities See Harmonic Mean.
| h1 | Scalar quantity |
| h2 | Scalar quantity |
Definition at line 363 of file MOM_lateral_boundary_diffusion.F90.
| subroutine, public mom_lateral_boundary_diffusion::lateral_boundary_diffusion | ( | type(ocean_grid_type), intent(inout) | G, |
| type(verticalgrid_type), intent(in) | GV, | ||
| type(unit_scale_type), intent(in) | US, | ||
| real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | h, | ||
| real, dimension( g %isdb: g %iedb, g %jsd: g %jed), intent(in) | Coef_x, | ||
| real, dimension( g %isd: g %ied, g %jsdb: g %jedb), intent(in) | Coef_y, | ||
| real, intent(in) | dt, | ||
| type(tracer_registry_type), pointer | Reg, | ||
| type(lateral_boundary_diffusion_cs), intent(in) | 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.
| [in,out] | g | Grid type |
| [in] | gv | ocean vertical grid structure |
| [in] | us | A dimensional unit scaling type |
| [in] | h | Layer thickness [H ~> m or kg m-2] |
| [in] | coef_x | dt * Kh * dy / dx at u-points [L2 ~> m2] |
| [in] | coef_y | dt * Kh * dx / dy at v-points [L2 ~> m2] |
| [in] | dt | Tracer time step * I_numitts (I_numitts in tracer_hordiff) |
| reg | Tracer registry | |
| [in] | cs | Control structure for this module |
Definition at line 137 of file MOM_lateral_boundary_diffusion.F90.
| logical function, public mom_lateral_boundary_diffusion::lateral_boundary_diffusion_init | ( | type(time_type), intent(in), target | Time, |
| type(ocean_grid_type), intent(in) | G, | ||
| type(param_file_type), intent(in) | param_file, | ||
| type(diag_ctrl), intent(inout), target | diag, | ||
| type(diabatic_cs), pointer | diabatic_CSp, | ||
| type(lateral_boundary_diffusion_cs), pointer | CS | ||
| ) |
Initialization routine that reads runtime parameters and sets up pointers to other control structures that might be needed for lateral boundary diffusion.
| [in] | time | Time structure |
| [in] | g | Grid structure |
| [in] | param_file | Parameter file structure |
| [in,out] | diag | Diagnostics control structure |
| diabatic_csp | KPP control structure needed to get BLD | |
| cs | Lateral boundary mixing control structure |
Definition at line 67 of file MOM_lateral_boundary_diffusion.F90.
| logical function, public mom_lateral_boundary_diffusion::near_boundary_unit_tests | ( | logical, intent(in) | verbose | ) |
Unit tests for near-boundary horizontal mixing.
| [in] | verbose | If true, output additional information for debugging unit tests |
Definition at line 806 of file MOM_lateral_boundary_diffusion.F90.
|
private |
Return true if output of unit tests for boundary_k_range does not match answers.
| 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 |
Definition at line 1131 of file MOM_lateral_boundary_diffusion.F90.
|
private |
Returns true if output of near-boundary unit tests does not match correct computed values and conditionally writes results to stream.
| [in] | verbose | If true, write results to stdout |
| [in] | test_name | Brief description of the unit test |
| [in] | nk | Number of layers |
| [in] | f_calc | Fluxes of the unitless tracer from the algorithm [s^-1] |
| [in] | f_ans | Fluxes of the unitless tracer calculated by hand [s^-1] |
Definition at line 1106 of file MOM_lateral_boundary_diffusion.F90.