mom_zanna_bolton module reference

Calculates Zanna and Bolton 2020 parameterization Implemented by Perezhogin P.A. Contact: pperezhogin@gmail.com .

More…

Data Types

zb2020_cs

Control structure for Zanna-Bolton-2020 parameterization.

Functions/Subroutines

zb2020_init()

Read parameters, allocate and precompute arrays, register diagnosicts used in Zanna_Bolton_2020().

zb2020_end()

Deallocate any variables allocated in ZB_2020_init.

zb2020_copy_gradient_and_thickness()

Save precomputed velocity gradients and thickness from the horizontal eddy viscosity module We save as much halo for velocity gradients as possible In symmetric (preferable) memory model: halo 2 for sh_xx and halo 1 for sh_xy and vort_xy We apply zero boundary conditions to velocity gradients which is required for filtering operations.

zb2020_lateral_stress()

Baroclinic Zanna-Bolton-2020 parameterization, see eq.

compute_c_diss()

Compute the attenuation parameter similarly to Klower2018, Juricke2019,2020: c_diss = 1/(1+(shear/(f*R_diss))) where shear = sqrt(sh_xx**2 + sh_xy**2) or shear = sqrt(sh_xx**2 + sh_xy**2 + vort_xy**2) In symmetric memory model, components of velocity gradient tensor should have halo 1 and zero boundary conditions.

compute_stress()

Compute stress tensor T = (Txx, Txy; Txy, Tyy) Which consists of the deviatoric and trace components, respectively: T = (-vort_xy * sh_xy, vort_xy * sh_xx; vort_xy * sh_xx, vort_xy * sh_xy) + 1/2 * (vort_xy^2 + sh_xy^2 + sh_xx^2, 0; 0, vort_xy^2 + sh_xy^2 + sh_xx^2) This stress tensor is multiplied by precomputed kappa=-CSamplitude * Garea: T -> T * kappa The sign of the stress tensor is such that (neglecting h): (du/dt, dv/dt) = div(T) In symmetric memory model: sh_xy and vort_xy should have halo 1 and zero B.C.; sh_xx should have halo 2 and zero B.C.

compute_stress_divergence()

Compute the divergence of subgrid stress weighted with thickness, i.e.

filter_velocity_gradients()

Filtering of the velocity gradients sh_xx, sh_xy, vort_xy.

filter_stress()

Filtering of the stress tensor Txx, Tyy, Txy.

filter_hq()

Wrapper for filter_3D function.

filter_3d()

Spatial lateral filter applied to 3D array.

compute_energy_source()

Computes the 3D energy source term for the ZB2020 scheme similarly to MOM_diagnostics.F90, specifically 1125 line.

Detailed Description

Calculates Zanna and Bolton 2020 parameterization Implemented by Perezhogin P.A. Contact: pperezhogin@gmail.com .

Type Documentation

type mom_zanna_bolton/zb2020_cs

Control structure for Zanna-Bolton-2020 parameterization.

Type fields:
  • % id_zb2020u [integer] :: Diagnostic handles.

  • % id_zb2020v [integer] :: Diagnostic handles.

  • % id_ke_zb2020 [integer] :: Diagnostic handles.

  • % id_txx [integer] :: Diagnostic handles.

  • % id_tyy [integer] :: Diagnostic handles.

  • % id_txy [integer] :: Diagnostic handles.

  • % id_cdiss [integer] :: Diagnostic handles.

  • % id_clock_module [integer] :: CPU time clock IDs.

  • % id_clock_copy [integer] :: CPU time clock IDs.

  • % id_clock_cdiss [integer] :: CPU time clock IDs.

  • % id_clock_stress [integer] :: CPU time clock IDs.

  • % id_clock_divergence [integer] :: CPU time clock IDs.

  • % id_clock_mpi [integer] :: CPU time clock IDs.

  • % id_clock_filter [integer] :: CPU time clock IDs.

  • % id_clock_post [integer] :: CPU time clock IDs.

  • % id_clock_source [integer] :: CPU time clock IDs.

  • % pass_tq [type(group_pass_type)] :: MPI group passes.

  • % pass_th [type(group_pass_type)] :: handles for halo passes of Txy and Txx, Tyy

  • % pass_xx [type(group_pass_type)] :: MPI group passes.

  • % pass_xy [type(group_pass_type)] :: handles for halo passes of sh_xx and sh_xy, vort_xy

  • % stress_halo [integer] :: The halo size in filter of the stress tensor.

  • % hpf_halo [integer] :: The halo size in filter of the velocity gradient.

  • % amplitude [real] :: The nondimensional scaling factor in ZB model, typically 0.1 - 10 [nondim].

  • % zb_type [integer] :: Select how to compute the trace part of ZB model: 0 - both deviatoric and trace components are computed 1 - only deviatoric component is computed 2 - only trace component is computed.

  • % zb_cons [integer] :: Select a discretization scheme for ZB model 0 - non-conservative scheme 1 - conservative scheme for deviatoric component.

  • % hpf_iter [integer] :: Number of sharpening passes for the Velocity Gradient (VG) components in ZB model.

  • % stress_iter [integer] :: Number of smoothing passes for the Stress tensor components in ZB model.

  • % klower_r_diss [real] :: Attenuation of the ZB parameterization in the regions of geostrophically-unbalanced flows (Klower 2018, Juricke2020,2019) Subgrid stress is multiplied by 1/(1+(shear/(f*R_diss))) R_diss=-1: attenuation is not used; typical value R_diss=1.0 [nondim].

  • % klower_shear [integer] :: Type of expression for shear in Klower formula 0: sqrt(sh_xx**2 + sh_xy**2) 1: sqrt(sh_xx**2 + sh_xy**2 + vort_xy**2)

  • % marching_halo [integer] :: The number of filter iterations per a single MPI exchange.

  • % sh_xx [real(:,:,:),allocatable] :: Horizontal tension (du/dx - dv/dy) in h (CENTER)

  • % sh_xy [real(:,:,:),allocatable] :: Horizontal shearing strain (du/dy + dv/dx) in q (CORNER)

  • % vort_xy [real(:,:,:),allocatable] :: Vertical vorticity (dv/dx - du/dy) in q (CORNER)

  • % hq [real(:,:,:),allocatable] :: Thickness in CORNER points [H ~> m or kg m-2].

  • % txx [real(:,:,:),allocatable] :: Subgrid stress xx component in h [L2 T-2 ~> m2 s-2].

  • % tyy [real(:,:,:),allocatable] :: Subgrid stress yy component in h [L2 T-2 ~> m2 s-2].

  • % txy [real(:,:,:),allocatable] :: Subgrid stress xy component in q [L2 T-2 ~> m2 s-2].

  • % kappa_h [real(:,:),allocatable] :: Scaling coefficient in h points [L2 ~> m2].

  • % kappa_q [real(:,:),allocatable] :: Scaling coefficient in q points [L2 ~> m2].

  • % icoriolis_h [real(:,:),allocatable] :: Inverse Coriolis parameter at h points [T ~> s].

  • % c_diss [real(:,:,:),allocatable] :: Attenuation parameter at h points.

  • % maskw_h [real(:,:),allocatable] :: Mask of land point at h points multiplied by filter weight [nondim].

  • % maskw_q [real(:,:),allocatable] :: Same mask but for q points [nondim].

  • % diag [type( diag_ctrl ),pointer] :: A type that regulates diagnostics output.

Function/Subroutine Documentation

subroutine mom_zanna_bolton/zb2020_init(Time, G, GV, US, param_file, diag, CS, use_ZB2020)

Read parameters, allocate and precompute arrays, register diagnosicts used in Zanna_Bolton_2020().

Parameters:
  • time :: [in] The current model time.

  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

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

  • param_file :: [in] Parameter file parser structure.

  • diag :: [inout] Diagnostics structure.

  • cs :: [inout] ZB2020 control structure.

  • use_zb2020 :: [out] If true, turns on ZB scheme.

Call to:

mom_diag_mediator::register_diag_field

Called from:

mom_hor_visc::hor_visc_init

subroutine mom_zanna_bolton/zb2020_end(CS)

Deallocate any variables allocated in ZB_2020_init.

Parameters:

cs :: [inout] ZB2020 control structure.

Called from:

mom_hor_visc::hor_visc_end

subroutine mom_zanna_bolton/zb2020_copy_gradient_and_thickness(sh_xx, sh_xy, vort_xy, hq, G, GV, CS, k)

Save precomputed velocity gradients and thickness from the horizontal eddy viscosity module We save as much halo for velocity gradients as possible In symmetric (preferable) memory model: halo 2 for sh_xx and halo 1 for sh_xy and vort_xy We apply zero boundary conditions to velocity gradients which is required for filtering operations.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

  • cs :: [inout] ZB2020 control structure.

  • sh_xy :: [in] horizontal shearing strain (du/dy + dv/dx)

  • vort_xy :: [in] Vertical vorticity (dv/dx - du/dy)

  • hq :: [in] harmonic mean of the harmonic means

  • sh_xx :: [in] horizontal tension (du/dx - dv/dy)

  • k :: [in] The vertical index of the layer to be passed.

Called from:

mom_hor_visc::horizontal_viscosity

subroutine mom_zanna_bolton/zb2020_lateral_stress(u, v, h, diffu, diffv, G, GV, CS, dx2h, dy2h, dx2q, dy2q)

Baroclinic Zanna-Bolton-2020 parameterization, see eq. 6 in https://laurezanna.github.io/files/Zanna-Bolton-2020.pdf We compute the lateral stress tensor according to ZB2020 model and update the acceleration due to eddy viscosity (diffu, diffv) as follows: diffu = diffu + ZB2020u diffv = diffv + ZB2020v.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

  • cs :: [inout] ZB2020 control structure.

  • u :: [in] The zonal velocity [L T-1 ~> m s-1].

  • v :: [in] The meridional velocity [L T-1 ~> m s-1].

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

  • diffu :: [inout] Zonal acceleration due to eddy viscosity.

  • diffv :: [inout] Meridional acceleration due to eddy viscosity.

  • dx2h :: [in] dx^2 at h points [L2 ~> m2]

  • dy2h :: [in] dy^2 at h points [L2 ~> m2]

  • dx2q :: [in] dx^2 at q points [L2 ~> m2]

  • dy2q :: [in] dy^2 at q points [L2 ~> m2]

Call to:

compute_c_diss compute_stress compute_stress_divergence filter_stress filter_velocity_gradients

Called from:

mom_hor_visc::horizontal_viscosity

subroutine mom_zanna_bolton/compute_c_diss(G, GV, CS)

Compute the attenuation parameter similarly to Klower2018, Juricke2019,2020: c_diss = 1/(1+(shear/(f*R_diss))) where shear = sqrt(sh_xx**2 + sh_xy**2) or shear = sqrt(sh_xx**2 + sh_xy**2 + vort_xy**2) In symmetric memory model, components of velocity gradient tensor should have halo 1 and zero boundary conditions. The result: c_diss having halo 1.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

  • cs :: [inout] ZB2020 control structure.

Called from:

zb2020_lateral_stress

subroutine mom_zanna_bolton/compute_stress(G, GV, CS)

Compute stress tensor T = (Txx, Txy; Txy, Tyy) Which consists of the deviatoric and trace components, respectively: T = (-vort_xy * sh_xy, vort_xy * sh_xx; vort_xy * sh_xx, vort_xy * sh_xy) + 1/2 * (vort_xy^2 + sh_xy^2 + sh_xx^2, 0; 0, vort_xy^2 + sh_xy^2 + sh_xx^2) This stress tensor is multiplied by precomputed kappa=-CSamplitude * Garea: T -> T * kappa The sign of the stress tensor is such that (neglecting h): (du/dt, dv/dt) = div(T) In symmetric memory model: sh_xy and vort_xy should have halo 1 and zero B.C.; sh_xx should have halo 2 and zero B.C. Result: Txx, Tyy, Txy with halo 1 and zero B.C.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

  • cs :: [inout] ZB2020 control structure.

Called from:

zb2020_lateral_stress

subroutine mom_zanna_bolton/compute_stress_divergence(u, v, h, diffu, diffv, dx2h, dy2h, dx2q, dy2q, G, GV, CS)

Compute the divergence of subgrid stress weighted with thickness, i.e. (fx,fy) = 1/h Div(h * [Txx, Txy; Txy, Tyy]) and update the acceleration due to eddy viscosity as diffu = diffu + dx; diffv = diffv + dy Optionally, before computing the divergence, we attenuate the stress according to the Klower formula. In symmetric memory model: Txx, Tyy, Txy, c_diss should have halo 1 with applied zero B.C.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

  • cs :: [in] ZB2020 control structure.

  • u :: [in] The zonal velocity [L T-1 ~> m s-1].

  • v :: [in] The meridional velocity [L T-1 ~> m s-1].

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

  • diffu :: [out] Zonal acceleration due to convergence of

  • diffv :: [out] Meridional acceleration due to convergence

  • dx2h :: [in] dx^2 at h points [L2 ~> m2]

  • dy2h :: [in] dy^2 at h points [L2 ~> m2]

  • dx2q :: [in] dx^2 at q points [L2 ~> m2]

  • dy2q :: [in] dy^2 at q points [L2 ~> m2]

Call to:

compute_energy_source

Called from:

zb2020_lateral_stress

subroutine mom_zanna_bolton/filter_velocity_gradients(G, GV, CS)

Filtering of the velocity gradients sh_xx, sh_xy, vort_xy. Here instead of smoothing we do sharpening, i.e. return (initial - smoothed) fields. The algorithm: marching halo with non-blocking grouped MPI exchanges. The input array sh_xx should have halo 2 with applied zero B.C. The arrays sh_xy and vort_xy should have halo 1 with applied B.C. The output have the same halo and B.C.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

  • cs :: [inout] ZB2020 control structure.

Call to:

filter_hq

Called from:

zb2020_lateral_stress

subroutine mom_zanna_bolton/filter_stress(G, GV, CS)

Filtering of the stress tensor Txx, Tyy, Txy. The algorithm: marching halo with non-blocking grouped MPI exchanges. The input arrays (Txx, Tyy, Txy) must have halo 1 with zero B.C. applied. The output have the same halo and B.C.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

  • cs :: [inout] ZB2020 control structure.

Call to:

filter_hq

Called from:

zb2020_lateral_stress

subroutine mom_zanna_bolton/filter_hq(G, GV, CS, current_halo, remaining_iterations, q, h)

Wrapper for filter_3D function. The border indices for q and h arrays are substituted.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

  • cs :: [in] ZB2020 control structure.

  • h :: [inout] Input/output array in h points [arbitrary]

  • q :: [inout] Input/output array in q points [arbitrary]

  • current_halo :: [inout] Currently available halo points

  • remaining_iterations :: [inout] The number of iterations to perform

Call to:

filter_3d

Called from:

filter_stress filter_velocity_gradients

subroutine mom_zanna_bolton/filter_3d(x, maskw, isd, ied, jsd, jed, is, ie, js, je, nz, current_halo, remaining_iterations, direction)

Spatial lateral filter applied to 3D array. The lateral filter is given by the convolutional kernel: [1 2 1] C = |2 4 2| * 1/16 [1 2 1] The fast algorithm decomposes the 2D filter into two 1D filters as follows: [1] C = |2| * [1 2 1] * 1/16 [1] The input array must have zero B.C. applied. B.C. is applied for output array. Note that maskw contains both land mask and 1/16 factor. Filter implements marching halo. The available halo is specified and as many filter iterations as possible and as needed are performed.

Parameters:
  • isd :: [in] Indices of array size

  • ied :: [in] Indices of array size

  • jsd :: [in] Indices of array size

  • jed :: [in] Indices of array size

  • is :: [in] Indices of owned points

  • ie :: [in] Indices of owned points

  • js :: [in] Indices of owned points

  • je :: [in] Indices of owned points

  • nz :: [in] Vertical array size

  • x :: [inout] Input/output array [arbitrary]

  • maskw :: [in] Mask array of land points divided by 16 [nondim]

  • current_halo :: [inout] Currently available halo points

  • remaining_iterations :: [inout] The number of iterations to perform

  • direction :: [in] The direction of the first 1D filter

Called from:

filter_hq

subroutine mom_zanna_bolton/compute_energy_source(u, v, h, fx, fy, G, GV, CS)

Computes the 3D energy source term for the ZB2020 scheme similarly to MOM_diagnostics.F90, specifically 1125 line.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

  • cs :: [in] ZB2020 control structure.

  • u :: [in] The zonal velocity [L T-1 ~> m s-1].

  • v :: [in] The meridional velocity [L T-1 ~> m s-1].

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

  • fx :: [in] Zonal acceleration due to convergence of

  • fy :: [in] Meridional acceleration due to convergence

Called from:

compute_stress_divergence