MOM6
MOM_lateral_mixing_coeffs.F90
1 !> Variable mixing coefficients
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_debugging, only : hchksum, uvchksum
7 use mom_error_handler, only : mom_error, fatal, warning, mom_mesg
8 use mom_diag_mediator, only : register_diag_field, safe_alloc_ptr, post_data
9 use mom_diag_mediator, only : diag_ctrl, time_type, query_averaging_enabled
10 use mom_domains, only : create_group_pass, do_group_pass
11 use mom_domains, only : group_pass_type, pass_var, pass_vector
14 use mom_isopycnal_slopes, only : calc_isoneutral_slopes
15 use mom_grid, only : ocean_grid_type
19 use mom_wave_speed, only : wave_speed, wave_speed_cs, wave_speed_init
20 use mom_open_boundary, only : ocean_obc_type, obc_none
21 
22 implicit none ; private
23 
24 #include <MOM_memory.h>
25 
26 !> Variable mixing coefficients
27 type, public :: varmix_cs
28  logical :: use_variable_mixing !< If true, use the variable mixing.
29  logical :: resoln_scaled_kh !< If true, scale away the Laplacian viscosity
30  !! when the deformation radius is well resolved.
31  logical :: resoln_scaled_khth !< If true, scale away the thickness diffusivity
32  !! when the deformation radius is well resolved.
33  logical :: depth_scaled_khth !< If true, KHTH is scaled away when the depth is
34  !! shallower than a reference depth.
35  logical :: resoln_scaled_khtr !< If true, scale away the tracer diffusivity
36  !! when the deformation radius is well resolved.
37  logical :: interpolate_res_fn !< If true, interpolate the resolution function
38  !! to the velocity points from the thickness
39  !! points; otherwise interpolate the wave
40  !! speed and calculate the resolution function
41  !! independently at each point.
42  logical :: use_stored_slopes !< If true, stores isopycnal slopes in this structure.
43  logical :: resoln_use_ebt !< If true, uses the equivalent barotropic wave speed instead
44  !! of first baroclinic wave for calculating the resolution fn.
45  logical :: khth_use_ebt_struct !< If true, uses the equivalent barotropic structure
46  !! as the vertical structure of thickness diffusivity.
47  logical :: calculate_cg1 !< If true, calls wave_speed() to calculate the first
48  !! baroclinic wave speed and populate CS%cg1.
49  !! This parameter is set depending on other parameters.
50  logical :: calculate_rd_dx !< If true, calculates Rd/dx and populate CS%Rd_dx_h.
51  !! This parameter is set depending on other parameters.
52  logical :: calculate_res_fns !< If true, calculate all the resolution factors.
53  !! This parameter is set depending on other parameters.
54  logical :: calculate_depth_fns !< If true, calculate all the depth factors.
55  !! This parameter is set depending on other parameters.
56  logical :: calculate_eady_growth_rate !< If true, calculate all the Eady growth rate.
57  !! This parameter is set depending on other parameters.
58  real, dimension(:,:), pointer :: &
59  sn_u => null(), & !< S*N at u-points [T-1 ~> s-1]
60  sn_v => null(), & !< S*N at v-points [T-1 ~> s-1]
61  l2u => null(), & !< Length scale^2 at u-points [L2 ~> m2]
62  l2v => null(), & !< Length scale^2 at v-points [L2 ~> m2]
63  cg1 => null(), & !< The first baroclinic gravity wave speed [L T-1 ~> m s-1].
64  res_fn_h => null(), & !< Non-dimensional function of the ratio the first baroclinic
65  !! deformation radius to the grid spacing at h points [nondim].
66  res_fn_q => null(), & !< Non-dimensional function of the ratio the first baroclinic
67  !! deformation radius to the grid spacing at q points [nondim].
68  res_fn_u => null(), & !< Non-dimensional function of the ratio the first baroclinic
69  !! deformation radius to the grid spacing at u points [nondim].
70  res_fn_v => null(), & !< Non-dimensional function of the ratio the first baroclinic
71  !! deformation radius to the grid spacing at v points [nondim].
72  depth_fn_u => null(), & !< Non-dimensional function of the ratio of the depth to
73  !! a reference depth (maximum 1) at u points [nondim]
74  depth_fn_v => null(), & !< Non-dimensional function of the ratio of the depth to
75  !! a reference depth (maximum 1) at v points [nondim]
76  beta_dx2_h => null(), & !< The magnitude of the gradient of the Coriolis parameter
77  !! times the grid spacing squared at h points [L T-1 ~> m s-1].
78  beta_dx2_q => null(), & !< The magnitude of the gradient of the Coriolis parameter
79  !! times the grid spacing squared at q points [L T-1 ~> m s-1].
80  beta_dx2_u => null(), & !< The magnitude of the gradient of the Coriolis parameter
81  !! times the grid spacing squared at u points [L T-1 ~> m s-1].
82  beta_dx2_v => null(), & !< The magnitude of the gradient of the Coriolis parameter
83  !! times the grid spacing squared at v points [L T-1 ~> m s-1].
84  f2_dx2_h => null(), & !< The Coriolis parameter squared times the grid
85  !! spacing squared at h [L2 T-2 ~> m2 s-2].
86  f2_dx2_q => null(), & !< The Coriolis parameter squared times the grid
87  !! spacing squared at q [L2 T-2 ~> m2 s-2].
88  f2_dx2_u => null(), & !< The Coriolis parameter squared times the grid
89  !! spacing squared at u [L2 T-2 ~> m2 s-2].
90  f2_dx2_v => null(), & !< The Coriolis parameter squared times the grid
91  !! spacing squared at v [L2 T-2 ~> m2 s-2].
92  rd_dx_h => null() !< Deformation radius over grid spacing [nondim]
93 
94  real, dimension(:,:,:), pointer :: &
95  slope_x => null(), & !< Zonal isopycnal slope [nondim]
96  slope_y => null(), & !< Meridional isopycnal slope [nondim]
97  ebt_struct => null() !< Vertical structure function to scale diffusivities with [nondim]
98  real allocable_, dimension(NIMEMB_PTR_,NJMEM_) :: &
99  laplac3_const_u !< Laplacian metric-dependent constants [L3 ~> m3]
100 
101  real allocable_, dimension(NIMEM_,NJMEMB_PTR_) :: &
102  laplac3_const_v !< Laplacian metric-dependent constants [L3 ~> m3]
103 
104  real allocable_, dimension(NIMEMB_PTR_,NJMEM_,NKMEM_) :: &
105  kh_u_qg !< QG Leith GM coefficient at u-points [L2 T-1 ~> m2 s-1]
106 
107  real allocable_, dimension(NIMEM_,NJMEMB_PTR_,NKMEM_) :: &
108  kh_v_qg !< QG Leith GM coefficient at v-points [L2 T-1 ~> m2 s-1]
109 
110  ! Parameters
111  logical :: use_visbeck !< Use Visbeck formulation for thickness diffusivity
112  integer :: varmix_ktop !< Top layer to start downward integrals
113  real :: visbeck_l_scale !< Fixed length scale in Visbeck formula
114  real :: res_coef_khth !< A non-dimensional number that determines the function
115  !! of resolution, used for thickness and tracer mixing, as:
116  !! F = 1 / (1 + (Res_coef_khth*Ld/dx)^Res_fn_power)
117  real :: res_coef_visc !< A non-dimensional number that determines the function
118  !! of resolution, used for lateral viscosity, as:
119  !! F = 1 / (1 + (Res_coef_visc*Ld/dx)^Res_fn_power)
120  real :: depth_scaled_khth_h0 !< The depth above which KHTH is linearly scaled away [Z ~> m]
121  real :: depth_scaled_khth_exp !< The exponent used in the depth dependent scaling function for KHTH [nondim]
122  real :: kappa_smooth !< A diffusivity for smoothing T/S in vanished layers [Z2 T-1 ~> m2 s-1]
123  integer :: res_fn_power_khth !< The power of dx/Ld in the KhTh resolution function. Any
124  !! positive integer power may be used, but even powers
125  !! and especially 2 are coded to be more efficient.
126  integer :: res_fn_power_visc !< The power of dx/Ld in the Kh resolution function. Any
127  !! positive integer power may be used, but even powers
128  !! and especially 2 are coded to be more efficient.
129  real :: visbeck_s_max !< Upper bound on slope used in Eady growth rate [nondim].
130 
131  ! Leith parameters
132  logical :: use_qg_leith_gm !< If true, uses the QG Leith viscosity as the GM coefficient
133  logical :: use_beta_in_qg_leith !< If true, includes the beta term in the QG Leith GM coefficient
134 
135  ! Diagnostics
136  !>@{
137  !! Diagnostic identifier
138  integer :: id_sn_u=-1, id_sn_v=-1, id_l2u=-1, id_l2v=-1, id_res_fn = -1
139  integer :: id_n2_u=-1, id_n2_v=-1, id_s2_u=-1, id_s2_v=-1
140  integer :: id_rd_dx=-1, id_kh_u_qg = -1, id_kh_v_qg = -1
141  type(diag_ctrl), pointer :: diag !< A structure that is used to regulate the
142  !! timing of diagnostic output.
143  !>@}
144 
145  type(wave_speed_cs), pointer :: wave_speed_csp => null() !< Wave speed control structure
146  type(group_pass_type) :: pass_cg1 !< For group halo pass
147  logical :: debug !< If true, write out checksums of data for debugging
148 end type varmix_cs
149 
150 public varmix_init, calc_slope_functions, calc_resoln_function
151 public calc_qg_leith_viscosity, calc_depth_function
152 
153 contains
154 
155 !> Calculates the non-dimensional depth functions.
156 subroutine calc_depth_function(G, CS)
157  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
158  type(varmix_cs), pointer :: cs !< Variable mixing coefficients
159 
160  ! Local variables
161  integer :: is, ie, js, je, isq, ieq, jsq, jeq
162  integer :: i, j
163  real :: h0 ! local variable for reference depth
164  real :: expo ! exponent used in the depth dependent scaling
165  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
166  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
167 
168  if (.not. associated(cs)) call mom_error(fatal, "calc_depth_function:"// &
169  "Module must be initialized before it is used.")
170  if (.not. cs%calculate_depth_fns) return
171  if (.not. associated(cs%Depth_fn_u)) call mom_error(fatal, &
172  "calc_depth_function: %Depth_fn_u is not associated with Depth_scaled_KhTh.")
173  if (.not. associated(cs%Depth_fn_v)) call mom_error(fatal, &
174  "calc_depth_function: %Depth_fn_v is not associated with Depth_scaled_KhTh.")
175 
176  h0 = cs%depth_scaled_khth_h0
177  expo = cs%depth_scaled_khth_exp
178 !$OMP do
179  do j=js,je ; do i=is-1,ieq
180  cs%Depth_fn_u(i,j) = (min(1.0, 0.5*(g%bathyT(i,j) + g%bathyT(i+1,j))/h0))**expo
181  enddo ; enddo
182 !$OMP do
183  do j=js-1,jeq ; do i=is,ie
184  cs%Depth_fn_v(i,j) = (min(1.0, 0.5*(g%bathyT(i,j) + g%bathyT(i,j+1))/h0))**expo
185  enddo ; enddo
186 
187 end subroutine calc_depth_function
188 
189 !> Calculates and stores the non-dimensional resolution functions
190 subroutine calc_resoln_function(h, tv, G, GV, US, CS)
191  type(ocean_grid_type), intent(inout) :: g !< Ocean grid structure
192  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
193  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables
194  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
195  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
196  type(varmix_cs), pointer :: cs !< Variable mixing coefficients
197 
198  ! Local variables
199  ! Depending on the power-function being used, dimensional rescaling may be limited, so some
200  ! of the following variables have units that depend on that power.
201  real :: cg1_q ! The gravity wave speed interpolated to q points [L T-1 ~> m s-1] or [m s-1].
202  real :: cg1_u ! The gravity wave speed interpolated to u points [L T-1 ~> m s-1] or [m s-1].
203  real :: cg1_v ! The gravity wave speed interpolated to v points [L T-1 ~> m s-1] or [m s-1].
204  real :: dx_term ! A term in the denominator [L2 T-2 ~> m2 s-2] or [m2 s-2]
205  integer :: power_2
206  integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz
207  integer :: i, j, k
208  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
209  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
210 
211  if (.not. associated(cs)) call mom_error(fatal, "calc_resoln_function:"// &
212  "Module must be initialized before it is used.")
213  if (cs%calculate_cg1) then
214  if (.not. associated(cs%cg1)) call mom_error(fatal, &
215  "calc_resoln_function: %cg1 is not associated with Resoln_scaled_Kh.")
216  if (cs%khth_use_ebt_struct) then
217  if (.not. associated(cs%ebt_struct)) call mom_error(fatal, &
218  "calc_resoln_function: %ebt_struct is not associated with RESOLN_USE_EBT.")
219  if (cs%Resoln_use_ebt) then
220  ! Both resolution fn and vertical structure are using EBT
221  call wave_speed(h, tv, g, gv, us, cs%cg1, cs%wave_speed_CSp, modal_structure=cs%ebt_struct)
222  else
223  ! Use EBT to get vertical structure first and then re-calculate cg1 using first baroclinic mode
224  call wave_speed(h, tv, g, gv, us, cs%cg1, cs%wave_speed_CSp, modal_structure=cs%ebt_struct, &
225  use_ebt_mode=.true.)
226  call wave_speed(h, tv, g, gv, us, cs%cg1, cs%wave_speed_CSp)
227  endif
228  call pass_var(cs%ebt_struct, g%Domain)
229  else
230  call wave_speed(h, tv, g, gv, us, cs%cg1, cs%wave_speed_CSp)
231  endif
232 
233  call create_group_pass(cs%pass_cg1, cs%cg1, g%Domain)
234  call do_group_pass(cs%pass_cg1, g%Domain)
235  endif
236 
237  ! Calculate and store the ratio between deformation radius and grid-spacing
238  ! at h-points [nondim].
239  if (cs%calculate_rd_dx) then
240  if (.not. associated(cs%Rd_dx_h)) call mom_error(fatal, &
241  "calc_resoln_function: %Rd_dx_h is not associated with calculate_rd_dx.")
242  !$OMP parallel do default(shared)
243  do j=js-1,je+1 ; do i=is-1,ie+1
244  cs%Rd_dx_h(i,j) = cs%cg1(i,j) / &
245  (sqrt(cs%f2_dx2_h(i,j) + cs%cg1(i,j)*cs%beta_dx2_h(i,j)))
246  enddo ; enddo
247  if (query_averaging_enabled(cs%diag)) then
248  if (cs%id_Rd_dx > 0) call post_data(cs%id_Rd_dx, cs%Rd_dx_h, cs%diag)
249  endif
250  endif
251 
252  if (.not. cs%calculate_res_fns) return
253 
254  if (.not. associated(cs%Res_fn_h)) call mom_error(fatal, &
255  "calc_resoln_function: %Res_fn_h is not associated with Resoln_scaled_Kh.")
256  if (.not. associated(cs%Res_fn_q)) call mom_error(fatal, &
257  "calc_resoln_function: %Res_fn_q is not associated with Resoln_scaled_Kh.")
258  if (.not. associated(cs%Res_fn_u)) call mom_error(fatal, &
259  "calc_resoln_function: %Res_fn_u is not associated with Resoln_scaled_Kh.")
260  if (.not. associated(cs%Res_fn_v)) call mom_error(fatal, &
261  "calc_resoln_function: %Res_fn_v is not associated with Resoln_scaled_Kh.")
262  if (.not. associated(cs%f2_dx2_h)) call mom_error(fatal, &
263  "calc_resoln_function: %f2_dx2_h is not associated with Resoln_scaled_Kh.")
264  if (.not. associated(cs%f2_dx2_q)) call mom_error(fatal, &
265  "calc_resoln_function: %f2_dx2_q is not associated with Resoln_scaled_Kh.")
266  if (.not. associated(cs%f2_dx2_u)) call mom_error(fatal, &
267  "calc_resoln_function: %f2_dx2_u is not associated with Resoln_scaled_Kh.")
268  if (.not. associated(cs%f2_dx2_v)) call mom_error(fatal, &
269  "calc_resoln_function: %f2_dx2_v is not associated with Resoln_scaled_Kh.")
270  if (.not. associated(cs%beta_dx2_h)) call mom_error(fatal, &
271  "calc_resoln_function: %beta_dx2_h is not associated with Resoln_scaled_Kh.")
272  if (.not. associated(cs%beta_dx2_q)) call mom_error(fatal, &
273  "calc_resoln_function: %beta_dx2_q is not associated with Resoln_scaled_Kh.")
274  if (.not. associated(cs%beta_dx2_u)) call mom_error(fatal, &
275  "calc_resoln_function: %beta_dx2_u is not associated with Resoln_scaled_Kh.")
276  if (.not. associated(cs%beta_dx2_v)) call mom_error(fatal, &
277  "calc_resoln_function: %beta_dx2_v is not associated with Resoln_scaled_Kh.")
278 
279  ! Do this calculation on the extent used in MOM_hor_visc.F90, and
280  ! MOM_tracer.F90 so that no halo update is needed.
281 
282 !$OMP parallel default(none) shared(is,ie,js,je,Ieq,Jeq,CS,US) &
283 !$OMP private(dx_term,cg1_q,power_2,cg1_u,cg1_v)
284  if (cs%Res_fn_power_visc >= 100) then
285 !$OMP do
286  do j=js-1,je+1 ; do i=is-1,ie+1
287  dx_term = cs%f2_dx2_h(i,j) + cs%cg1(i,j)*cs%beta_dx2_h(i,j)
288  if ((cs%Res_coef_visc * cs%cg1(i,j))**2 > dx_term) then
289  cs%Res_fn_h(i,j) = 0.0
290  else
291  cs%Res_fn_h(i,j) = 1.0
292  endif
293  enddo ; enddo
294 !$OMP do
295  do j=js-1,jeq ; do i=is-1,ieq
296  cg1_q = 0.25 * ((cs%cg1(i,j) + cs%cg1(i+1,j+1)) + (cs%cg1(i+1,j) + cs%cg1(i,j+1)))
297  dx_term = cs%f2_dx2_q(i,j) + cg1_q * cs%beta_dx2_q(i,j)
298  if ((cs%Res_coef_visc * cg1_q)**2 > dx_term) then
299  cs%Res_fn_q(i,j) = 0.0
300  else
301  cs%Res_fn_q(i,j) = 1.0
302  endif
303  enddo ; enddo
304  elseif (cs%Res_fn_power_visc == 2) then
305 !$OMP do
306  do j=js-1,je+1 ; do i=is-1,ie+1
307  dx_term = cs%f2_dx2_h(i,j) + cs%cg1(i,j)*cs%beta_dx2_h(i,j)
308  cs%Res_fn_h(i,j) = dx_term / (dx_term + (cs%Res_coef_visc * cs%cg1(i,j))**2)
309  enddo ; enddo
310 !$OMP do
311  do j=js-1,jeq ; do i=is-1,ieq
312  cg1_q = 0.25 * ((cs%cg1(i,j) + cs%cg1(i+1,j+1)) + (cs%cg1(i+1,j) + cs%cg1(i,j+1)))
313  dx_term = cs%f2_dx2_q(i,j) + cg1_q * cs%beta_dx2_q(i,j)
314  cs%Res_fn_q(i,j) = dx_term / (dx_term + (cs%Res_coef_visc * cg1_q)**2)
315  enddo ; enddo
316  elseif (mod(cs%Res_fn_power_visc, 2) == 0) then
317  power_2 = cs%Res_fn_power_visc / 2
318 !$OMP do
319  do j=js-1,je+1 ; do i=is-1,ie+1
320  dx_term = (us%L_T_to_m_s**2*(cs%f2_dx2_h(i,j) + cs%cg1(i,j)*cs%beta_dx2_h(i,j)))**power_2
321  cs%Res_fn_h(i,j) = dx_term / &
322  (dx_term + (cs%Res_coef_visc * us%L_T_to_m_s*cs%cg1(i,j))**cs%Res_fn_power_visc)
323  enddo ; enddo
324 !$OMP do
325  do j=js-1,jeq ; do i=is-1,ieq
326  cg1_q = 0.25 * ((cs%cg1(i,j) + cs%cg1(i+1,j+1)) + (cs%cg1(i+1,j) + cs%cg1(i,j+1)))
327  dx_term = (us%L_T_to_m_s**2*(cs%f2_dx2_q(i,j) + cg1_q * cs%beta_dx2_q(i,j)))**power_2
328  cs%Res_fn_q(i,j) = dx_term / &
329  (dx_term + (cs%Res_coef_visc * us%L_T_to_m_s*cg1_q)**cs%Res_fn_power_visc)
330  enddo ; enddo
331  else
332 !$OMP do
333  do j=js-1,je+1 ; do i=is-1,ie+1
334  dx_term = (us%L_T_to_m_s*sqrt(cs%f2_dx2_h(i,j) + &
335  cs%cg1(i,j)*cs%beta_dx2_h(i,j)))**cs%Res_fn_power_visc
336  cs%Res_fn_h(i,j) = dx_term / &
337  (dx_term + (cs%Res_coef_visc * us%L_T_to_m_s*cs%cg1(i,j))**cs%Res_fn_power_visc)
338  enddo ; enddo
339 !$OMP do
340  do j=js-1,jeq ; do i=is-1,ieq
341  cg1_q = 0.25 * ((cs%cg1(i,j) + cs%cg1(i+1,j+1)) + (cs%cg1(i+1,j) + cs%cg1(i,j+1)))
342  dx_term = (us%L_T_to_m_s*sqrt(cs%f2_dx2_q(i,j) + &
343  cg1_q * cs%beta_dx2_q(i,j)))**cs%Res_fn_power_visc
344  cs%Res_fn_q(i,j) = dx_term / &
345  (dx_term + (cs%Res_coef_visc * us%L_T_to_m_s*cg1_q)**cs%Res_fn_power_visc)
346  enddo ; enddo
347  endif
348 
349  if (cs%interpolate_Res_fn) then
350  do j=js,je ; do i=is-1,ieq
351  cs%Res_fn_u(i,j) = 0.5*(cs%Res_fn_h(i,j) + cs%Res_fn_h(i+1,j))
352  enddo ; enddo
353  do j=js-1,jeq ; do i=is,ie
354  cs%Res_fn_v(i,j) = 0.5*(cs%Res_fn_h(i,j) + cs%Res_fn_h(i,j+1))
355  enddo ; enddo
356  else ! .not.CS%interpolate_Res_fn
357  if (cs%Res_fn_power_khth >= 100) then
358 !$OMP do
359  do j=js,je ; do i=is-1,ieq
360  cg1_u = 0.5 * (cs%cg1(i,j) + cs%cg1(i+1,j))
361  dx_term = cs%f2_dx2_u(i,j) + cg1_u * cs%beta_dx2_u(i,j)
362  if ((cs%Res_coef_khth * cg1_u)**2 > dx_term) then
363  cs%Res_fn_u(i,j) = 0.0
364  else
365  cs%Res_fn_u(i,j) = 1.0
366  endif
367  enddo ; enddo
368 !$OMP do
369  do j=js-1,jeq ; do i=is,ie
370  cg1_v = 0.5 * (cs%cg1(i,j) + cs%cg1(i,j+1))
371  dx_term = cs%f2_dx2_v(i,j) + cg1_v * cs%beta_dx2_v(i,j)
372  if ((cs%Res_coef_khth * cg1_v)**2 > dx_term) then
373  cs%Res_fn_v(i,j) = 0.0
374  else
375  cs%Res_fn_v(i,j) = 1.0
376  endif
377  enddo ; enddo
378  elseif (cs%Res_fn_power_khth == 2) then
379 !$OMP do
380  do j=js,je ; do i=is-1,ieq
381  cg1_u = 0.5 * (cs%cg1(i,j) + cs%cg1(i+1,j))
382  dx_term = cs%f2_dx2_u(i,j) + cg1_u * cs%beta_dx2_u(i,j)
383  cs%Res_fn_u(i,j) = dx_term / (dx_term + (cs%Res_coef_khth * cg1_u)**2)
384  enddo ; enddo
385 !$OMP do
386  do j=js-1,jeq ; do i=is,ie
387  cg1_v = 0.5 * (cs%cg1(i,j) + cs%cg1(i,j+1))
388  dx_term = cs%f2_dx2_v(i,j) + cg1_v * cs%beta_dx2_v(i,j)
389  cs%Res_fn_v(i,j) = dx_term / (dx_term + (cs%Res_coef_khth * cg1_v)**2)
390  enddo ; enddo
391  elseif (mod(cs%Res_fn_power_khth, 2) == 0) then
392  power_2 = cs%Res_fn_power_khth / 2
393 !$OMP do
394  do j=js,je ; do i=is-1,ieq
395  cg1_u = 0.5 * (cs%cg1(i,j) + cs%cg1(i+1,j))
396  dx_term = (us%L_T_to_m_s**2 * (cs%f2_dx2_u(i,j) + cg1_u * cs%beta_dx2_u(i,j)))**power_2
397  cs%Res_fn_u(i,j) = dx_term / &
398  (dx_term + (cs%Res_coef_khth * us%L_T_to_m_s*cg1_u)**cs%Res_fn_power_khth)
399  enddo ; enddo
400 !$OMP do
401  do j=js-1,jeq ; do i=is,ie
402  cg1_v = 0.5 * (cs%cg1(i,j) + cs%cg1(i,j+1))
403  dx_term = (us%L_T_to_m_s**2 * (cs%f2_dx2_v(i,j) + cg1_v * cs%beta_dx2_v(i,j)))**power_2
404  cs%Res_fn_v(i,j) = dx_term / &
405  (dx_term + (cs%Res_coef_khth * us%L_T_to_m_s*cg1_v)**cs%Res_fn_power_khth)
406  enddo ; enddo
407  else
408 !$OMP do
409  do j=js,je ; do i=is-1,ieq
410  cg1_u = 0.5 * (cs%cg1(i,j) + cs%cg1(i+1,j))
411  dx_term = (us%L_T_to_m_s*sqrt(cs%f2_dx2_u(i,j) + &
412  cg1_u * cs%beta_dx2_u(i,j)))**cs%Res_fn_power_khth
413  cs%Res_fn_u(i,j) = dx_term / &
414  (dx_term + (cs%Res_coef_khth * us%L_T_to_m_s*cg1_u)**cs%Res_fn_power_khth)
415  enddo ; enddo
416 !$OMP do
417  do j=js-1,jeq ; do i=is,ie
418  cg1_v = 0.5 * (cs%cg1(i,j) + cs%cg1(i,j+1))
419  dx_term = (us%L_T_to_m_s*sqrt(cs%f2_dx2_v(i,j) + &
420  cg1_v * cs%beta_dx2_v(i,j)))**cs%Res_fn_power_khth
421  cs%Res_fn_v(i,j) = dx_term / &
422  (dx_term + (cs%Res_coef_khth * us%L_T_to_m_s*cg1_v)**cs%Res_fn_power_khth)
423  enddo ; enddo
424  endif
425  endif
426 !$OMP end parallel
427 
428  if (query_averaging_enabled(cs%diag)) then
429  if (cs%id_Res_fn > 0) call post_data(cs%id_Res_fn, cs%Res_fn_h, cs%diag)
430  endif
431 
432 end subroutine calc_resoln_function
433 
434 !> Calculates and stores functions of isopycnal slopes, e.g. Sx, Sy, S*N, mostly used in the Visbeck et al.
435 !! style scaling of diffusivity
436 subroutine calc_slope_functions(h, tv, dt, G, GV, US, CS, OBC)
437  type(ocean_grid_type), intent(inout) :: g !< Ocean grid structure
438  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
439  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
440  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
441  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables
442  real, intent(in) :: dt !< Time increment [T ~> s]
443  type(varmix_cs), pointer :: cs !< Variable mixing coefficients
444  type(ocean_obc_type), optional, pointer :: obc !< Open boundaries control structure.
445  ! Local variables
446  real, dimension(SZI_(G), SZJ_(G), SZK_(G)+1) :: &
447  e ! The interface heights relative to mean sea level [Z ~> m].
448  real, dimension(SZIB_(G), SZJ_(G), SZK_(G)+1) :: n2_u ! Square of Brunt-Vaisala freq at u-points [T-2 ~> s-2]
449  real, dimension(SZI_(G), SZJB_(G), SZK_(G)+1) :: n2_v ! Square of Brunt-Vaisala freq at v-points [T-2 ~> s-2]
450 
451  if (.not. associated(cs)) call mom_error(fatal, "MOM_lateral_mixing_coeffs.F90, calc_slope_functions:"//&
452  "Module must be initialized before it is used.")
453 
454  if (cs%calculate_Eady_growth_rate) then
455  call find_eta(h, tv, g, gv, us, e, halo_size=2)
456  if (cs%use_stored_slopes) then
457  call calc_isoneutral_slopes(g, gv, us, h, e, tv, dt*cs%kappa_smooth, &
458  cs%slope_x, cs%slope_y, n2_u, n2_v, 1, obc=obc)
459  call calc_visbeck_coeffs(h, cs%slope_x, cs%slope_y, n2_u, n2_v, g, gv, us, cs, obc=obc)
460 ! call calc_slope_functions_using_just_e(h, G, CS, e, .false.)
461  else
462  !call calc_isoneutral_slopes(G, GV, h, e, tv, dt*CS%kappa_smooth, CS%slope_x, CS%slope_y)
463  call calc_slope_functions_using_just_e(h, g, gv, us, cs, e, .true., obc=obc)
464  endif
465  endif
466 
467  if (query_averaging_enabled(cs%diag)) then
468  if (cs%id_SN_u > 0) call post_data(cs%id_SN_u, cs%SN_u, cs%diag)
469  if (cs%id_SN_v > 0) call post_data(cs%id_SN_v, cs%SN_v, cs%diag)
470  if (cs%id_L2u > 0) call post_data(cs%id_L2u, cs%L2u, cs%diag)
471  if (cs%id_L2v > 0) call post_data(cs%id_L2v, cs%L2v, cs%diag)
472  if (cs%calculate_Eady_growth_rate .and. cs%use_stored_slopes) then
473  if (cs%id_N2_u > 0) call post_data(cs%id_N2_u, n2_u, cs%diag)
474  if (cs%id_N2_v > 0) call post_data(cs%id_N2_v, n2_v, cs%diag)
475  endif
476  endif
477 
478 end subroutine calc_slope_functions
479 
480 !> Calculates factors used when setting diffusivity coefficients similar to Visbeck et al.
481 subroutine calc_visbeck_coeffs(h, slope_x, slope_y, N2_u, N2_v, G, GV, US, CS, OBC)
482  type(ocean_grid_type), intent(inout) :: G !< Ocean grid structure
483  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
484  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
485  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), intent(in) :: slope_x !< Zonal isoneutral slope
486  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), intent(in) :: N2_u !< Buoyancy (Brunt-Vaisala) frequency
487  !! at u-points [T-2 ~> s-2]
488  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), intent(in) :: slope_y !< Meridional isoneutral slope
489  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), intent(in) :: N2_v !< Buoyancy (Brunt-Vaisala) frequency
490  !! at v-points [T-2 ~> s-2]
491  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
492  type(varmix_cs), pointer :: CS !< Variable mixing coefficients
493  type(ocean_obc_type), optional, pointer :: OBC !< Open boundaries control structure.
494 
495  ! Local variables
496  real :: S2 ! Interface slope squared [nondim]
497  real :: N2 ! Positive buoyancy frequency or zero [T-2 ~> s-2]
498  real :: Hup, Hdn ! Thickness from above, below [H ~> m or kg m-2]
499  real :: H_geom ! The geometric mean of Hup*Hdn [H ~> m or kg m-2].
500  integer :: is, ie, js, je, nz
501  integer :: i, j, k, kb_max
502  integer :: l_seg
503  real :: S2max, wNE, wSE, wSW, wNW
504  real :: H_u(SZIB_(G)), H_v(SZI_(G))
505  real :: S2_u(SZIB_(G), SZJ_(G))
506  real :: S2_v(SZI_(G), SZJB_(G))
507  logical :: local_open_u_BC, local_open_v_BC
508 
509  if (.not. associated(cs)) call mom_error(fatal, "calc_slope_function:"// &
510  "Module must be initialized before it is used.")
511  if (.not. cs%calculate_Eady_growth_rate) return
512  if (.not. associated(cs%SN_u)) call mom_error(fatal, "calc_slope_function:"// &
513  "%SN_u is not associated with use_variable_mixing.")
514  if (.not. associated(cs%SN_v)) call mom_error(fatal, "calc_slope_function:"// &
515  "%SN_v is not associated with use_variable_mixing.")
516 
517  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
518 
519  local_open_u_bc = .false.
520  local_open_v_bc = .false.
521  if (present(obc)) then ; if (associated(obc)) then
522  local_open_u_bc = obc%open_u_BCs_exist_globally
523  local_open_v_bc = obc%open_v_BCs_exist_globally
524  endif ; endif
525 
526  s2max = cs%Visbeck_S_max**2
527 
528  !$OMP parallel do default(shared)
529  do j=js-1,je+1 ; do i=is-1,ie+1
530  cs%SN_u(i,j) = 0.0
531  cs%SN_v(i,j) = 0.0
532  enddo ; enddo
533 
534  ! To set the length scale based on the deformation radius, use wave_speed to
535  ! calculate the first-mode gravity wave speed and then blend the equatorial
536  ! and midlatitude deformation radii, using calc_resoln_function as a template.
537 
538  !$OMP parallel do default(shared) private(S2,H_u,Hdn,Hup,H_geom,N2,wNE,wSE,wSW,wNW)
539  do j = js,je
540  do i=is-1,ie
541  cs%SN_u(i,j) = 0. ; h_u(i) = 0. ; s2_u(i,j) = 0.
542  enddo
543  do k=2,nz ; do i=is-1,ie
544  hdn = sqrt( h(i,j,k) * h(i+1,j,k) )
545  hup = sqrt( h(i,j,k-1) * h(i+1,j,k-1) )
546  h_geom = sqrt( hdn * hup )
547  !H_geom = H_geom * sqrt(N2) ! WKB-ish
548  !H_geom = H_geom * N2 ! WKB-ish
549  wse = g%mask2dCv(i+1,j-1) * ( (h(i+1,j,k)*h(i+1,j-1,k)) * (h(i+1,j,k-1)*h(i+1,j-1,k-1)) )
550  wnw = g%mask2dCv(i ,j ) * ( (h(i ,j,k)*h(i ,j+1,k)) * (h(i ,j,k-1)*h(i ,j+1,k-1)) )
551  wne = g%mask2dCv(i+1,j ) * ( (h(i+1,j,k)*h(i+1,j+1,k)) * (h(i+1,j,k-1)*h(i+1,j+1,k-1)) )
552  wsw = g%mask2dCv(i ,j-1) * ( (h(i ,j,k)*h(i ,j-1,k)) * (h(i ,j,k-1)*h(i ,j-1,k-1)) )
553  s2 = slope_x(i,j,k)**2 + &
554  ((wnw*slope_y(i,j,k)**2 + wse*slope_y(i+1,j-1,k)**2) + &
555  (wne*slope_y(i+1,j,k)**2 + wsw*slope_y(i,j-1,k)**2) ) / &
556  ( ((wse+wnw) + (wne+wsw)) + gv%H_subroundoff**4 )
557  if (s2max>0.) s2 = s2 * s2max / (s2 + s2max) ! Limit S2
558 
559  n2 = max(0., n2_u(i,j,k))
560  cs%SN_u(i,j) = cs%SN_u(i,j) + sqrt( s2*n2 )*h_geom
561  s2_u(i,j) = s2_u(i,j) + s2*h_geom
562  h_u(i) = h_u(i) + h_geom
563  enddo ; enddo
564  do i=is-1,ie
565  if (h_u(i)>0.) then
566  cs%SN_u(i,j) = g%mask2dCu(i,j) * cs%SN_u(i,j) / h_u(i)
567  s2_u(i,j) = g%mask2dCu(i,j) * s2_u(i,j) / h_u(i)
568  else
569  cs%SN_u(i,j) = 0.
570  endif
571  if (local_open_u_bc) then
572  l_seg = obc%segnum_u(i,j)
573 
574  if (l_seg /= obc_none) then
575  if (obc%segment(l_seg)%open) then
576  cs%SN_u(i,j) = 0.
577  endif
578  endif
579  endif
580  enddo
581  enddo
582 
583  !$OMP parallel do default(shared) private(S2,H_v,Hdn,Hup,H_geom,N2,wNE,wSE,wSW,wNW)
584  do j = js-1,je
585  do i=is,ie
586  cs%SN_v(i,j) = 0. ; h_v(i) = 0. ; s2_v(i,j) = 0.
587  enddo
588  do k=2,nz ; do i=is,ie
589  hdn = sqrt( h(i,j,k) * h(i,j+1,k) )
590  hup = sqrt( h(i,j,k-1) * h(i,j+1,k-1) )
591  h_geom = sqrt( hdn * hup )
592  !H_geom = H_geom * sqrt(N2) ! WKB-ish
593  !H_geom = H_geom * N2 ! WKB-ish
594  wse = g%mask2dCu(i,j) * ( (h(i,j ,k)*h(i+1,j ,k)) * (h(i,j ,k-1)*h(i+1,j ,k-1)) )
595  wnw = g%mask2dCu(i-1,j+1) * ( (h(i,j+1,k)*h(i-1,j+1,k)) * (h(i,j+1,k-1)*h(i-1,j+1,k-1)) )
596  wne = g%mask2dCu(i,j+1) * ( (h(i,j+1,k)*h(i+1,j+1,k)) * (h(i,j+1,k-1)*h(i+1,j+1,k-1)) )
597  wsw = g%mask2dCu(i-1,j) * ( (h(i,j ,k)*h(i-1,j ,k)) * (h(i,j ,k-1)*h(i-1,j ,k-1)) )
598  s2 = slope_y(i,j,k)**2 + &
599  ((wse*slope_x(i,j,k)**2 + wnw*slope_x(i-1,j+1,k)**2) + &
600  (wne*slope_x(i,j+1,k)**2 + wsw*slope_x(i-1,j,k)**2) ) / &
601  ( ((wse+wnw) + (wne+wsw)) + gv%H_subroundoff**4 )
602  if (s2max>0.) s2 = s2 * s2max / (s2 + s2max) ! Limit S2
603 
604  n2 = max(0., n2_v(i,j,k))
605  cs%SN_v(i,j) = cs%SN_v(i,j) + sqrt( s2*n2 )*h_geom
606  s2_v(i,j) = s2_v(i,j) + s2*h_geom
607  h_v(i) = h_v(i) + h_geom
608  enddo ; enddo
609  do i=is,ie
610  if (h_v(i)>0.) then
611  cs%SN_v(i,j) = g%mask2dCv(i,j) * cs%SN_v(i,j) / h_v(i)
612  s2_v(i,j) = g%mask2dCv(i,j) * s2_v(i,j) / h_v(i)
613  else
614  cs%SN_v(i,j) = 0.
615  endif
616  if (local_open_v_bc) then
617  l_seg = obc%segnum_v(i,j)
618 
619  if (l_seg /= obc_none) then
620  if (obc%segment(obc%segnum_v(i,j))%open) then
621  cs%SN_v(i,j) = 0.
622  endif
623  endif
624  endif
625  enddo
626  enddo
627 
628 ! Offer diagnostic fields for averaging.
629  if (query_averaging_enabled(cs%diag)) then
630  if (cs%id_S2_u > 0) call post_data(cs%id_S2_u, s2_u, cs%diag)
631  if (cs%id_S2_v > 0) call post_data(cs%id_S2_v, s2_v, cs%diag)
632  endif
633 
634  if (cs%debug) then
635  call uvchksum("calc_Visbeck_coeffs slope_[xy]", slope_x, slope_y, g%HI, haloshift=1)
636  call uvchksum("calc_Visbeck_coeffs N2_u, N2_v", n2_u, n2_v, g%HI, &
637  scale=us%s_to_T**2, scalar_pair=.true.)
638  call uvchksum("calc_Visbeck_coeffs SN_[uv]", cs%SN_u, cs%SN_v, g%HI, &
639  scale=us%s_to_T, scalar_pair=.true.)
640  endif
641 
642 end subroutine calc_visbeck_coeffs
643 
644 !> The original calc_slope_function() that calculated slopes using
645 !! interface positions only, not accounting for density variations.
646 subroutine calc_slope_functions_using_just_e(h, G, GV, US, CS, e, calculate_slopes, OBC)
647  type(ocean_grid_type), intent(inout) :: G !< Ocean grid structure
648  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
649  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
650  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
651  type(varmix_cs), pointer :: CS !< Variable mixing coefficients
652  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(in) :: e !< Interface position [Z ~> m]
653  logical, intent(in) :: calculate_slopes !< If true, calculate slopes internally
654  !! otherwise use slopes stored in CS
655  type(ocean_obc_type), optional, pointer :: OBC !< Open boundaries control structure.
656  ! Local variables
657  real :: E_x(SZIB_(G), SZJ_(G)) ! X-slope of interface at u points [nondim] (for diagnostics)
658  real :: E_y(SZI_(G), SZJB_(G)) ! Y-slope of interface at v points [nondim] (for diagnostics)
659  real :: H_cutoff ! Local estimate of a minimum thickness for masking [H ~> m or kg m-2]
660  real :: h_neglect ! A thickness that is so small it is usually lost
661  ! in roundoff and can be neglected [H ~> m or kg m-2].
662  real :: S2 ! Interface slope squared [nondim]
663  real :: N2 ! Brunt-Vaisala frequency squared [T-2 ~> s-2]
664  real :: Hup, Hdn ! Thickness from above, below [H ~> m or kg m-2]
665  real :: H_geom ! The geometric mean of Hup*Hdn [H ~> m or kg m-2].
666  real :: one_meter ! One meter in thickness units [H ~> m or kg m-2].
667  integer :: is, ie, js, je, nz
668  integer :: i, j, k, kb_max
669  integer :: l_seg
670  real :: S2N2_u_local(SZIB_(G), SZJ_(G),SZK_(G))
671  real :: S2N2_v_local(SZI_(G), SZJB_(G),SZK_(G))
672  logical :: local_open_u_BC, local_open_v_BC
673 
674  if (.not. associated(cs)) call mom_error(fatal, "calc_slope_function:"// &
675  "Module must be initialized before it is used.")
676  if (.not. cs%calculate_Eady_growth_rate) return
677  if (.not. associated(cs%SN_u)) call mom_error(fatal, "calc_slope_function:"// &
678  "%SN_u is not associated with use_variable_mixing.")
679  if (.not. associated(cs%SN_v)) call mom_error(fatal, "calc_slope_function:"// &
680  "%SN_v is not associated with use_variable_mixing.")
681 
682  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
683 
684  local_open_u_bc = .false.
685  local_open_v_bc = .false.
686  if (present(obc)) then ; if (associated(obc)) then
687  local_open_u_bc = obc%open_u_BCs_exist_globally
688  local_open_v_bc = obc%open_v_BCs_exist_globally
689  endif ; endif
690 
691  one_meter = 1.0 * gv%m_to_H
692  h_neglect = gv%H_subroundoff
693  h_cutoff = real(2*nz) * (gv%Angstrom_H + h_neglect)
694 
695  ! To set the length scale based on the deformation radius, use wave_speed to
696  ! calculate the first-mode gravity wave speed and then blend the equatorial
697  ! and midlatitude deformation radii, using calc_resoln_function as a template.
698 
699  !$OMP parallel do default(shared) private(E_x,E_y,S2,Hdn,Hup,H_geom,N2)
700  do k=nz,cs%VarMix_Ktop,-1
701 
702  if (calculate_slopes) then
703  ! Calculate the interface slopes E_x and E_y and u- and v- points respectively
704  do j=js-1,je+1 ; do i=is-1,ie
705  e_x(i,j) = us%Z_to_L*(e(i+1,j,k)-e(i,j,k))*g%IdxCu(i,j)
706  ! Mask slopes where interface intersects topography
707  if (min(h(i,j,k),h(i+1,j,k)) < h_cutoff) e_x(i,j) = 0.
708  enddo ; enddo
709  do j=js-1,je ; do i=is-1,ie+1
710  e_y(i,j) = us%Z_to_L*(e(i,j+1,k)-e(i,j,k))*g%IdyCv(i,j)
711  ! Mask slopes where interface intersects topography
712  if (min(h(i,j,k),h(i,j+1,k)) < h_cutoff) e_y(i,j) = 0.
713  enddo ; enddo
714  else
715  do j=js-1,je+1 ; do i=is-1,ie
716  e_x(i,j) = cs%slope_x(i,j,k)
717  if (min(h(i,j,k),h(i+1,j,k)) < h_cutoff) e_x(i,j) = 0.
718  enddo ; enddo
719  do j=js-1,je ; do i=is-1,ie+1
720  e_y(i,j) = cs%slope_y(i,j,k)
721  if (min(h(i,j,k),h(i,j+1,k)) < h_cutoff) e_y(i,j) = 0.
722  enddo ; enddo
723  endif
724 
725  ! Calculate N*S*h from this layer and add to the sum
726  do j=js,je ; do i=is-1,ie
727  s2 = ( e_x(i,j)**2 + 0.25*( &
728  (e_y(i,j)**2+e_y(i+1,j-1)**2)+(e_y(i+1,j)**2+e_y(i,j-1)**2) ) )
729  hdn = 2.*h(i,j,k)*h(i,j,k-1) / (h(i,j,k) + h(i,j,k-1) + h_neglect)
730  hup = 2.*h(i+1,j,k)*h(i+1,j,k-1) / (h(i+1,j,k) + h(i+1,j,k-1) + h_neglect)
731  h_geom = sqrt(hdn*hup)
732  n2 = gv%g_prime(k)*us%L_to_Z**2 / (gv%H_to_Z * max(hdn,hup,one_meter))
733  if (min(h(i,j,k-1), h(i+1,j,k-1), h(i,j,k), h(i+1,j,k)) < h_cutoff) &
734  s2 = 0.0
735  s2n2_u_local(i,j,k) = (h_geom * gv%H_to_Z) * s2 * n2
736  enddo ; enddo
737  do j=js-1,je ; do i=is,ie
738  s2 = ( e_y(i,j)**2 + 0.25*( &
739  (e_x(i,j)**2+e_x(i-1,j+1)**2)+(e_x(i,j+1)**2+e_x(i-1,j)**2) ) )
740  hdn = 2.*h(i,j,k)*h(i,j,k-1) / (h(i,j,k) + h(i,j,k-1) + h_neglect)
741  hup = 2.*h(i,j+1,k)*h(i,j+1,k-1) / (h(i,j+1,k) + h(i,j+1,k-1) + h_neglect)
742  h_geom = sqrt(hdn*hup)
743  n2 = gv%g_prime(k)*us%L_to_Z**2 / (gv%H_to_Z * max(hdn,hup,one_meter))
744  if (min(h(i,j,k-1), h(i,j+1,k-1), h(i,j,k), h(i,j+1,k)) < h_cutoff) &
745  s2 = 0.0
746  s2n2_v_local(i,j,k) = (h_geom * gv%H_to_Z) * s2 * n2
747  enddo ; enddo
748 
749  enddo ! k
750  !$OMP parallel do default(shared)
751  do j=js,je
752  do i=is-1,ie ; cs%SN_u(i,j) = 0.0 ; enddo
753  do k=nz,cs%VarMix_Ktop,-1 ; do i=is-1,ie
754  cs%SN_u(i,j) = cs%SN_u(i,j) + s2n2_u_local(i,j,k)
755  enddo ; enddo
756  ! SN above contains S^2*N^2*H, convert to vertical average of S*N
757  do i=is-1,ie
758  !SN_u(I,j) = sqrt( SN_u(I,j) / ( max(G%bathyT(I,j), G%bathyT(I+1,j)) + GV%Angstrom_Z ) ))
759  !The code below behaves better than the line above. Not sure why? AJA
760  if ( min(g%bathyT(i,j), g%bathyT(i+1,j)) > h_cutoff*gv%H_to_Z ) then
761  cs%SN_u(i,j) = g%mask2dCu(i,j) * sqrt( cs%SN_u(i,j) / &
762  (max(g%bathyT(i,j), g%bathyT(i+1,j))) )
763  else
764  cs%SN_u(i,j) = 0.0
765  endif
766  if (local_open_u_bc) then
767  l_seg = obc%segnum_u(i,j)
768 
769  if (l_seg /= obc_none) then
770  if (obc%segment(l_seg)%open) then
771  cs%SN_u(i,j) = 0.
772  endif
773  endif
774  endif
775  enddo
776  enddo
777  !$OMP parallel do default(shared)
778  do j=js-1,je
779  do i=is,ie ; cs%SN_v(i,j) = 0.0 ; enddo
780  do k=nz,cs%VarMix_Ktop,-1 ; do i=is,ie
781  cs%SN_v(i,j) = cs%SN_v(i,j) + s2n2_v_local(i,j,k)
782  enddo ; enddo
783  do i=is,ie
784  !SN_v(i,J) = sqrt( SN_v(i,J) / ( max(G%bathyT(i,J), G%bathyT(i,J+1)) + GV%Angstrom_Z ) ))
785  !The code below behaves better than the line above. Not sure why? AJA
786  if ( min(g%bathyT(i,j), g%bathyT(i+1,j)) > h_cutoff*gv%H_to_Z ) then
787  cs%SN_v(i,j) = g%mask2dCv(i,j) * sqrt( cs%SN_v(i,j) / &
788  (max(g%bathyT(i,j), g%bathyT(i,j+1))) )
789  else
790  cs%SN_v(i,j) = 0.0
791  endif
792  if (local_open_v_bc) then
793  l_seg = obc%segnum_v(i,j)
794 
795  if (l_seg /= obc_none) then
796  if (obc%segment(obc%segnum_v(i,j))%open) then
797  cs%SN_v(i,j) = 0.
798  endif
799  endif
800  endif
801  enddo
802  enddo
803 
804 end subroutine calc_slope_functions_using_just_e
805 
806 !> Calculates the Leith Laplacian and bi-harmonic viscosity coefficients
807 subroutine calc_qg_leith_viscosity(CS, G, GV, US, h, k, div_xx_dx, div_xx_dy, vort_xy_dx, vort_xy_dy)
808  type(varmix_cs), pointer :: cs !< Variable mixing coefficients
809  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
810  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
811  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
812  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
813  integer, intent(in) :: k !< Layer for which to calculate vorticity magnitude
814  real, dimension(SZIB_(G),SZJ_(G)), intent(in) :: div_xx_dx !< x-derivative of horizontal divergence
815  !! (d/dx(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]
816  real, dimension(SZI_(G),SZJB_(G)), intent(in) :: div_xx_dy !< y-derivative of horizontal divergence
817  !! (d/dy(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]
818  real, dimension(SZI_(G),SZJB_(G)), intent(inout) :: vort_xy_dx !< x-derivative of vertical vorticity
819  !! (d/dx(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]
820  real, dimension(SZIB_(G),SZJ_(G)), intent(inout) :: vort_xy_dy !< y-derivative of vertical vorticity
821  !! (d/dy(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]
822  ! Local variables
823  real, dimension(SZI_(G),SZJB_(G)) :: &
824  dslopey_dz, & ! z-derivative of y-slope at v-points [Z-1 ~> m-1]
825  h_at_v, & ! Thickness at v-points [H ~> m or kg m-2]
826  beta_v, & ! Beta at v-points [T-1 L-1 ~> s-1 m-1]
827  grad_vort_mag_v, & ! Magnitude of vorticity gradient at v-points [T-1 L-1 ~> s-1 m-1]
828  grad_div_mag_v ! Magnitude of divergence gradient at v-points [T-1 L-1 ~> s-1 m-1]
829 
830  real, dimension(SZIB_(G),SZJ_(G)) :: &
831  dslopex_dz, & ! z-derivative of x-slope at u-points [Z-1 ~> m-1]
832  h_at_u, & ! Thickness at u-points [H ~> m or kg m-2]
833  beta_u, & ! Beta at u-points [T-1 L-1 ~> s-1 m-1]
834  grad_vort_mag_u, & ! Magnitude of vorticity gradient at u-points [T-1 L-1 ~> s-1 m-1]
835  grad_div_mag_u ! Magnitude of divergence gradient at u-points [T-1 L-1 ~> s-1 m-1]
836  real :: h_at_slope_above ! The thickness above [H ~> m or kg m-2]
837  real :: h_at_slope_below ! The thickness below [H ~> m or kg m-2]
838  real :: ih ! The inverse of a combination of thicknesses [H-1 ~> m-1 or m2 kg-1]
839  real :: f ! A copy of the Coriolis parameter [T-1 ~> s-1]
840  integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq,nz
841  real :: inv_pi3
842 
843  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
844  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
845  nz = g%ke
846 
847  inv_pi3 = 1.0/((4.0*atan(1.0))**3)
848 
849  if ((k > 1) .and. (k < nz)) then
850 
851  do j=js-1,je+1 ; do i=is-2,ieq+1
852  h_at_slope_above = 2. * ( h(i,j,k-1) * h(i+1,j,k-1) ) * ( h(i,j,k) * h(i+1,j,k) ) / &
853  ( ( h(i,j,k-1) * h(i+1,j,k-1) ) * ( h(i,j,k) + h(i+1,j,k) ) &
854  + ( h(i,j,k) * h(i+1,j,k) ) * ( h(i,j,k-1) + h(i+1,j,k-1) ) + gv%H_subroundoff )
855  h_at_slope_below = 2. * ( h(i,j,k) * h(i+1,j,k) ) * ( h(i,j,k+1) * h(i+1,j,k+1) ) / &
856  ( ( h(i,j,k) * h(i+1,j,k) ) * ( h(i,j,k+1) + h(i+1,j,k+1) ) &
857  + ( h(i,j,k+1) * h(i+1,j,k+1) ) * ( h(i,j,k) + h(i+1,j,k) ) + gv%H_subroundoff )
858  ih = 1./ ( ( h_at_slope_above + h_at_slope_below + gv%H_subroundoff ) * gv%H_to_Z )
859  dslopex_dz(i,j) = 2. * ( cs%slope_x(i,j,k) - cs%slope_x(i,j,k+1) ) * ih
860  h_at_u(i,j) = 2. * ( h_at_slope_above * h_at_slope_below ) * ih
861  enddo ; enddo
862 
863  do j=js-2,jeq+1 ; do i=is-1,ie+1
864  h_at_slope_above = 2. * ( h(i,j,k-1) * h(i,j+1,k-1) ) * ( h(i,j,k) * h(i,j+1,k) ) / &
865  ( ( h(i,j,k-1) * h(i,j+1,k-1) ) * ( h(i,j,k) + h(i,j+1,k) ) &
866  + ( h(i,j,k) * h(i,j+1,k) ) * ( h(i,j,k-1) + h(i,j+1,k-1) ) + gv%H_subroundoff )
867  h_at_slope_below = 2. * ( h(i,j,k) * h(i,j+1,k) ) * ( h(i,j,k+1) * h(i,j+1,k+1) ) / &
868  ( ( h(i,j,k) * h(i,j+1,k) ) * ( h(i,j,k+1) + h(i,j+1,k+1) ) &
869  + ( h(i,j,k+1) * h(i,j+1,k+1) ) * ( h(i,j,k) + h(i,j+1,k) ) + gv%H_subroundoff )
870  ih = 1./ ( ( h_at_slope_above + h_at_slope_below + gv%H_subroundoff ) * gv%H_to_Z )
871  dslopey_dz(i,j) = 2. * ( cs%slope_y(i,j,k) - cs%slope_y(i,j,k+1) ) * ih
872  h_at_v(i,j) = 2. * ( h_at_slope_above * h_at_slope_below ) * ih
873  enddo ; enddo
874 
875  do j=js-1,je ; do i=is-1,ieq+1
876  f = 0.5 * ( g%CoriolisBu(i,j) + g%CoriolisBu(i-1,j) )
877  vort_xy_dx(i,j) = vort_xy_dx(i,j) - f * us%L_to_Z * &
878  ( ( h_at_u(i,j) * dslopex_dz(i,j) + h_at_u(i-1,j+1) * dslopex_dz(i-1,j+1) ) &
879  + ( h_at_u(i-1,j) * dslopex_dz(i-1,j) + h_at_u(i,j+1) * dslopex_dz(i,j+1) ) ) / &
880  ( ( h_at_u(i,j) + h_at_u(i-1,j+1) ) + ( h_at_u(i-1,j) + h_at_u(i,j+1) ) + gv%H_subroundoff)
881  enddo ; enddo
882 
883  do j=js-1,jeq+1 ; do i=is-1,ie
884  f = 0.5 * ( g%CoriolisBu(i,j) + g%CoriolisBu(i,j-1) )
885  vort_xy_dy(i,j) = vort_xy_dy(i,j) - f * us%L_to_Z * &
886  ( ( h_at_v(i,j) * dslopey_dz(i,j) + h_at_v(i+1,j-1) * dslopey_dz(i+1,j-1) ) &
887  + ( h_at_v(i,j-1) * dslopey_dz(i,j-1) + h_at_v(i+1,j) * dslopey_dz(i+1,j) ) ) / &
888  ( ( h_at_v(i,j) + h_at_v(i+1,j-1) ) + ( h_at_v(i,j-1) + h_at_v(i+1,j) ) + gv%H_subroundoff)
889  enddo ; enddo
890  endif ! k > 1
891 
892  if (cs%use_QG_Leith_GM) then
893 
894  do j=js,je ; do i=is-1,ieq
895  grad_vort_mag_u(i,j) = sqrt(vort_xy_dy(i,j)**2 + (0.25*((vort_xy_dx(i,j) + vort_xy_dx(i+1,j-1)) &
896  + (vort_xy_dx(i+1,j) + vort_xy_dx(i,j-1))))**2)
897  grad_div_mag_u(i,j) = sqrt(div_xx_dx(i,j)**2 + (0.25*((div_xx_dy(i,j) + div_xx_dy(i+1,j-1)) &
898  + (div_xx_dy(i+1,j) + div_xx_dy(i,j-1))))**2)
899  if (cs%use_beta_in_QG_Leith) then
900  beta_u(i,j) = sqrt((0.5*(g%dF_dx(i,j)+g%dF_dx(i+1,j))**2) + &
901  (0.5*(g%dF_dy(i,j)+g%dF_dy(i+1,j))**2))
902  cs%KH_u_QG(i,j,k) = min(grad_vort_mag_u(i,j) + grad_div_mag_u(i,j), 3.0*beta_u(i,j)) * &
903  cs%Laplac3_const_u(i,j) * inv_pi3
904  else
905  cs%KH_u_QG(i,j,k) = (grad_vort_mag_u(i,j) + grad_div_mag_u(i,j)) * &
906  cs%Laplac3_const_u(i,j) * inv_pi3
907  endif
908  enddo ; enddo
909 
910  do j=js-1,jeq ; do i=is,ie
911  grad_vort_mag_v(i,j) = sqrt(vort_xy_dx(i,j)**2 + (0.25*((vort_xy_dy(i,j) + vort_xy_dy(i-1,j+1)) &
912  + (vort_xy_dy(i,j+1) + vort_xy_dy(i-1,j))))**2)
913  grad_div_mag_v(i,j) = sqrt(div_xx_dy(i,j)**2 + (0.25*((div_xx_dx(i,j) + div_xx_dx(i-1,j+1)) &
914  + (div_xx_dx(i,j+1) + div_xx_dx(i-1,j))))**2)
915  if (cs%use_beta_in_QG_Leith) then
916  beta_v(i,j) = sqrt((0.5*(g%dF_dx(i,j)+g%dF_dx(i,j+1))**2) + &
917  (0.5*(g%dF_dy(i,j)+g%dF_dy(i,j+1))**2))
918  cs%KH_v_QG(i,j,k) = min(grad_vort_mag_v(i,j) + grad_div_mag_v(i,j), 3.0*beta_v(i,j)) * &
919  cs%Laplac3_const_v(i,j) * inv_pi3
920  else
921  cs%KH_v_QG(i,j,k) = (grad_vort_mag_v(i,j) + grad_div_mag_v(i,j)) * &
922  cs%Laplac3_const_v(i,j) * inv_pi3
923  endif
924  enddo ; enddo
925  ! post diagnostics
926 
927  if (k==nz) then
928  if (cs%id_KH_v_QG > 0) call post_data(cs%id_KH_v_QG, cs%KH_v_QG, cs%diag)
929  if (cs%id_KH_u_QG > 0) call post_data(cs%id_KH_u_QG, cs%KH_u_QG, cs%diag)
930  endif
931  endif
932 
933 end subroutine calc_qg_leith_viscosity
934 
935 !> Initializes the variables mixing coefficients container
936 subroutine varmix_init(Time, G, GV, US, param_file, diag, CS)
937  type(time_type), intent(in) :: time !< Current model time
938  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
939  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
940  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
941  type(param_file_type), intent(in) :: param_file !< Parameter file handles
942  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure
943  type(varmix_cs), pointer :: cs !< Variable mixing coefficients
944  ! Local variables
945  real :: khtr_slope_cff, khth_slope_cff, oneortwo
946  real :: n2_filter_depth ! A depth below which stratification is treated as monotonic when
947  ! calculating the first-mode wave speed [Z ~> m]
948  real :: khtr_passivity_coeff
949  real :: absurdly_small_freq ! A miniscule frequency that is used to avoid division by 0 [T-1 ~> s-1]. The
950  ! default value is roughly (pi / (the age of the universe)).
951  logical :: gill_equatorial_ld, use_fgnv_streamfn, use_meke, in_use
952  logical :: default_2018_answers, remap_answers_2018
953  real :: mle_front_length
954  real :: leith_lap_const ! The non-dimensional coefficient in the Leith viscosity
955  real :: grid_sp_u2, grid_sp_v2 ! Intermediate quantities for Leith metrics [L2 ~> m2]
956  real :: grid_sp_u3, grid_sp_v3 ! Intermediate quantities for Leith metrics [L3 ~> m3]
957  real :: wave_speed_min ! A floor in the first mode speed below which 0 is returned [L T-1 ~> m s-1]
958  real :: wave_speed_tol ! The fractional tolerance for finding the wave speeds [nondim]
959  logical :: better_speed_est ! If true, use a more robust estimate of the first
960  ! mode wave speed as the starting point for iterations.
961 ! This include declares and sets the variable "version".
962 # include "version_variable.h"
963  character(len=40) :: mdl = "MOM_lateral_mixing_coeffs" ! This module's name.
964  integer :: is, ie, js, je, isq, ieq, jsq, jeq, i, j
965  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
966  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
967  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
968  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
969  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
970 
971  if (associated(cs)) then
972  call mom_error(warning, "VarMix_init called with an associated "// &
973  "control structure.")
974  return
975  endif
976 
977  allocate(cs)
978  in_use = .false. ! Set to true to avoid deallocating
979  cs%diag => diag ! Diagnostics pointer
980  cs%calculate_cg1 = .false.
981  cs%calculate_Rd_dx = .false.
982  cs%calculate_res_fns = .false.
983  cs%calculate_Eady_growth_rate = .false.
984  cs%calculate_depth_fns = .false.
985  ! Read all relevant parameters and write them to the model log.
986  call log_version(param_file, mdl, version, "")
987  call get_param(param_file, mdl, "USE_VARIABLE_MIXING", cs%use_variable_mixing,&
988  "If true, the variable mixing code will be called. This "//&
989  "allows diagnostics to be created even if the scheme is "//&
990  "not used. If KHTR_SLOPE_CFF>0 or KhTh_Slope_Cff>0, "//&
991  "this is set to true regardless of what is in the "//&
992  "parameter file.", default=.false.)
993  call get_param(param_file, mdl, "USE_VISBECK", cs%use_Visbeck,&
994  "If true, use the Visbeck et al. (1997) formulation for \n"//&
995  "thickness diffusivity.", default=.false.)
996  call get_param(param_file, mdl, "RESOLN_SCALED_KH", cs%Resoln_scaled_Kh, &
997  "If true, the Laplacian lateral viscosity is scaled away "//&
998  "when the first baroclinic deformation radius is well "//&
999  "resolved.", default=.false.)
1000  call get_param(param_file, mdl, "DEPTH_SCALED_KHTH", cs%Depth_scaled_KhTh, &
1001  "If true, KHTH is scaled away when the depth is shallower"//&
1002  "than a reference depth: KHTH = MIN(1,H/H0)**N * KHTH, "//&
1003  "where H0 is a reference depth, controlled via DEPTH_SCALED_KHTH_H0, "//&
1004  "and the exponent (N) is controlled via DEPTH_SCALED_KHTH_EXP.",&
1005  default=.false.)
1006  call get_param(param_file, mdl, "RESOLN_SCALED_KHTH", cs%Resoln_scaled_KhTh, &
1007  "If true, the interface depth diffusivity is scaled away "//&
1008  "when the first baroclinic deformation radius is well "//&
1009  "resolved.", default=.false.)
1010  call get_param(param_file, mdl, "RESOLN_SCALED_KHTR", cs%Resoln_scaled_KhTr, &
1011  "If true, the epipycnal tracer diffusivity is scaled "//&
1012  "away when the first baroclinic deformation radius is "//&
1013  "well resolved.", default=.false.)
1014  call get_param(param_file, mdl, "RESOLN_USE_EBT", cs%Resoln_use_ebt, &
1015  "If true, uses the equivalent barotropic wave speed instead "//&
1016  "of first baroclinic wave for calculating the resolution fn.",&
1017  default=.false.)
1018  call get_param(param_file, mdl, "KHTH_USE_EBT_STRUCT", cs%khth_use_ebt_struct, &
1019  "If true, uses the equivalent barotropic structure "//&
1020  "as the vertical structure of thickness diffusivity.",&
1021  default=.false.)
1022  call get_param(param_file, mdl, "KHTH_SLOPE_CFF", khth_slope_cff, &
1023  "The nondimensional coefficient in the Visbeck formula "//&
1024  "for the interface depth diffusivity", units="nondim", &
1025  default=0.0)
1026  call get_param(param_file, mdl, "KHTR_SLOPE_CFF", khtr_slope_cff, &
1027  "The nondimensional coefficient in the Visbeck formula "//&
1028  "for the epipycnal tracer diffusivity", units="nondim", &
1029  default=0.0)
1030  call get_param(param_file, mdl, "USE_STORED_SLOPES", cs%use_stored_slopes,&
1031  "If true, the isopycnal slopes are calculated once and "//&
1032  "stored for re-use. This uses more memory but avoids calling "//&
1033  "the equation of state more times than should be necessary.", &
1034  default=.false.)
1035  call get_param(param_file, mdl, "VERY_SMALL_FREQUENCY", absurdly_small_freq, &
1036  "A miniscule frequency that is used to avoid division by 0. The default "//&
1037  "value is roughly (pi / (the age of the universe)).", &
1038  default=1.0e-17, units="s-1", scale=us%T_to_s)
1039  call get_param(param_file, mdl, "KHTH_USE_FGNV_STREAMFUNCTION", use_fgnv_streamfn, &
1040  default=.false., do_not_log=.true.)
1041  cs%calculate_cg1 = cs%calculate_cg1 .or. use_fgnv_streamfn
1042  call get_param(param_file, mdl, "USE_MEKE", use_meke, &
1043  default=.false., do_not_log=.true.)
1044  cs%calculate_Rd_dx = cs%calculate_Rd_dx .or. use_meke
1045  cs%calculate_Eady_growth_rate = cs%calculate_Eady_growth_rate .or. use_meke
1046  call get_param(param_file, mdl, "KHTR_PASSIVITY_COEFF", khtr_passivity_coeff, &
1047  default=0., do_not_log=.true.)
1048  cs%calculate_Rd_dx = cs%calculate_Rd_dx .or. (khtr_passivity_coeff>0.)
1049  call get_param(param_file, mdl, "MLE_FRONT_LENGTH", mle_front_length, &
1050  default=0., do_not_log=.true.)
1051  cs%calculate_Rd_dx = cs%calculate_Rd_dx .or. (mle_front_length>0.)
1052 
1053  call get_param(param_file, mdl, "DEBUG", cs%debug, default=.false., do_not_log=.true.)
1054 
1055 
1056  if (cs%Resoln_use_ebt .or. cs%khth_use_ebt_struct) then
1057  in_use = .true.
1058  call get_param(param_file, mdl, "RESOLN_N2_FILTER_DEPTH", n2_filter_depth, &
1059  "The depth below which N2 is monotonized to avoid stratification "//&
1060  "artifacts from altering the equivalent barotropic mode structure.",&
1061  units="m", default=2000., scale=us%m_to_Z)
1062  allocate(cs%ebt_struct(isd:ied,jsd:jed,g%ke)) ; cs%ebt_struct(:,:,:) = 0.0
1063  endif
1064 
1065  if (khtr_slope_cff>0. .or. khth_slope_cff>0.) then
1066  cs%calculate_Eady_growth_rate = .true.
1067  call get_param(param_file, mdl, "VISBECK_MAX_SLOPE", cs%Visbeck_S_max, &
1068  "If non-zero, is an upper bound on slopes used in the "//&
1069  "Visbeck formula for diffusivity. This does not affect the "//&
1070  "isopycnal slope calculation used within thickness diffusion.", &
1071  units="nondim", default=0.0)
1072  endif
1073 
1074  if (cs%use_stored_slopes) then
1075  in_use = .true.
1076  allocate(cs%slope_x(isdb:iedb,jsd:jed,g%ke+1)) ; cs%slope_x(:,:,:) = 0.0
1077  allocate(cs%slope_y(isd:ied,jsdb:jedb,g%ke+1)) ; cs%slope_y(:,:,:) = 0.0
1078  call get_param(param_file, mdl, "KD_SMOOTH", cs%kappa_smooth, &
1079  "A diapycnal diffusivity that is used to interpolate "//&
1080  "more sensible values of T & S into thin layers.", &
1081  units="m2 s-1", default=1.0e-6, scale=us%m_to_Z**2*us%T_to_s)
1082  endif
1083 
1084  if (cs%calculate_Eady_growth_rate) then
1085  in_use = .true.
1086  allocate(cs%SN_u(isdb:iedb,jsd:jed)) ; cs%SN_u(:,:) = 0.0
1087  allocate(cs%SN_v(isd:ied,jsdb:jedb)) ; cs%SN_v(:,:) = 0.0
1088  cs%id_SN_u = register_diag_field('ocean_model', 'SN_u', diag%axesCu1, time, &
1089  'Inverse eddy time-scale, S*N, at u-points', 's-1', conversion=us%s_to_T)
1090  cs%id_SN_v = register_diag_field('ocean_model', 'SN_v', diag%axesCv1, time, &
1091  'Inverse eddy time-scale, S*N, at v-points', 's-1', conversion=us%s_to_T)
1092  call get_param(param_file, mdl, "VARMIX_KTOP", cs%VarMix_Ktop, &
1093  "The layer number at which to start vertical integration "//&
1094  "of S*N for purposes of finding the Eady growth rate.", &
1095  units="nondim", default=2)
1096  endif
1097 
1098  if (khtr_slope_cff>0. .or. khth_slope_cff>0.) then
1099  in_use = .true.
1100  call get_param(param_file, mdl, "VISBECK_L_SCALE", cs%Visbeck_L_scale, &
1101  "The fixed length scale in the Visbeck formula.", units="m", &
1102  default=0.0)
1103  allocate(cs%L2u(isdb:iedb,jsd:jed)) ; cs%L2u(:,:) = 0.0
1104  allocate(cs%L2v(isd:ied,jsdb:jedb)) ; cs%L2v(:,:) = 0.0
1105  if (cs%Visbeck_L_scale<0) then
1106  do j=js,je ; do i=is-1,ieq
1107  cs%L2u(i,j) = cs%Visbeck_L_scale**2 * g%areaCu(i,j)
1108  enddo; enddo
1109  do j=js-1,jeq ; do i=is,ie
1110  cs%L2v(i,j) = cs%Visbeck_L_scale**2 * g%areaCv(i,j)
1111  enddo; enddo
1112  else
1113  cs%L2u(:,:) = us%m_to_L**2*cs%Visbeck_L_scale**2
1114  cs%L2v(:,:) = us%m_to_L**2*cs%Visbeck_L_scale**2
1115  endif
1116 
1117  cs%id_L2u = register_diag_field('ocean_model', 'L2u', diag%axesCu1, time, &
1118  'Length scale squared for mixing coefficient, at u-points', &
1119  'm2', conversion=us%L_to_m**2)
1120  cs%id_L2v = register_diag_field('ocean_model', 'L2v', diag%axesCv1, time, &
1121  'Length scale squared for mixing coefficient, at v-points', &
1122  'm2', conversion=us%L_to_m**2)
1123  endif
1124 
1125  if (cs%calculate_Eady_growth_rate .and. cs%use_stored_slopes) then
1126  cs%id_N2_u = register_diag_field('ocean_model', 'N2_u', diag%axesCui, time, &
1127  'Square of Brunt-Vaisala frequency, N^2, at u-points, as used in Visbeck et al.', &
1128  's-2', conversion=us%s_to_T**2)
1129  cs%id_N2_v = register_diag_field('ocean_model', 'N2_v', diag%axesCvi, time, &
1130  'Square of Brunt-Vaisala frequency, N^2, at v-points, as used in Visbeck et al.', &
1131  's-2', conversion=us%s_to_T**2)
1132  endif
1133  if (cs%use_stored_slopes) then
1134  cs%id_S2_u = register_diag_field('ocean_model', 'S2_u', diag%axesCu1, time, &
1135  'Depth average square of slope magnitude, S^2, at u-points, as used in Visbeck et al.', 'nondim')
1136  cs%id_S2_v = register_diag_field('ocean_model', 'S2_v', diag%axesCv1, time, &
1137  'Depth average square of slope magnitude, S^2, at v-points, as used in Visbeck et al.', 'nondim')
1138  endif
1139 
1140  oneortwo = 1.0
1141  if (cs%Resoln_scaled_Kh .or. cs%Resoln_scaled_KhTh .or. cs%Resoln_scaled_KhTr) then
1142  cs%calculate_Rd_dx = .true.
1143  cs%calculate_res_fns = .true.
1144  allocate(cs%Res_fn_h(isd:ied,jsd:jed)) ; cs%Res_fn_h(:,:) = 0.0
1145  allocate(cs%Res_fn_q(isdb:iedb,jsdb:jedb)) ; cs%Res_fn_q(:,:) = 0.0
1146  allocate(cs%Res_fn_u(isdb:iedb,jsd:jed)) ; cs%Res_fn_u(:,:) = 0.0
1147  allocate(cs%Res_fn_v(isd:ied,jsdb:jedb)) ; cs%Res_fn_v(:,:) = 0.0
1148  allocate(cs%beta_dx2_q(isdb:iedb,jsdb:jedb)) ; cs%beta_dx2_q(:,:) = 0.0
1149  allocate(cs%beta_dx2_u(isdb:iedb,jsd:jed)) ; cs%beta_dx2_u(:,:) = 0.0
1150  allocate(cs%beta_dx2_v(isd:ied,jsdb:jedb)) ; cs%beta_dx2_v(:,:) = 0.0
1151  allocate(cs%f2_dx2_q(isdb:iedb,jsdb:jedb)) ; cs%f2_dx2_q(:,:) = 0.0
1152  allocate(cs%f2_dx2_u(isdb:iedb,jsd:jed)) ; cs%f2_dx2_u(:,:) = 0.0
1153  allocate(cs%f2_dx2_v(isd:ied,jsdb:jedb)) ; cs%f2_dx2_v(:,:) = 0.0
1154 
1155  cs%id_Res_fn = register_diag_field('ocean_model', 'Res_fn', diag%axesT1, time, &
1156  'Resolution function for scaling diffusivities', 'nondim')
1157 
1158  call get_param(param_file, mdl, "KH_RES_SCALE_COEF", cs%Res_coef_khth, &
1159  "A coefficient that determines how KhTh is scaled away if "//&
1160  "RESOLN_SCALED_... is true, as "//&
1161  "F = 1 / (1 + (KH_RES_SCALE_COEF*Rd/dx)^KH_RES_FN_POWER).", &
1162  units="nondim", default=1.0)
1163  call get_param(param_file, mdl, "KH_RES_FN_POWER", cs%Res_fn_power_khth, &
1164  "The power of dx/Ld in the Kh resolution function. Any "//&
1165  "positive integer may be used, although even integers "//&
1166  "are more efficient to calculate. Setting this greater "//&
1167  "than 100 results in a step-function being used.", &
1168  units="nondim", default=2)
1169  call get_param(param_file, mdl, "VISC_RES_SCALE_COEF", cs%Res_coef_visc, &
1170  "A coefficient that determines how Kh is scaled away if "//&
1171  "RESOLN_SCALED_... is true, as "//&
1172  "F = 1 / (1 + (KH_RES_SCALE_COEF*Rd/dx)^KH_RES_FN_POWER). "//&
1173  "This function affects lateral viscosity, Kh, and not KhTh.", &
1174  units="nondim", default=cs%Res_coef_khth)
1175  call get_param(param_file, mdl, "VISC_RES_FN_POWER", cs%Res_fn_power_visc, &
1176  "The power of dx/Ld in the Kh resolution function. Any "//&
1177  "positive integer may be used, although even integers "//&
1178  "are more efficient to calculate. Setting this greater "//&
1179  "than 100 results in a step-function being used. "//&
1180  "This function affects lateral viscosity, Kh, and not KhTh.", &
1181  units="nondim", default=cs%Res_fn_power_khth)
1182  call get_param(param_file, mdl, "INTERPOLATE_RES_FN", cs%interpolate_Res_fn, &
1183  "If true, interpolate the resolution function to the "//&
1184  "velocity points from the thickness points; otherwise "//&
1185  "interpolate the wave speed and calculate the resolution "//&
1186  "function independently at each point.", default=.false.)
1187  if (cs%interpolate_Res_fn) then
1188  if (cs%Res_coef_visc /= cs%Res_coef_khth) call mom_error(fatal, &
1189  "MOM_lateral_mixing_coeffs.F90, VarMix_init:"//&
1190  "When INTERPOLATE_RES_FN=True, VISC_RES_FN_POWER must equal KH_RES_SCALE_COEF.")
1191  if (cs%Res_fn_power_visc /= cs%Res_fn_power_khth) call mom_error(fatal, &
1192  "MOM_lateral_mixing_coeffs.F90, VarMix_init:"//&
1193  "When INTERPOLATE_RES_FN=True, VISC_RES_FN_POWER must equal KH_RES_FN_POWER.")
1194  endif
1195  call get_param(param_file, mdl, "GILL_EQUATORIAL_LD", gill_equatorial_ld, &
1196  "If true, uses Gill's definition of the baroclinic "//&
1197  "equatorial deformation radius, otherwise, if false, use "//&
1198  "Pedlosky's definition. These definitions differ by a factor "//&
1199  "of 2 in front of the beta term in the denominator. Gill's "//&
1200  "is the more appropriate definition.", default=.true.)
1201  if (gill_equatorial_ld) then
1202  oneortwo = 2.0
1203  endif
1204 
1205  do j=js-1,jeq ; do i=is-1,ieq
1206  cs%f2_dx2_q(i,j) = (g%dxBu(i,j)**2 + g%dyBu(i,j)**2) * &
1207  max(g%CoriolisBu(i,j)**2, absurdly_small_freq**2)
1208  cs%beta_dx2_q(i,j) = oneortwo * ((g%dxBu(i,j))**2 + (g%dyBu(i,j))**2) * (sqrt(0.5 * &
1209  ( (((g%CoriolisBu(i,j)-g%CoriolisBu(i-1,j)) * g%IdxCv(i,j))**2 + &
1210  ((g%CoriolisBu(i+1,j)-g%CoriolisBu(i,j)) * g%IdxCv(i+1,j))**2) + &
1211  (((g%CoriolisBu(i,j)-g%CoriolisBu(i,j-1)) * g%IdyCu(i,j))**2 + &
1212  ((g%CoriolisBu(i,j+1)-g%CoriolisBu(i,j)) * g%IdyCu(i,j+1))**2) ) ))
1213  enddo ; enddo
1214 
1215  do j=js,je ; do i=is-1,ieq
1216  cs%f2_dx2_u(i,j) = (g%dxCu(i,j)**2 + g%dyCu(i,j)**2) * &
1217  max(0.5* (g%CoriolisBu(i,j)**2+g%CoriolisBu(i,j-1)**2), absurdly_small_freq**2)
1218  cs%beta_dx2_u(i,j) = oneortwo * ((g%dxCu(i,j))**2 + (g%dyCu(i,j))**2) * (sqrt( &
1219  0.25*( (((g%CoriolisBu(i,j-1)-g%CoriolisBu(i-1,j-1)) * g%IdxCv(i,j-1))**2 + &
1220  ((g%CoriolisBu(i+1,j)-g%CoriolisBu(i,j)) * g%IdxCv(i+1,j))**2) + &
1221  (((g%CoriolisBu(i+1,j-1)-g%CoriolisBu(i,j-1)) * g%IdxCv(i+1,j-1))**2 + &
1222  ((g%CoriolisBu(i,j)-g%CoriolisBu(i-1,j)) * g%IdxCv(i,j))**2) ) + &
1223  ((g%CoriolisBu(i,j)-g%CoriolisBu(i,j-1)) * g%IdyCu(i,j))**2 ))
1224  enddo ; enddo
1225 
1226  do j=js-1,jeq ; do i=is,ie
1227  cs%f2_dx2_v(i,j) = ((g%dxCv(i,j))**2 + (g%dyCv(i,j))**2) * &
1228  max(0.5*(g%CoriolisBu(i,j)**2+g%CoriolisBu(i-1,j)**2), absurdly_small_freq**2)
1229  cs%beta_dx2_v(i,j) = oneortwo * ((g%dxCv(i,j))**2 + (g%dyCv(i,j))**2) * (sqrt( &
1230  ((g%CoriolisBu(i,j)-g%CoriolisBu(i-1,j)) * g%IdxCv(i,j))**2 + &
1231  0.25*( (((g%CoriolisBu(i,j)-g%CoriolisBu(i,j-1)) * g%IdyCu(i,j))**2 + &
1232  ((g%CoriolisBu(i-1,j+1)-g%CoriolisBu(i-1,j)) * g%IdyCu(i-1,j+1))**2) + &
1233  (((g%CoriolisBu(i,j+1)-g%CoriolisBu(i,j)) * g%IdyCu(i,j+1))**2 + &
1234  ((g%CoriolisBu(i-1,j)-g%CoriolisBu(i-1,j-1)) * g%IdyCu(i-1,j))**2) ) ))
1235  enddo ; enddo
1236 
1237  endif
1238 
1239  if (cs%Depth_scaled_KhTh) then
1240  cs%calculate_depth_fns = .true.
1241  allocate(cs%Depth_fn_u(isdb:iedb,jsd:jed)) ; cs%Depth_fn_u(:,:) = 0.0
1242  allocate(cs%Depth_fn_v(isd:ied,jsdb:jedb)) ; cs%Depth_fn_v(:,:) = 0.0
1243  call get_param(param_file, mdl, "DEPTH_SCALED_KHTH_H0", cs%depth_scaled_khth_h0, &
1244  "The depth above which KHTH is scaled away.",&
1245  units="m", default=1000.)
1246  call get_param(param_file, mdl, "DEPTH_SCALED_KHTH_EXP", cs%depth_scaled_khth_exp, &
1247  "The exponent used in the depth dependent scaling function for KHTH.",&
1248  units="nondim", default=3.0)
1249  endif
1250 
1251  ! Resolution %Rd_dx_h
1252  cs%id_Rd_dx = register_diag_field('ocean_model', 'Rd_dx', diag%axesT1, time, &
1253  'Ratio between deformation radius and grid spacing', 'm m-1')
1254  cs%calculate_Rd_dx = cs%calculate_Rd_dx .or. (cs%id_Rd_dx>0)
1255 
1256  if (cs%calculate_Rd_dx) then
1257  cs%calculate_cg1 = .true. ! We will need %cg1
1258  allocate(cs%Rd_dx_h(isd:ied,jsd:jed)) ; cs%Rd_dx_h(:,:) = 0.0
1259  allocate(cs%beta_dx2_h(isd:ied,jsd:jed)); cs%beta_dx2_h(:,:) = 0.0
1260  allocate(cs%f2_dx2_h(isd:ied,jsd:jed)) ; cs%f2_dx2_h(:,:) = 0.0
1261  do j=js-1,je+1 ; do i=is-1,ie+1
1262  cs%f2_dx2_h(i,j) = (g%dxT(i,j)**2 + g%dyT(i,j)**2) * &
1263  max(0.25 * ((g%CoriolisBu(i,j)**2 + g%CoriolisBu(i-1,j-1)**2) + &
1264  (g%CoriolisBu(i-1,j)**2 + g%CoriolisBu(i,j-1)**2)), &
1265  absurdly_small_freq**2)
1266  cs%beta_dx2_h(i,j) = oneortwo * ((g%dxT(i,j))**2 + (g%dyT(i,j))**2) * (sqrt(0.5 * &
1267  ( (((g%CoriolisBu(i,j)-g%CoriolisBu(i-1,j)) * g%IdxCv(i,j))**2 + &
1268  ((g%CoriolisBu(i,j-1)-g%CoriolisBu(i-1,j-1)) * g%IdxCv(i,j-1))**2) + &
1269  (((g%CoriolisBu(i,j)-g%CoriolisBu(i,j-1)) * g%IdyCu(i,j))**2 + &
1270  ((g%CoriolisBu(i-1,j)-g%CoriolisBu(i-1,j-1)) * g%IdyCu(i-1,j))**2) ) ))
1271  enddo ; enddo
1272  endif
1273 
1274  if (cs%calculate_cg1) then
1275  in_use = .true.
1276  allocate(cs%cg1(isd:ied,jsd:jed)) ; cs%cg1(:,:) = 0.0
1277  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
1278  "This sets the default value for the various _2018_ANSWERS parameters.", &
1279  default=.false.)
1280  call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", remap_answers_2018, &
1281  "If true, use the order of arithmetic and expressions that recover the "//&
1282  "answers from the end of 2018. Otherwise, use updated and more robust "//&
1283  "forms of the same expressions.", default=default_2018_answers)
1284  call get_param(param_file, mdl, "INTERNAL_WAVE_SPEED_TOL", wave_speed_tol, &
1285  "The fractional tolerance for finding the wave speeds.", &
1286  units="nondim", default=0.001)
1287  !### Set defaults so that wave_speed_min*wave_speed_tol >= 1e-9 m s-1
1288  call get_param(param_file, mdl, "INTERNAL_WAVE_SPEED_MIN", wave_speed_min, &
1289  "A floor in the first mode speed below which 0 used instead.", &
1290  units="m s-1", default=0.0, scale=us%m_s_to_L_T)
1291  call get_param(param_file, mdl, "INTERNAL_WAVE_SPEED_BETTER_EST", better_speed_est, &
1292  "If true, use a more robust estimate of the first mode wave speed as the "//&
1293  "starting point for iterations.", default=.false.) !### Change the default.
1294  call wave_speed_init(cs%wave_speed_CSp, use_ebt_mode=cs%Resoln_use_ebt, &
1295  mono_n2_depth=n2_filter_depth, remap_answers_2018=remap_answers_2018, &
1296  better_speed_est=better_speed_est, min_speed=wave_speed_min, &
1297  wave_speed_tol=wave_speed_tol)
1298  endif
1299 
1300  ! Leith parameters
1301  call get_param(param_file, mdl, "USE_QG_LEITH_GM", cs%use_QG_Leith_GM, &
1302  "If true, use the QG Leith viscosity as the GM coefficient.", &
1303  default=.false.)
1304 
1305  if (cs%Use_QG_Leith_GM) then
1306  call get_param(param_file, mdl, "LEITH_LAP_CONST", leith_lap_const, &
1307  "The nondimensional Laplacian Leith constant, \n"//&
1308  "often set to 1.0", units="nondim", default=0.0)
1309 
1310  call get_param(param_file, mdl, "USE_BETA_IN_LEITH", cs%use_beta_in_QG_Leith, &
1311  "If true, include the beta term in the Leith nonlinear eddy viscosity.", &
1312  default=.true.)
1313 
1314  alloc_(cs%Laplac3_const_u(isdb:iedb,jsd:jed)) ; cs%Laplac3_const_u(:,:) = 0.0
1315  alloc_(cs%Laplac3_const_v(isd:ied,jsdb:jedb)) ; cs%Laplac3_const_v(:,:) = 0.0
1316  alloc_(cs%KH_u_QG(isdb:iedb,jsd:jed,g%ke)) ; cs%KH_u_QG(:,:,:) = 0.0
1317  alloc_(cs%KH_v_QG(isd:ied,jsdb:jedb,g%ke)) ; cs%KH_v_QG(:,:,:) = 0.0
1318  ! register diagnostics
1319 
1320  cs%id_KH_u_QG = register_diag_field('ocean_model', 'KH_u_QG', diag%axesCuL, time, &
1321  'Horizontal viscosity from Leith QG, at u-points', 'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
1322  cs%id_KH_v_QG = register_diag_field('ocean_model', 'KH_v_QG', diag%axesCvL, time, &
1323  'Horizontal viscosity from Leith QG, at v-points', 'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
1324 
1325  do j=jsq,jeq+1 ; do i=is-1,ieq
1326  ! Static factors in the Leith schemes
1327  grid_sp_u2 = g%dyCu(i,j)*g%dxCu(i,j)
1328  grid_sp_u3 = grid_sp_u2*sqrt(grid_sp_u2)
1329  cs%Laplac3_const_u(i,j) = leith_lap_const * grid_sp_u3
1330  enddo ; enddo
1331  do j=js-1,jeq ; do i=isq,ieq+1
1332  ! Static factors in the Leith schemes
1333  grid_sp_v2 = g%dyCv(i,j)*g%dxCv(i,j)
1334  grid_sp_v3 = grid_sp_v2*sqrt(grid_sp_v2)
1335  cs%Laplac3_const_v(i,j) = leith_lap_const * grid_sp_v3
1336  enddo ; enddo
1337 
1338  if (.not. cs%use_stored_slopes) call mom_error(fatal, &
1339  "MOM_lateral_mixing_coeffs.F90, VarMix_init:"//&
1340  "USE_STORED_SLOPES must be True when using QG Leith.")
1341  endif
1342 
1343  ! If nothing is being stored in this class then deallocate
1344  if (in_use) then
1345  cs%use_variable_mixing = .true.
1346  else
1347  deallocate(cs)
1348  return
1349  endif
1350 
1351 end subroutine varmix_init
1352 
1353 !> \namespace mom_lateral_mixing_coeffs
1354 !!
1355 !! This module provides a container for various factors used in prescribing diffusivities, that are
1356 !! a function of the state (in particular the stratification and isoneutral slopes).
1357 !!
1358 !! \section section_Resolution_Function The resolution function
1359 !!
1360 !! The resolution function is expressed in terms of the ratio of grid-spacing to deformation radius.
1361 !! The square of the resolution parameter is
1362 !!
1363 !! \f[
1364 !! R^2 = \frac{L_d^2}{\Delta^2} = \frac{ c_g^2 }{ f^2 \Delta^2 + c_g \beta \Delta^2 }
1365 !! \f]
1366 !!
1367 !! where the grid spacing is calculated as
1368 !!
1369 !! \f[
1370 !! \Delta^2 = \Delta x^2 + \Delta y^2 .
1371 !! \f]
1372 !!
1373 !! \todo Check this reference to Bob on/off paper.
1374 !! The resolution function used in scaling diffusivities (Hallberg, 2010) is
1375 !!
1376 !! \f[
1377 !! r(\Delta,L_d) = \frac{1}{1+(\alpha R)^p}
1378 !! \f]
1379 !!
1380 !! The resolution function can be applied independently to thickness diffusion (module mom_thickness_diffuse),
1381 !! tracer diffusion (mom_tracer_hordiff) lateral viscosity (mom_hor_visc).
1382 !!
1383 !! Robert Hallberg, 2013: Using a resolution function to regulate parameterizations of oceanic mesoscale eddy effects.
1384 !! Ocean Modelling, 71, pp 92-103. http://dx.doi.org/10.1016/j.ocemod.2013.08.007
1385 !!
1386 !! | Symbol | Module parameter |
1387 !! | ------ | --------------- |
1388 !! | - | <code>USE_VARIABLE_MIXING</code> |
1389 !! | - | <code>RESOLN_SCALED_KH</code> |
1390 !! | - | <code>RESOLN_SCALED_KHTH</code> |
1391 !! | - | <code>RESOLN_SCALED_KHTR</code> |
1392 !! | \f$ \alpha \f$ | <code>KH_RES_SCALE_COEF</code> (for thickness and tracer diffusivity) |
1393 !! | \f$ p \f$ | <code>KH_RES_FN_POWER</code> (for thickness and tracer diffusivity) |
1394 !! | \f$ \alpha \f$ | <code>VISC_RES_SCALE_COEF</code> (for lateral viscosity) |
1395 !! | \f$ p \f$ | <code>VISC_RES_FN_POWER</code> (for lateral viscosity) |
1396 !! | - | <code>GILL_EQUATORIAL_LD</code> |
1397 !!
1398 !!
1399 !!
1400 !! \section section_Vicbeck Visbeck diffusivity
1401 !!
1402 !! This module also calculates factors used in setting the thickness diffusivity similar to a Visbeck et al., 1997,
1403 !! scheme. The factors are combined in mom_thickness_diffuse::thickness_diffuse() but calculated in this module.
1404 !!
1405 !! \f[
1406 !! \kappa_h = \alpha_s L_s^2 S N
1407 !! \f]
1408 !!
1409 !! where \f$S\f$ is the magnitude of the isoneutral slope and \f$N\f$ is the Brunt-Vaisala frequency.
1410 !!
1411 !! Visbeck, Marshall, Haine and Spall, 1997: Specification of Eddy Transfer Coefficients in Coarse-Resolution
1412 !! Ocean Circulation Models. J. Phys. Oceanogr. http://dx.doi.org/10.1175/1520-0485(1997)027%3C0381:SOETCI%3E2.0.CO;2
1413 !!
1414 !! | Symbol | Module parameter |
1415 !! | ------ | --------------- |
1416 !! | - | <code>USE_VARIABLE_MIXING</code> |
1417 !! | \f$ \alpha_s \f$ | <code>KHTH_SLOPE_CFF</code> (for mom_thickness_diffuse module)|
1418 !! | \f$ \alpha_s \f$ | <code>KHTR_SLOPE_CFF</code> (for mom_tracer_hordiff module)|
1419 !! | \f$ L_{s} \f$ | <code>VISBECK_L_SCALE</code> |
1420 !! | \f$ S_{max} \f$ | <code>VISBECK_MAX_SLOPE</code> |
1421 !!
1422 !!
1423 !! \section section_vertical_structure_khth Vertical structure function for KhTh
1424 !!
1425 !! The thickness diffusivity can be prescribed a vertical distribution with the shape of the equivalent barotropic
1426 !! velocity mode. The structure function is stored in the control structure for thie module (varmix_cs) but is
1427 !! calculated using subroutines in mom_wave_speed.
1428 !!
1429 !! | Symbol | Module parameter |
1430 !! | ------ | --------------- |
1431 !! | - | <code>KHTH_USE_EBT_STRUCT</code> |
1432 
1433 end module mom_lateral_mixing_coeffs
Calculates the heights of sruface or all interfaces from layer thicknesses.
Control structure for MOM_wave_speed.
Functions for calculating interface heights, including free surface height.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Routines for calculating baroclinic wave speeds.
The MOM6 facility to parse input files for runtime parameters.
Variable mixing coefficients.
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:59
Make a diagnostic available for averaging or output.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Routines for error handling and I/O management.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Set up a group of halo updates.
Definition: MOM_domains.F90:84
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Controls where open boundary conditions are applied.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Do a halo update on an array.
Definition: MOM_domains.F90:54
Calculations of isoneutral slopes and stratification.
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.