MOM6
MOM_set_diffusivity.F90
1 !> Calculate vertical diffusivity from all mixing processes
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_bkgnd_mixing, only : calculate_bkgnd_mixing, bkgnd_mixing_init, bkgnd_mixing_cs
7 use mom_bkgnd_mixing, only : bkgnd_mixing_end
8 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
9 use mom_cpu_clock, only : clock_module_driver, clock_module, clock_routine
10 use mom_cvmix_ddiff, only : cvmix_ddiff_init, cvmix_ddiff_end, cvmix_ddiff_cs
11 use mom_cvmix_ddiff, only : compute_ddiff_coeffs
12 use mom_cvmix_shear, only : calculate_cvmix_shear, cvmix_shear_init, cvmix_shear_cs
13 use mom_cvmix_shear, only : cvmix_shear_end
14 use mom_diag_mediator, only : diag_ctrl, time_type
15 use mom_diag_mediator, only : post_data, register_diag_field
16 use mom_debugging, only : hchksum, uvchksum, bchksum, hchksum_pair
17 use mom_eos, only : calculate_density, calculate_density_derivs, eos_domain
18 use mom_error_handler, only : mom_error, is_root_pe, fatal, warning, note
19 use mom_error_handler, only : calltree_showquery
20 use mom_error_handler, only : calltree_enter, calltree_leave, calltree_waypoint
22 use mom_forcing_type, only : forcing, optics_type
23 use mom_full_convection, only : full_convection
24 use mom_grid, only : ocean_grid_type
25 use mom_internal_tides, only : int_tide_cs, get_lowmode_loss
26 use mom_intrinsic_functions, only : invcosh
27 use mom_io, only : slasher, mom_read_data
28 use mom_isopycnal_slopes, only : vert_fill_ts
29 use mom_kappa_shear, only : calculate_kappa_shear, kappa_shear_init, kappa_shear_cs
30 use mom_kappa_shear, only : calc_kappa_shear_vertex, kappa_shear_at_vertex
32 use mom_open_boundary, only : obc_direction_e, obc_direction_w, obc_direction_n, obc_direction_s
33 use mom_string_functions, only : uppercase
34 use mom_tidal_mixing, only : tidal_mixing_cs, calculate_tidal_mixing, tidal_mixing_h_amp
35 use mom_tidal_mixing, only : setup_tidal_diagnostics, post_tidal_diagnostics
36 use mom_tidal_mixing, only : tidal_mixing_init, tidal_mixing_end
40 use user_change_diffusivity, only : user_change_diff, user_change_diff_init
41 use user_change_diffusivity, only : user_change_diff_end, user_change_diff_cs
42 
43 implicit none ; private
44 
45 #include <MOM_memory.h>
46 
47 public set_diffusivity
48 public set_bbl_tke
49 public set_diffusivity_init
50 public set_diffusivity_end
51 
52 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
53 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
54 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
55 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
56 
57 !> This control structure contains parameters for MOM_set_diffusivity.
58 type, public :: set_diffusivity_cs ; private
59  logical :: debug !< If true, write verbose checksums for debugging.
60 
61  logical :: bulkmixedlayer !< If true, a refined bulk mixed layer is used with
62  !! GV%nk_rho_varies variable density mixed & buffer layers.
63  real :: fluxri_max !< The flux Richardson number where the stratification is
64  !! large enough that N2 > omega2 [nondim]. The full expression
65  !! for the Flux Richardson number is usually
66  !! FLUX_RI_MAX*N2/(N2+OMEGA2). The default is 0.2.
67  logical :: bottomdraglaw !< If true, the bottom stress is calculated with a
68  !! drag law c_drag*|u|*u.
69  logical :: bbl_mixing_as_max !< If true, take the maximum of the diffusivity
70  !! from the BBL mixing and the other diffusivities.
71  !! Otherwise, diffusivities from the BBL_mixing is
72  !! added.
73  logical :: use_lotw_bbl_diffusivity !< If true, use simpler/less precise, BBL diffusivity.
74  logical :: lotw_bbl_use_omega !< If true, use simpler/less precise, BBL diffusivity.
75  real :: bbl_effic !< efficiency with which the energy extracted
76  !! by bottom drag drives BBL diffusion [nondim]
77  real :: cdrag !< quadratic drag coefficient [nondim]
78  real :: imax_decay !< inverse of a maximum decay scale for
79  !! bottom-drag driven turbulence [Z-1 ~> m-1].
80  real :: kv !< The interior vertical viscosity [Z2 T-1 ~> m2 s-1].
81  real :: kd !< interior diapycnal diffusivity [Z2 T-1 ~> m2 s-1].
82  real :: kd_min !< minimum diapycnal diffusivity [Z2 T-1 ~> m2 s-1].
83  real :: kd_max !< maximum increment for diapycnal diffusivity [Z2 T-1 ~> m2 s-1].
84  !! Set to a negative value to have no limit.
85  real :: kd_add !< uniform diffusivity added everywhere without
86  !! filtering or scaling [Z2 T-1 ~> m2 s-1].
87  real :: kd_smooth !< Vertical diffusivity used to interpolate more
88  !! sensible values of T & S into thin layers [Z2 T-1 ~> m2 s-1].
89  type(diag_ctrl), pointer :: diag => null() !< structure to regulate diagnostic output timing
90 
91  logical :: limit_dissipation !< If enabled, dissipation is limited to be larger
92  !! than the following:
93  real :: dissip_min !< Minimum dissipation [R Z2 T-3 ~> W m-3]
94  real :: dissip_n0 !< Coefficient a in minimum dissipation = a+b*N [R Z2 T-3 ~> W m-3]
95  real :: dissip_n1 !< Coefficient b in minimum dissipation = a+b*N [R Z2 T-2 ~> J m-3]
96  real :: dissip_n2 !< Coefficient c in minimum dissipation = c*N2 [R Z2 T-1 ~> J s m-3]
97  real :: dissip_kd_min !< Minimum Kd [Z2 T-1 ~> m2 s-1], with dissipation Rho0*Kd_min*N^2
98 
99  real :: omega !< Earth's rotation frequency [T-1 ~> s-1]
100  logical :: ml_radiation !< allow a fraction of TKE available from wind work
101  !! to penetrate below mixed layer base with a vertical
102  !! decay scale determined by the minimum of
103  !! (1) The depth of the mixed layer, or
104  !! (2) An Ekman length scale.
105  !! Energy available to drive mixing below the mixed layer is
106  !! given by E = ML_RAD_COEFF*MSTAR*USTAR**3. Optionally, if
107  !! ML_rad_TKE_decay is true, this is further reduced by a factor
108  !! of exp(-h_ML*Idecay_len_TkE), where Idecay_len_TKE is
109  !! calculated the same way as in the mixed layer code.
110  !! The diapycnal diffusivity is KD(k) = E/(N2(k)+OMEGA2),
111  !! where N2 is the squared buoyancy frequency [T-2 ~> s-2] and OMEGA2
112  !! is the rotation rate of the earth squared.
113  real :: ml_rad_kd_max !< Maximum diapycnal diffusivity due to turbulence
114  !! radiated from the base of the mixed layer [Z2 T-1 ~> m2 s-1].
115  real :: ml_rad_efold_coeff !< non-dim coefficient to scale penetration depth
116  real :: ml_rad_coeff !< coefficient, which scales MSTAR*USTAR^3 to
117  !! obtain energy available for mixing below
118  !! mixed layer base [nondim]
119  logical :: ml_rad_bug !< If true use code with a bug that reduces the energy available
120  !! in the transition layer by a factor of the inverse of the energy
121  !! deposition lenthscale (in m).
122  logical :: ml_rad_tke_decay !< If true, apply same exponential decay
123  !! to ML_rad as applied to the other surface
124  !! sources of TKE in the mixed layer code.
125  real :: ustar_min !< A minimum value of ustar to avoid numerical
126  !! problems [Z T-1 ~> m s-1]. If the value is small enough,
127  !! this parameter should not affect the solution.
128  real :: tke_decay !< ratio of natural Ekman depth to TKE decay scale [nondim]
129  real :: mstar !< ratio of friction velocity cubed to
130  !! TKE input to the mixed layer [nondim]
131  logical :: ml_use_omega !< If true, use absolute rotation rate instead
132  !! of the vertical component of rotation when
133  !! setting the decay scale for mixed layer turbulence.
134  real :: ml_omega_frac !< When setting the decay scale for turbulence, use
135  !! this fraction of the absolute rotation rate blended
136  !! with the local value of f, as f^2 ~= (1-of)*f^2 + of*4*omega^2.
137  logical :: user_change_diff !< If true, call user-defined code to change diffusivity.
138  logical :: usekappashear !< If true, use the kappa_shear module to find the
139  !! shear-driven diapycnal diffusivity.
140  logical :: vertex_shear !< If true, do the calculations of the shear-driven mixing
141  !! at the cell vertices (i.e., the vorticity points).
142  logical :: use_cvmix_shear !< If true, use one of the CVMix modules to find
143  !! shear-driven diapycnal diffusivity.
144  logical :: double_diffusion !< If true, enable double-diffusive mixing using an old method.
145  logical :: use_cvmix_ddiff !< If true, enable double-diffusive mixing via CVMix.
146  logical :: use_tidal_mixing !< If true, activate tidal mixing diffusivity.
147  logical :: simple_tke_to_kd !< If true, uses a simple estimate of Kd/TKE that
148  !! does not rely on a layer-formulation.
149  real :: max_rrho_salt_fingers !< max density ratio for salt fingering
150  real :: max_salt_diff_salt_fingers !< max salt diffusivity for salt fingers [Z2 T-1 ~> m2 s-1]
151  real :: kv_molecular !< molecular visc for double diff convect [Z2 T-1 ~> m2 s-1]
152 
153  logical :: answers_2018 !< If true, use the order of arithmetic and expressions that recover the
154  !! answers from the end of 2018. Otherwise, use updated and more robust
155  !! forms of the same expressions.
156 
157  character(len=200) :: inputdir !< The directory in which input files are found
158  type(user_change_diff_cs), pointer :: user_change_diff_csp => null() !< Control structure for a child module
159  type(kappa_shear_cs), pointer :: kappashear_csp => null() !< Control structure for a child module
160  type(cvmix_shear_cs), pointer :: cvmix_shear_csp => null() !< Control structure for a child module
161  type(cvmix_ddiff_cs), pointer :: cvmix_ddiff_csp => null() !< Control structure for a child module
162  type(bkgnd_mixing_cs), pointer :: bkgnd_mixing_csp => null() !< Control structure for a child module
163  type(int_tide_cs), pointer :: int_tide_csp => null() !< Control structure for a child module
164  type(tidal_mixing_cs), pointer :: tidal_mixing_csp => null() !< Control structure for a child module
165 
166  !>@{ Diagnostic IDs
167  integer :: id_maxtke = -1, id_tke_to_kd = -1, id_kd_user = -1
168  integer :: id_kd_layer = -1, id_kd_bbl = -1, id_n2 = -1
169  integer :: id_kd_work = -1, id_kt_extra = -1, id_ks_extra = -1, id_r_rho = -1
170  integer :: id_kd_bkgnd = -1, id_kv_bkgnd = -1
171  !>@}
172 
173 end type set_diffusivity_cs
174 
175 !> This structure has memory for used in calculating diagnostics of diffusivity
177  real, pointer, dimension(:,:,:) :: &
178  n2_3d => null(), & !< squared buoyancy frequency at interfaces [T-2 ~> s-2]
179  kd_user => null(), & !< user-added diffusivity at interfaces [Z2 T-1 ~> m2 s-1]
180  kd_bbl => null(), & !< BBL diffusivity at interfaces [Z2 T-1 ~> m2 s-1]
181  kd_work => null(), & !< layer integrated work by diapycnal mixing [R Z3 T-3 ~> W m-2]
182  maxtke => null(), & !< energy required to entrain to h_max [Z3 T-3 ~> m3 s-3]
183  kd_bkgnd => null(), & !< Background diffusivity at interfaces [Z2 T-1 ~> m2 s-1]
184  kv_bkgnd => null(), & !< Viscosity from ackground diffusivity at interfaces [Z2 T-1 ~> m2 s-1]
185  kt_extra => null(), & !< double diffusion diffusivity for temp [Z2 T-1 ~> m2 s-1].
186  ks_extra => null(), & !< double diffusion diffusivity for saln [Z2 T-1 ~> m2 s-1].
187  drho_rat => null() !< The density difference ratio used in double diffusion [nondim].
188  real, pointer, dimension(:,:,:) :: tke_to_kd => null()
189  !< conversion rate (~1.0 / (G_Earth + dRho_lay)) between TKE
190  !! dissipated within a layer and Kd in that layer
191  !! [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
192 
193 end type diffusivity_diags
194 
195 !>@{ CPU time clocks
196 integer :: id_clock_kappashear, id_clock_cvmix_ddiff
197 !>@}
198 
199 contains
200 
201 !> Sets the interior vertical diffusion of scalars due to the following processes:
202 !! 1. Shear-driven mixing: two options, Jackson et at. and KPP interior;
203 !! 2. Background mixing via CVMix (Bryan-Lewis profile) or the scheme described by
204 !! Harrison & Hallberg, JPO 2008;
205 !! 3. Double-diffusion, old method and new method via CVMix;
206 !! 4. Tidal mixing: many options available, see MOM_tidal_mixing.F90;
207 !! In addition, this subroutine has the option to set the interior vertical
208 !! viscosity associated with processes 1,2 and 4 listed above, which is stored in
209 !! visc%Kv_slow. Vertical viscosity due to shear-driven mixing is passed via
210 !! visc%Kv_shear
211 subroutine set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, optics, visc, dt, &
212  G, GV, US, CS, Kd_lay, Kd_int, Kd_extra_T, Kd_extra_S)
213  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
214  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
215  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
216  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
217  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
218  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
219  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
220  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
221  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
222  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
223  intent(in) :: u_h !< Zonal velocity interpolated to h points [L T-1 ~> m s-1].
224  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
225  intent(in) :: v_h !< Meridional velocity interpolated to h points [L T-1 ~> m s-1].
226  type(thermo_var_ptrs), intent(inout) :: tv !< Structure with pointers to thermodynamic
227  !! fields. Out is for tv%TempxPmE.
228  type(forcing), intent(in) :: fluxes !< A structure of thermodynamic surface fluxes
229  type(optics_type), pointer :: optics !< A structure describing the optical
230  !! properties of the ocean.
231  type(vertvisc_type), intent(inout) :: visc !< Structure containing vertical viscosities, bottom
232  !! boundary layer properies, and related fields.
233  real, intent(in) :: dt !< Time increment [T ~> s].
234  type(set_diffusivity_cs), pointer :: CS !< Module control structure.
235  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
236  optional, intent(out) :: Kd_lay !< Diapycnal diffusivity of each layer [Z2 T-1 ~> m2 s-1].
237  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
238  optional, intent(out) :: Kd_int !< Diapycnal diffusivity at each interface [Z2 T-1 ~> m2 s-1].
239  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
240  optional, intent(out) :: Kd_extra_T !< The extra diffusivity at interfaces of
241  !! temperature due to double diffusion relative to
242  !! the diffusivity of density [Z2 T-1 ~> m2 s-1].
243  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
244  optional, intent(out) :: Kd_extra_S !< The extra diffusivity at interfaces of
245  !! salinity due to double diffusion relative to
246  !! the diffusivity of density [Z2 T-1 ~> m2 s-1].
247 
248  ! local variables
249  real, dimension(SZI_(G)) :: &
250  N2_bot ! bottom squared buoyancy frequency [T-2 ~> s-2]
251 
252  type(diffusivity_diags) :: dd ! structure with arrays of available diags
253 
254  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
255  T_f, S_f ! Temperature and salinity [degC] and [ppt] with properties in massless layers
256  ! filled vertically by diffusion or the properties after full convective adjustment.
257 
258  real, dimension(SZI_(G),SZK_(G)) :: &
259  N2_lay, & !< Squared buoyancy frequency associated with layers [T-2 ~> s-2]
260  Kd_lay_2d, & !< The layer diffusivities [Z2 T-1 ~> m2 s-1]
261  maxTKE, & !< Energy required to entrain to h_max [Z3 T-3 ~> m3 s-3]
262  TKE_to_Kd !< Conversion rate (~1.0 / (G_Earth + dRho_lay)) between
263  !< TKE dissipated within a layer and Kd in that layer
264  !< [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
265 
266  real, dimension(SZI_(G),SZK_(G)+1) :: &
267  N2_int, & !< squared buoyancy frequency associated at interfaces [T-2 ~> s-2]
268  Kd_int_2d, & !< The interface diffusivities [Z2 T-1 ~> m2 s-1]
269  Kv_bkgnd, & !< The background diffusion related interface viscosities [Z2 T-1 ~> m2 s-1]
270  dRho_int, & !< Locally referenced potential density difference across interfaces [R ~> kg m-3]
271  KT_extra, & !< Double difusion diffusivity of temperature [Z2 T-1 ~> m2 s-1]
272  KS_extra !< Double difusion diffusivity of salinity [Z2 T-1 ~> m2 s-1]
273 
274  real :: dissip ! local variable for dissipation calculations [Z2 R T-3 ~> W m-3]
275  real :: Omega2 ! squared absolute rotation rate [T-2 ~> s-2]
276 
277  logical :: use_EOS ! If true, compute density from T/S using equation of state.
278  logical :: TKE_to_Kd_used ! If true, TKE_to_Kd and maxTKE need to be calculated.
279  integer :: kb(szi_(g)) ! The index of the lightest layer denser than the
280  ! buffer layer, or -1 without a bulk mixed layer.
281  logical :: showCallTree ! If true, show the call tree.
282 
283  integer :: i, j, k, is, ie, js, je, nz, isd, ied, jsd, jed
284 
285  real :: kappa_dt_fill ! diffusivity times a timestep used to fill massless layers [Z2 ~> m2]
286 
287  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
288  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
289  showcalltree = calltree_showquery()
290  if (showcalltree) call calltree_enter("set_diffusivity(), MOM_set_diffusivity.F90")
291 
292  if (.not.associated(cs)) call mom_error(fatal,"set_diffusivity: "//&
293  "Module must be initialized before it is used.")
294 
295  if (cs%answers_2018) then
296  ! These hard-coded dimensional parameters are being replaced.
297  kappa_dt_fill = us%m_to_Z**2 * 1.e-3 * 7200.
298  else
299  kappa_dt_fill = cs%Kd_smooth * dt
300  endif
301  omega2 = cs%omega * cs%omega
302 
303  use_eos = associated(tv%eqn_of_state)
304 
305  if ((cs%use_CVMix_ddiff .or. cs%double_diffusion) .and. &
306  .not.(present(kd_extra_t) .and. present(kd_extra_s))) &
307  call mom_error(fatal, "set_diffusivity: both Kd_extra_T and Kd_extra_S must be present "//&
308  "when USE_CVMIX_DDIFF or DOUBLE_DIFFUSION are true.")
309 
310  tke_to_kd_used = (cs%use_tidal_mixing .or. cs%ML_radiation .or. &
311  (cs%bottomdraglaw .and. .not.cs%use_LOTW_BBL_diffusivity))
312 
313  ! Set Kd_lay, Kd_int and Kv_slow to constant values, mostly to fill the halos.
314  if (present(kd_lay)) kd_lay(:,:,:) = cs%Kd
315  if (present(kd_int)) kd_int(:,:,:) = cs%Kd
316  if (present(kd_extra_t)) kd_extra_t(:,:,:) = 0.0
317  if (present(kd_extra_s)) kd_extra_s(:,:,:) = 0.0
318  if (associated(visc%Kv_slow)) visc%Kv_slow(:,:,:) = cs%Kv
319 
320  ! Set up arrays for diagnostics.
321 
322  if (cs%id_N2 > 0) then
323  allocate(dd%N2_3d(isd:ied,jsd:jed,nz+1)) ; dd%N2_3d(:,:,:) = 0.0
324  endif
325  if (cs%id_Kd_user > 0) then
326  allocate(dd%Kd_user(isd:ied,jsd:jed,nz+1)) ; dd%Kd_user(:,:,:) = 0.0
327  endif
328  if (cs%id_Kd_work > 0) then
329  allocate(dd%Kd_work(isd:ied,jsd:jed,nz)) ; dd%Kd_work(:,:,:) = 0.0
330  endif
331  if (cs%id_maxTKE > 0) then
332  allocate(dd%maxTKE(isd:ied,jsd:jed,nz)) ; dd%maxTKE(:,:,:) = 0.0
333  endif
334  if (cs%id_TKE_to_Kd > 0) then
335  allocate(dd%TKE_to_Kd(isd:ied,jsd:jed,nz)) ; dd%TKE_to_Kd(:,:,:) = 0.0
336  endif
337  if ((cs%double_diffusion) .and. (cs%id_KT_extra > 0)) then
338  allocate(dd%KT_extra(isd:ied,jsd:jed,nz+1)) ; dd%KT_extra(:,:,:) = 0.0
339  endif
340  if ((cs%double_diffusion) .and. (cs%id_KS_extra > 0)) then
341  allocate(dd%KS_extra(isd:ied,jsd:jed,nz+1)) ; dd%KS_extra(:,:,:) = 0.0
342  endif
343  if (cs%id_R_rho > 0) then
344  allocate(dd%drho_rat(isd:ied,jsd:jed,nz+1)) ; dd%drho_rat(:,:,:) = 0.0
345  endif
346  if (cs%id_Kd_BBL > 0) then
347  allocate(dd%Kd_BBL(isd:ied,jsd:jed,nz+1)) ; dd%Kd_BBL(:,:,:) = 0.0
348  endif
349 
350  if (cs%id_Kd_bkgnd > 0) then
351  allocate(dd%Kd_bkgnd(isd:ied,jsd:jed,nz+1)) ; dd%Kd_bkgnd(:,:,:) = 0.
352  endif
353  if (cs%id_Kv_bkgnd > 0) then
354  allocate(dd%Kv_bkgnd(isd:ied,jsd:jed,nz+1)) ; dd%Kv_bkgnd(:,:,:) = 0.
355  endif
356 
357  ! set up arrays for tidal mixing diagnostics
358  call setup_tidal_diagnostics(g, cs%tidal_mixing_CSp)
359 
360  if (cs%useKappaShear) then
361  if (cs%debug) then
362  call hchksum_pair("before calc_KS [uv]_h", u_h, v_h, g%HI, scale=us%L_T_to_m_s)
363  endif
364  call cpu_clock_begin(id_clock_kappashear)
365  if (cs%Vertex_shear) then
366  call full_convection(g, gv, us, h, tv, t_f, s_f, fluxes%p_surf, &
367  (gv%Z_to_H**2)*kappa_dt_fill, halo=1)
368 
369  call calc_kappa_shear_vertex(u, v, h, t_f, s_f, tv, fluxes%p_surf, visc%Kd_shear, &
370  visc%TKE_turb, visc%Kv_shear_Bu, dt, g, gv, us, cs%kappaShear_CSp)
371  if (associated(visc%Kv_shear)) visc%Kv_shear(:,:,:) = 0.0 ! needed for other parameterizations
372  if (cs%debug) then
373  call hchksum(visc%Kd_shear, "after calc_KS_vert visc%Kd_shear", g%HI, scale=us%Z2_T_to_m2_s)
374  call bchksum(visc%Kv_shear_Bu, "after calc_KS_vert visc%Kv_shear_Bu", g%HI, scale=us%Z2_T_to_m2_s)
375  call bchksum(visc%TKE_turb, "after calc_KS_vert visc%TKE_turb", g%HI, scale=us%Z_to_m**2*us%s_to_T**2)
376  endif
377  else
378  ! Changes: visc%Kd_shear ; Sets: visc%Kv_shear and visc%TKE_turb
379  call calculate_kappa_shear(u_h, v_h, h, tv, fluxes%p_surf, visc%Kd_shear, visc%TKE_turb, &
380  visc%Kv_shear, dt, g, gv, us, cs%kappaShear_CSp)
381  if (cs%debug) then
382  call hchksum(visc%Kd_shear, "after calc_KS visc%Kd_shear", g%HI, scale=us%Z2_T_to_m2_s)
383  call hchksum(visc%Kv_shear, "after calc_KS visc%Kv_shear", g%HI, scale=us%Z2_T_to_m2_s)
384  call hchksum(visc%TKE_turb, "after calc_KS visc%TKE_turb", g%HI, scale=us%Z_to_m**2*us%s_to_T**2)
385  endif
386  endif
387  call cpu_clock_end(id_clock_kappashear)
388  if (showcalltree) call calltree_waypoint("done with calculate_kappa_shear (set_diffusivity)")
389  elseif (cs%use_CVMix_shear) then
390  !NOTE{BGR}: this needs to be cleaned up. It works in 1D case, but has not been tested outside.
391  call calculate_cvmix_shear(u_h, v_h, h, tv, visc%Kd_shear, visc%Kv_shear, g, gv, us, cs%CVMix_shear_CSp)
392  if (cs%debug) then
393  call hchksum(visc%Kd_shear, "after CVMix_shear visc%Kd_shear", g%HI, scale=us%Z2_T_to_m2_s)
394  call hchksum(visc%Kv_shear, "after CVMix_shear visc%Kv_shear", g%HI, scale=us%Z2_T_to_m2_s)
395  endif
396  elseif (associated(visc%Kv_shear)) then
397  visc%Kv_shear(:,:,:) = 0.0 ! needed if calculate_kappa_shear is not enabled
398  endif
399 
400  ! Smooth the properties through massless layers.
401  if (use_eos) then
402  if (cs%debug) then
403  call hchksum(tv%T, "before vert_fill_TS tv%T",g%HI)
404  call hchksum(tv%S, "before vert_fill_TS tv%S",g%HI)
405  call hchksum(h, "before vert_fill_TS h",g%HI, scale=gv%H_to_m)
406  endif
407  call vert_fill_ts(h, tv%T, tv%S, kappa_dt_fill, t_f, s_f, g, gv, larger_h_denom=.true.)
408  if (cs%debug) then
409  call hchksum(tv%T, "after vert_fill_TS tv%T",g%HI)
410  call hchksum(tv%S, "after vert_fill_TS tv%S",g%HI)
411  call hchksum(h, "after vert_fill_TS h",g%HI, scale=gv%H_to_m)
412  endif
413  endif
414 
415  ! Calculate the diffusivities, Kd_lay and Kd_int, for each layer and interface. This would
416  ! be an appropriate place to add a depth-dependent parameterization or another explicit
417  ! parameterization of Kd.
418 
419  !$OMP parallel do default(shared) private(dRho_int,N2_lay,Kd_lay_2d,Kd_int_2d,Kv_bkgnd,N2_int,&
420  !$OMP N2_bot,KT_extra,KS_extra,TKE_to_Kd,maxTKE,dissip,kb)
421  do j=js,je
422 
423  ! Set up variables related to the stratification.
424  call find_n2(h, tv, t_f, s_f, fluxes, j, g, gv, us, cs, drho_int, n2_lay, n2_int, n2_bot)
425 
426  if (associated(dd%N2_3d)) then
427  do k=1,nz+1 ; do i=is,ie ; dd%N2_3d(i,j,k) = n2_int(i,k) ; enddo ; enddo
428  endif
429 
430  ! Add background mixing
431  call calculate_bkgnd_mixing(h, tv, n2_lay, kd_lay_2d, kd_int_2d, kv_bkgnd, j, g, gv, us, cs%bkgnd_mixing_csp)
432  ! Update Kv and 3-d diffusivity diagnostics.
433  if (associated(visc%Kv_slow)) then ; do k=1,nz+1 ; do i=is,ie
434  visc%Kv_slow(i,j,k) = visc%Kv_slow(i,j,k) + kv_bkgnd(i,k)
435  enddo ; enddo ; endif
436  if (cs%id_Kv_bkgnd > 0) then ; do k=1,nz+1 ; do i=is,ie
437  dd%Kv_bkgnd(i,j,k) = kv_bkgnd(i,k)
438  enddo ; enddo ; endif
439  if (cs%id_Kd_bkgnd > 0) then ; do k=1,nz+1 ; do i=is,ie
440  dd%Kd_bkgnd(i,j,k) = kd_int_2d(i,k)
441  enddo ; enddo ; endif
442 
443  ! Double-diffusion (old method)
444  if (cs%double_diffusion) then
445  call double_diffusion(tv, h, t_f, s_f, j, g, gv, us, cs, kt_extra, ks_extra)
446  ! One of Kd_extra_T and Kd_extra_S is always 0. Kd_extra_S is positive for salt fingering.
447  ! Kd_extra_T is positive for double diffusive convection.
448  do k=2,nz ; do i=is,ie
449  if (ks_extra(i,k) > kt_extra(i,k)) then ! salt fingering
450  kd_lay_2d(i,k-1) = kd_lay_2d(i,k-1) + 0.5 * kt_extra(i,k)
451  kd_lay_2d(i,k) = kd_lay_2d(i,k) + 0.5 * kt_extra(i,k)
452  kd_extra_s(i,j,k) = (ks_extra(i,k) - kt_extra(i,k))
453  kd_extra_t(i,j,k) = 0.0
454  elseif (kt_extra(i,k) > 0.0) then ! double-diffusive convection
455  kd_lay_2d(i,k-1) = kd_lay_2d(i,k-1) + 0.5 * ks_extra(i,k)
456  kd_lay_2d(i,k) = kd_lay_2d(i,k) + 0.5 * ks_extra(i,k)
457  kd_extra_t(i,j,k) = (kt_extra(i,k) - ks_extra(i,k))
458  kd_extra_s(i,j,k) = 0.0
459  else ! There is no double diffusion at this interface.
460  kd_extra_t(i,j,k) = 0.0
461  kd_extra_s(i,j,k) = 0.0
462  endif
463  enddo ; enddo
464  if (associated(dd%KT_extra)) then ; do k=1,nz+1 ; do i=is,ie
465  dd%KT_extra(i,j,k) = kt_extra(i,k)
466  enddo ; enddo ; endif
467 
468  if (associated(dd%KS_extra)) then ; do k=1,nz+1 ; do i=is,ie
469  dd%KS_extra(i,j,k) = ks_extra(i,k)
470  enddo ; enddo ; endif
471  endif
472 
473  ! Apply double diffusion via CVMix
474  ! GMM, we need to pass HBL to compute_ddiff_coeffs, but it is not yet available.
475  if (cs%use_CVMix_ddiff) then
476  call cpu_clock_begin(id_clock_cvmix_ddiff)
477  if (associated(dd%drho_rat)) then
478  call compute_ddiff_coeffs(h, tv, g, gv, us, j, kd_extra_t, kd_extra_s, &
479  cs%CVMix_ddiff_csp, dd%drho_rat)
480  else
481  call compute_ddiff_coeffs(h, tv, g, gv, us, j, kd_extra_t, kd_extra_s, cs%CVMix_ddiff_csp)
482  endif
483  call cpu_clock_end(id_clock_cvmix_ddiff)
484  endif
485 
486  ! Calculate conversion ratios from TKE to layer diffusivities.
487  if (tke_to_kd_used) then
488  call find_tke_to_kd(h, tv, drho_int, n2_lay, j, dt, g, gv, us, cs, tke_to_kd, maxtke, kb)
489  if (associated(dd%maxTKE)) then ; do k=1,nz ; do i=is,ie
490  dd%maxTKE(i,j,k) = maxtke(i,k)
491  enddo ; enddo ; endif
492  if (associated(dd%TKE_to_Kd)) then ; do k=1,nz ; do i=is,ie
493  dd%TKE_to_Kd(i,j,k) = tke_to_kd(i,k)
494  enddo ; enddo ; endif
495  endif
496 
497  ! Add the input turbulent diffusivity.
498  if (cs%useKappaShear .or. cs%use_CVMix_shear) then
499  if (present(kd_int)) then
500  do k=2,nz ; do i=is,ie
501  kd_int_2d(i,k) = visc%Kd_shear(i,j,k) + 0.5 * (kd_lay_2d(i,k-1) + kd_lay_2d(i,k))
502  enddo ; enddo
503  do i=is,ie
504  kd_int_2d(i,1) = visc%Kd_shear(i,j,1) ! This isn't actually used. It could be 0.
505  kd_int_2d(i,nz+1) = 0.0
506  enddo
507  endif
508  do k=1,nz ; do i=is,ie
509  kd_lay_2d(i,k) = kd_lay_2d(i,k) + 0.5 * (visc%Kd_shear(i,j,k) + visc%Kd_shear(i,j,k+1))
510  enddo ; enddo
511  else
512  if (present(kd_int)) then
513  do i=is,ie
514  kd_int_2d(i,1) = kd_lay_2d(i,1) ; kd_int_2d(i,nz+1) = 0.0
515  enddo
516  do k=2,nz ; do i=is,ie
517  kd_int_2d(i,k) = 0.5 * (kd_lay_2d(i,k-1) + kd_lay_2d(i,k))
518  enddo ; enddo
519  endif
520  endif
521 
522  if (present(kd_int)) then
523  ! Add the ML_Rad diffusivity.
524  if (cs%ML_radiation) &
525  call add_mlrad_diffusivity(h, fluxes, j, g, gv, us, cs, tke_to_kd, kd_lay_2d, kd_int_2d)
526 
527  ! Add the Nikurashin and / or tidal bottom-driven mixing
528  if (cs%use_tidal_mixing) &
529  call calculate_tidal_mixing(h, n2_bot, j, tke_to_kd, maxtke, g, gv, us, cs%tidal_mixing_CSp, &
530  n2_lay, n2_int, kd_lay_2d, kd_int_2d, cs%Kd_max, visc%Kv_slow)
531 
532  ! This adds the diffusion sustained by the energy extracted from the flow by the bottom drag.
533  if (cs%bottomdraglaw .and. (cs%BBL_effic>0.0)) then
534  if (cs%use_LOTW_BBL_diffusivity) then
535  call add_lotw_bbl_diffusivity(h, u, v, tv, fluxes, visc, j, n2_int, g, gv, us, cs, &
536  dd%Kd_BBL, kd_lay_2d, kd_int_2d)
537  else
538  call add_drag_diffusivity(h, u, v, tv, fluxes, visc, j, tke_to_kd, &
539  maxtke, kb, g, gv, us, cs, kd_lay_2d, kd_int_2d, dd%Kd_BBL)
540  endif
541  endif
542 
543  if (cs%limit_dissipation) then
544  ! This calculates the dissipation ONLY from Kd calculated in this routine
545  ! dissip has units of W/m3 (= kg/m3 * m2/s * 1/s2)
546  ! 1) a global constant,
547  ! 2) a dissipation proportional to N (aka Gargett) and
548  ! 3) dissipation corresponding to a (nearly) constant diffusivity.
549  do k=2,nz ; do i=is,ie
550  dissip = max( cs%dissip_min, & ! Const. floor on dissip.
551  cs%dissip_N0 + cs%dissip_N1 * sqrt(n2_int(i,k)), & ! Floor aka Gargett
552  cs%dissip_N2 * n2_int(i,k)) ! Floor of Kd_min*rho0/F_Ri
553  kd_int_2d(i,k) = max(kd_int_2d(i,k) , & ! Apply floor to Kd
554  dissip * (cs%FluxRi_max / (gv%Rho0 * (n2_int(i,k) + omega2))))
555  enddo ; enddo
556  endif
557 
558  ! Optionally add a uniform diffusivity at the interfaces.
559  if (cs%Kd_add > 0.0) then ; do k=1,nz+1 ; do i=is,ie
560  kd_int_2d(i,k) = kd_int_2d(i,k) + cs%Kd_add
561  enddo ; enddo ; endif
562 
563  ! Copy the 2-d slices into the 3-d array that is exported.
564  do k=1,nz+1 ; do i=is,ie
565  kd_int(i,j,k) = kd_int_2d(i,k)
566  enddo ; enddo
567 
568  else ! Kd_int is not present.
569 
570  ! Add the ML_Rad diffusivity.
571  if (cs%ML_radiation) &
572  call add_mlrad_diffusivity(h, fluxes, j, g, gv, us, cs, tke_to_kd, kd_lay_2d)
573 
574  ! Add the Nikurashin and / or tidal bottom-driven mixing
575  if (cs%use_tidal_mixing) &
576  call calculate_tidal_mixing(h, n2_bot, j, tke_to_kd, maxtke, g, gv, us, cs%tidal_mixing_CSp, &
577  n2_lay, n2_int, kd_lay_2d, kd_max=cs%Kd_max, kv=visc%Kv_slow)
578 
579  ! This adds the diffusion sustained by the energy extracted from the flow by the bottom drag.
580  if (cs%bottomdraglaw .and. (cs%BBL_effic>0.0)) then
581  if (cs%use_LOTW_BBL_diffusivity) then
582  call add_lotw_bbl_diffusivity(h, u, v, tv, fluxes, visc, j, n2_int, g, gv, us, cs, &
583  dd%Kd_BBL, kd_lay_2d)
584  else
585  call add_drag_diffusivity(h, u, v, tv, fluxes, visc, j, tke_to_kd, &
586  maxtke, kb, g, gv, us, cs, kd_lay_2d, kd_bbl=dd%Kd_BBL)
587  endif
588  endif
589 
590  endif
591 
592  if (cs%limit_dissipation) then
593  ! This calculates the layer dissipation ONLY from Kd calculated in this routine
594  ! dissip has units of W/m3 (= kg/m3 * m2/s * 1/s2)
595  ! 1) a global constant,
596  ! 2) a dissipation proportional to N (aka Gargett) and
597  ! 3) dissipation corresponding to a (nearly) constant diffusivity.
598  do k=2,nz-1 ; do i=is,ie
599  dissip = max( cs%dissip_min, & ! Const. floor on dissip.
600  cs%dissip_N0 + cs%dissip_N1 * sqrt(n2_lay(i,k)), & ! Floor aka Gargett
601  cs%dissip_N2 * n2_lay(i,k)) ! Floor of Kd_min*rho0/F_Ri
602  kd_lay_2d(i,k) = max(kd_lay_2d(i,k) , & ! Apply floor to Kd
603  dissip * (cs%FluxRi_max / (gv%Rho0 * (n2_lay(i,k) + omega2))))
604  enddo ; enddo
605  endif
606 
607  if (associated(dd%Kd_work)) then
608  do k=1,nz ; do i=is,ie
609  dd%Kd_Work(i,j,k) = gv%Rho0 * kd_lay_2d(i,k) * n2_lay(i,k) * &
610  gv%H_to_Z*h(i,j,k) ! Watt m-2 s = kg s-3
611  enddo ; enddo
612  endif
613 
614  ! Optionally add a uniform diffusivity to the layers.
615  if ((cs%Kd_add > 0.0) .and. (present(kd_lay))) then
616  do k=1,nz ; do i=is,ie
617  kd_lay_2d(i,k) = kd_lay_2d(i,k) + cs%Kd_add
618  enddo ; enddo
619  endif
620 
621  ! Copy the 2-d slices into the 3-d array that is exported; this was done above for Kd_int.
622  if (present(kd_lay)) then ; do k=1,nz ; do i=is,ie
623  kd_lay(i,j,k) = kd_lay_2d(i,k)
624  enddo ; enddo ; endif
625  enddo ! j-loop
626 
627  if (cs%user_change_diff) then
628  call user_change_diff(h, tv, g, gv, us, cs%user_change_diff_CSp, kd_lay, kd_int, &
629  t_f, s_f, dd%Kd_user)
630  endif
631 
632  if (cs%debug) then
633  if (present(kd_lay)) call hchksum(kd_lay, "Kd_lay", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
634 
635  if (cs%useKappaShear) call hchksum(visc%Kd_shear, "Turbulent Kd", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
636 
637  if (cs%use_CVMix_ddiff) then
638  call hchksum(kd_extra_t, "MOM_set_diffusivity: Kd_extra_T", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
639  call hchksum(kd_extra_s, "MOM_set_diffusivity: Kd_extra_S", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
640  endif
641 
642  if (associated(visc%kv_bbl_u) .and. associated(visc%kv_bbl_v)) then
643  call uvchksum("BBL Kv_bbl_[uv]", visc%kv_bbl_u, visc%kv_bbl_v, g%HI, &
644  haloshift=0, symmetric=.true., scale=us%Z2_T_to_m2_s, &
645  scalar_pair=.true.)
646  endif
647 
648  if (associated(visc%bbl_thick_u) .and. associated(visc%bbl_thick_v)) then
649  call uvchksum("BBL bbl_thick_[uv]", visc%bbl_thick_u, visc%bbl_thick_v, &
650  g%HI, haloshift=0, symmetric=.true., scale=us%Z_to_m, &
651  scalar_pair=.true.)
652  endif
653 
654  if (associated(visc%Ray_u) .and. associated(visc%Ray_v)) then
655  call uvchksum("Ray_[uv]", visc%Ray_u, visc%Ray_v, g%HI, 0, symmetric=.true., scale=us%Z_to_m*us%s_to_T)
656  endif
657 
658  endif
659 
660  ! post diagnostics
661  if (present(kd_lay) .and. (cs%id_Kd_layer > 0)) call post_data(cs%id_Kd_layer, kd_lay, cs%diag)
662 
663  ! background mixing
664  if (cs%id_Kd_bkgnd > 0) call post_data(cs%id_Kd_bkgnd, dd%Kd_bkgnd, cs%diag)
665  if (cs%id_Kv_bkgnd > 0) call post_data(cs%id_Kv_bkgnd, dd%Kv_bkgnd, cs%diag)
666 
667  ! tidal mixing
668  call post_tidal_diagnostics(g, gv, h, cs%tidal_mixing_CSp)
669  if (cs%id_N2 > 0) call post_data(cs%id_N2, dd%N2_3d, cs%diag)
670  if (cs%id_Kd_Work > 0) call post_data(cs%id_Kd_Work, dd%Kd_Work, cs%diag)
671  if (cs%id_maxTKE > 0) call post_data(cs%id_maxTKE, dd%maxTKE, cs%diag)
672  if (cs%id_TKE_to_Kd > 0) call post_data(cs%id_TKE_to_Kd, dd%TKE_to_Kd, cs%diag)
673 
674  if (cs%id_Kd_user > 0) call post_data(cs%id_Kd_user, dd%Kd_user, cs%diag)
675 
676  ! double diffusive mixing
677  if (cs%double_diffusion) then
678  if (cs%id_KT_extra > 0) call post_data(cs%id_KT_extra, dd%KT_extra, cs%diag)
679  if (cs%id_KS_extra > 0) call post_data(cs%id_KS_extra, dd%KS_extra, cs%diag)
680  elseif (cs%use_CVMix_ddiff) then
681  if (cs%id_KT_extra > 0) call post_data(cs%id_KT_extra, kd_extra_t, cs%diag)
682  if (cs%id_KS_extra > 0) call post_data(cs%id_KS_extra, kd_extra_s, cs%diag)
683  if (cs%id_R_rho > 0) call post_data(cs%id_R_rho, dd%drho_rat, cs%diag)
684  endif
685  if (cs%id_Kd_BBL > 0) call post_data(cs%id_Kd_BBL, dd%Kd_BBL, cs%diag)
686 
687  if (associated(dd%N2_3d)) deallocate(dd%N2_3d)
688  if (associated(dd%Kd_work)) deallocate(dd%Kd_work)
689  if (associated(dd%Kd_user)) deallocate(dd%Kd_user)
690  if (associated(dd%maxTKE)) deallocate(dd%maxTKE)
691  if (associated(dd%TKE_to_Kd)) deallocate(dd%TKE_to_Kd)
692  if (associated(dd%KT_extra)) deallocate(dd%KT_extra)
693  if (associated(dd%KS_extra)) deallocate(dd%KS_extra)
694  if (associated(dd%drho_rat)) deallocate(dd%drho_rat)
695  if (associated(dd%Kd_BBL)) deallocate(dd%Kd_BBL)
696 
697  if (showcalltree) call calltree_leave("set_diffusivity()")
698 
699 end subroutine set_diffusivity
700 
701 !> Convert turbulent kinetic energy to diffusivity
702 subroutine find_tke_to_kd(h, tv, dRho_int, N2_lay, j, dt, G, GV, US, CS, &
703  TKE_to_Kd, maxTKE, kb)
704  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
705  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
706  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
707  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
708  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
709  type(thermo_var_ptrs), intent(in) :: tv !< Structure containing pointers to any available
710  !! thermodynamic fields.
711  real, dimension(SZI_(G),SZK_(G)+1), intent(in) :: dRho_int !< Change in locally referenced potential density
712  !! across each interface [R ~> kg m-3].
713  real, dimension(SZI_(G),SZK_(G)), intent(in) :: N2_lay !< The squared buoyancy frequency of the
714  !! layers [T-2 ~> s-2].
715  integer, intent(in) :: j !< j-index of row to work on
716  real, intent(in) :: dt !< Time increment [T ~> s].
717  type(set_diffusivity_cs), pointer :: CS !< Diffusivity control structure
718  real, dimension(SZI_(G),SZK_(G)), intent(out) :: TKE_to_Kd !< The conversion rate between the
719  !! TKE dissipated within a layer and the
720  !! diapycnal diffusivity witin that layer,
721  !! usually (~Rho_0 / (G_Earth * dRho_lay))
722  !! [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
723  real, dimension(SZI_(G),SZK_(G)), intent(out) :: maxTKE !< The energy required to for a layer to entrain
724  !! to its maximum realizable thickness [Z3 T-3 ~> m3 s-3]
725  integer, dimension(SZI_(G)), intent(out) :: kb !< Index of lightest layer denser than the buffer
726  !! layer, or -1 without a bulk mixed layer.
727  ! Local variables
728  real, dimension(SZI_(G),SZK_(G)) :: &
729  ds_dsp1, & ! coordinate variable (sigma-2) difference across an
730  ! interface divided by the difference across the interface
731  ! below it [nondim]
732  dsp1_ds, & ! inverse coordinate variable (sigma-2) difference
733  ! across an interface times the difference across the
734  ! interface above it [nondim]
735  rho_0, & ! Layer potential densities relative to surface pressure [R ~> kg m-3]
736  maxent ! maxEnt is the maximum value of entrainment from below (with
737  ! compensating entrainment from above to keep the layer
738  ! density from changing) that will not deplete all of the
739  ! layers above or below a layer within a timestep [Z ~> m].
740  real, dimension(SZI_(G)) :: &
741  htot, & ! total thickness above or below a layer, or the
742  ! integrated thickness in the BBL [Z ~> m].
743  mfkb, & ! total thickness in the mixed and buffer layers times ds_dsp1 [Z ~> m].
744  p_ref, & ! array of tv%P_Ref pressures [R L2 T-2 ~> Pa]
745  rcv_kmb, & ! coordinate density in the lowest buffer layer [R ~> kg m-3]
746  p_0 ! An array of 0 pressures [R L2 T-2 ~> Pa]
747 
748  real :: dh_max ! maximum amount of entrainment a layer could
749  ! undergo before entraining all fluid in the layers
750  ! above or below [Z ~> m].
751  real :: dRho_lay ! density change across a layer [R ~> kg m-3]
752  real :: Omega2 ! rotation rate squared [T-2 ~> s-2]
753  real :: G_Rho0 ! gravitation accel divided by Bouss ref density [Z T-2 R-1 ~> m4 s-2 kg-1]
754  real :: G_IRho0 ! Alternate calculation of G_Rho0 for reproducibility [Z T-2 R-1 ~> m4 s-2 kg-1]
755  real :: I_Rho0 ! inverse of Boussinesq reference density [R-1 ~> m3 kg-1]
756  real :: I_dt ! 1/dt [T-1 ~> s-1]
757  real :: H_neglect ! negligibly small thickness [H ~> m or kg m-2]
758  real :: hN2pO2 ! h (N^2 + Omega^2), in [Z T-2 ~> m s-2].
759  logical :: do_i(szi_(g))
760 
761  integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
762  integer :: i, k, is, ie, nz, i_rem, kmb, kb_min
763  is = g%isc ; ie = g%iec ; nz = g%ke
764 
765  i_dt = 1.0 / dt
766  omega2 = cs%omega**2
767  h_neglect = gv%H_subroundoff
768  g_rho0 = (us%L_to_Z**2 * gv%g_Earth) / (gv%Rho0)
769  if (cs%answers_2018) then
770  i_rho0 = 1.0 / (gv%Rho0)
771  g_irho0 = (us%L_to_Z**2 * gv%g_Earth) * i_rho0
772  else
773  g_irho0 = g_rho0
774  endif
775 
776  ! Simple but coordinate-independent estimate of Kd/TKE
777  if (cs%simple_TKE_to_Kd) then
778  do k=1,nz ; do i=is,ie
779  hn2po2 = (gv%H_to_Z * h(i,j,k)) * (n2_lay(i,k) + omega2) ! Units of Z T-2.
780  if (hn2po2>0.) then
781  tke_to_kd(i,k) = 1.0 / hn2po2 ! Units of T2 Z-1.
782  else; tke_to_kd(i,k) = 0.; endif
783  ! The maximum TKE conversion we allow is really a statement
784  ! about the upper diffusivity we allow. Kd_max must be set.
785  maxtke(i,k) = hn2po2 * cs%Kd_max ! Units of Z3 T-3.
786  enddo ; enddo
787  kb(is:ie) = -1 ! kb should not be used by any code in non-layered mode -AJA
788  return
789  endif
790 
791  ! Determine kb - the index of the shallowest active interior layer.
792  if (cs%bulkmixedlayer) then
793  kmb = gv%nk_rho_varies
794  do i=is,ie ; p_0(i) = 0.0 ; p_ref(i) = tv%P_Ref ; enddo
795  eosdom(:) = eos_domain(g%HI)
796  do k=1,nz
797  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_0, rho_0(:,k), tv%eqn_of_state, eosdom)
798  enddo
799  call calculate_density(tv%T(:,j,kmb), tv%S(:,j,kmb), p_ref, rcv_kmb, tv%eqn_of_state, eosdom)
800 
801  kb_min = kmb+1
802  do i=is,ie
803  ! Determine the next denser layer than the buffer layer in the
804  ! coordinate density (sigma-2).
805  do k=kmb+1,nz-1 ; if (rcv_kmb(i) <= gv%Rlay(k)) exit ; enddo
806  kb(i) = k
807 
808  ! Backtrack, in case there are massive layers above that are stable
809  ! in sigma-0.
810  do k=kb(i)-1,kmb+1,-1
811  if (rho_0(i,kmb) > rho_0(i,k)) exit
812  if (h(i,j,k)>2.0*gv%Angstrom_H) kb(i) = k
813  enddo
814  enddo
815 
816  call set_density_ratios(h, tv, kb, g, gv, us, cs, j, ds_dsp1, rho_0)
817  else ! not bulkmixedlayer
818  kb_min = 2 ; kmb = 0
819  do i=is,ie ; kb(i) = 1 ; enddo
820  call set_density_ratios(h, tv, kb, g, gv, us, cs, j, ds_dsp1)
821  endif
822 
823  ! Determine maxEnt - the maximum permitted entrainment from below by each
824  ! interior layer.
825  do k=2,nz-1 ; do i=is,ie
826  dsp1_ds(i,k) = 1.0 / ds_dsp1(i,k)
827  enddo ; enddo
828  do i=is,ie ; dsp1_ds(i,nz) = 0.0 ; enddo
829 
830  if (cs%bulkmixedlayer) then
831  kmb = gv%nk_rho_varies
832  do i=is,ie
833  htot(i) = gv%H_to_Z*h(i,j,kmb)
834  mfkb(i) = 0.0
835  if (kb(i) < nz) &
836  mfkb(i) = ds_dsp1(i,kb(i)) * (gv%H_to_Z*(h(i,j,kmb) - gv%Angstrom_H))
837  enddo
838  do k=1,kmb-1 ; do i=is,ie
839  htot(i) = htot(i) + gv%H_to_Z*h(i,j,k)
840  mfkb(i) = mfkb(i) + ds_dsp1(i,k+1)*(gv%H_to_Z*(h(i,j,k) - gv%Angstrom_H))
841  enddo ; enddo
842  else
843  do i=is,i
844  maxent(i,1) = 0.0 ; htot(i) = gv%H_to_Z*(h(i,j,1) - gv%Angstrom_H)
845  enddo
846  endif
847  do k=kb_min,nz-1 ; do i=is,ie
848  if (k == kb(i)) then
849  maxent(i,kb(i)) = mfkb(i)
850  elseif (k > kb(i)) then
851  if (cs%answers_2018) then
852  maxent(i,k) = (1.0/dsp1_ds(i,k))*(maxent(i,k-1) + htot(i))
853  else
854  maxent(i,k) = ds_dsp1(i,k)*(maxent(i,k-1) + htot(i))
855  endif
856  htot(i) = htot(i) + gv%H_to_Z*(h(i,j,k) - gv%Angstrom_H)
857  endif
858  enddo ; enddo
859 
860  do i=is,ie
861  htot(i) = gv%H_to_Z*(h(i,j,nz) - gv%Angstrom_H) ; maxent(i,nz) = 0.0
862  do_i(i) = (g%mask2dT(i,j) > 0.5)
863  enddo
864  do k=nz-1,kb_min,-1
865  i_rem = 0
866  do i=is,ie ; if (do_i(i)) then
867  if (k<kb(i)) then ; do_i(i) = .false. ; cycle ; endif
868  i_rem = i_rem + 1 ! Count the i-rows that are still being worked on.
869  maxent(i,k) = min(maxent(i,k),dsp1_ds(i,k+1)*maxent(i,k+1) + htot(i))
870  htot(i) = htot(i) + gv%H_to_Z*(h(i,j,k) - gv%Angstrom_H)
871  endif ; enddo
872  if (i_rem == 0) exit
873  enddo ! k-loop
874 
875  ! Now set maxTKE and TKE_to_Kd.
876  do i=is,ie
877  maxtke(i,1) = 0.0 ; tke_to_kd(i,1) = 0.0
878  maxtke(i,nz) = 0.0 ; tke_to_kd(i,nz) = 0.0
879  enddo
880  do k=2,kmb ; do i=is,ie
881  maxtke(i,k) = 0.0
882  tke_to_kd(i,k) = 1.0 / ((n2_lay(i,k) + omega2) * &
883  (gv%H_to_Z*(h(i,j,k) + h_neglect)))
884  enddo ; enddo
885  do k=kmb+1,kb_min-1 ; do i=is,ie
886  ! These are the properties in the deeper mixed and buffer layers, and
887  ! should perhaps be revisited.
888  maxtke(i,k) = 0.0 ; tke_to_kd(i,k) = 0.0
889  enddo ; enddo
890  do k=kb_min,nz-1 ; do i=is,ie
891  if (k<kb(i)) then
892  maxtke(i,k) = 0.0
893  tke_to_kd(i,k) = 0.0
894  else
895  ! maxTKE is found by determining the kappa that gives maxEnt.
896  ! kappa_max = I_dt * dRho_int(i,K+1) * maxEnt(i,k) * &
897  ! (GV%H_to_Z*h(i,j,k) + dh_max) / dRho_lay
898  ! maxTKE(i,k) = (GV%g_Earth*US%L_to_Z**2) * dRho_lay * kappa_max
899  ! dRho_int should already be non-negative, so the max is redundant?
900  dh_max = maxent(i,k) * (1.0 + dsp1_ds(i,k))
901  drho_lay = 0.5 * max(drho_int(i,k) + drho_int(i,k+1), 0.0)
902  maxtke(i,k) = i_dt * (g_irho0 * &
903  (0.5*max(drho_int(i,k+1) + dsp1_ds(i,k)*drho_int(i,k), 0.0))) * &
904  ((gv%H_to_Z*h(i,j,k) + dh_max) * maxent(i,k))
905  ! TKE_to_Kd should be rho_InSitu / G_Earth * (delta rho_InSitu)
906  ! The omega^2 term in TKE_to_Kd is due to a rescaling of the efficiency of turbulent
907  ! mixing by a factor of N^2 / (N^2 + Omega^2), as proposed by Melet et al., 2013?
908  tke_to_kd(i,k) = 1.0 / (g_rho0 * drho_lay + &
909  cs%omega**2 * gv%H_to_Z*(h(i,j,k) + h_neglect))
910  endif
911  enddo ; enddo
912 
913 end subroutine find_tke_to_kd
914 
915 !> Calculate Brunt-Vaisala frequency, N^2.
916 subroutine find_n2(h, tv, T_f, S_f, fluxes, j, G, GV, US, CS, dRho_int, &
917  N2_lay, N2_int, N2_bot)
918  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
919  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
920  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
921  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
922  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
923  type(thermo_var_ptrs), intent(in) :: tv !< Structure containing pointers to any available
924  !! thermodynamic fields.
925  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
926  intent(in) :: T_f !< layer temperature with the values in massless layers
927  !! filled vertically by diffusion [degC].
928  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
929  intent(in) :: S_f !< Layer salinities with values in massless
930  !! layers filled vertically by diffusion [ppt].
931  type(forcing), intent(in) :: fluxes !< A structure of thermodynamic surface fluxes
932  integer, intent(in) :: j !< j-index of row to work on
933  type(set_diffusivity_cs), pointer :: CS !< Diffusivity control structure
934  real, dimension(SZI_(G),SZK_(G)+1), &
935  intent(out) :: dRho_int !< Change in locally referenced potential density
936  !! across each interface [R ~> kg m-3].
937  real, dimension(SZI_(G),SZK_(G)+1), &
938  intent(out) :: N2_int !< The squared buoyancy frequency at the interfaces [T-2 ~> s-2].
939  real, dimension(SZI_(G),SZK_(G)), &
940  intent(out) :: N2_lay !< The squared buoyancy frequency of the layers [T-2 ~> s-2].
941  real, dimension(SZI_(G)), intent(out) :: N2_bot !< The near-bottom squared buoyancy frequency [T-2 ~> s-2].
942  ! Local variables
943  real, dimension(SZI_(G),SZK_(G)+1) :: &
944  dRho_int_unfilt, & ! unfiltered density differences across interfaces [R ~> kg m-3]
945  dRho_dT, & ! partial derivative of density wrt temp [R degC-1 ~> kg m-3 degC-1]
946  dRho_dS ! partial derivative of density wrt saln [R ppt-1 ~> kg m-3 ppt-1]
947 
948  real, dimension(SZI_(G)) :: &
949  pres, & ! pressure at each interface [R L2 T-2 ~> Pa]
950  Temp_int, & ! temperature at each interface [degC]
951  Salin_int, & ! salinity at each interface [ppt]
952  drho_bot, & ! A density difference [R ~> kg m-3]
953  h_amp, & ! The topographic roughness amplitude [Z ~> m].
954  hb, & ! The thickness of the bottom layer [Z ~> m].
955  z_from_bot ! The hieght above the bottom [Z ~> m].
956 
957  real :: Rml_base ! density of the deepest variable density layer
958  real :: dz_int ! thickness associated with an interface [Z ~> m].
959  real :: G_Rho0 ! gravitation acceleration divided by Bouss reference density
960  ! times some unit conversion factors [Z T-2 R-1 ~> m4 s-2 kg-1].
961  real :: H_neglect ! negligibly small thickness, in the same units as h.
962 
963  logical :: do_i(szi_(g)), do_any
964  integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
965  integer :: i, k, is, ie, nz
966 
967  is = g%isc ; ie = g%iec ; nz = g%ke
968  g_rho0 = (us%L_to_Z**2 * gv%g_Earth) / (gv%Rho0)
969  h_neglect = gv%H_subroundoff
970 
971  ! Find the (limited) density jump across each interface.
972  do i=is,ie
973  drho_int(i,1) = 0.0 ; drho_int(i,nz+1) = 0.0
974  drho_int_unfilt(i,1) = 0.0 ; drho_int_unfilt(i,nz+1) = 0.0
975  enddo
976  if (associated(tv%eqn_of_state)) then
977  if (associated(fluxes%p_surf)) then
978  do i=is,ie ; pres(i) = fluxes%p_surf(i,j) ; enddo
979  else
980  do i=is,ie ; pres(i) = 0.0 ; enddo
981  endif
982  eosdom(:) = eos_domain(g%HI)
983  do k=2,nz
984  do i=is,ie
985  pres(i) = pres(i) + (gv%g_Earth*gv%H_to_RZ)*h(i,j,k-1)
986  temp_int(i) = 0.5 * (t_f(i,j,k) + t_f(i,j,k-1))
987  salin_int(i) = 0.5 * (s_f(i,j,k) + s_f(i,j,k-1))
988  enddo
989  call calculate_density_derivs(temp_int, salin_int, pres, drho_dt(:,k), drho_ds(:,k), &
990  tv%eqn_of_state, eosdom)
991  do i=is,ie
992  drho_int(i,k) = max(drho_dt(i,k)*(t_f(i,j,k) - t_f(i,j,k-1)) + &
993  drho_ds(i,k)*(s_f(i,j,k) - s_f(i,j,k-1)), 0.0)
994  drho_int_unfilt(i,k) = max(drho_dt(i,k)*(tv%T(i,j,k) - tv%T(i,j,k-1)) + &
995  drho_ds(i,k)*(tv%S(i,j,k) - tv%S(i,j,k-1)), 0.0)
996  enddo
997  enddo
998  else
999  do k=2,nz ; do i=is,ie
1000  drho_int(i,k) = gv%Rlay(k) - gv%Rlay(k-1)
1001  enddo ; enddo
1002  endif
1003 
1004  ! Set the buoyancy frequencies.
1005  do k=1,nz ; do i=is,ie
1006  n2_lay(i,k) = g_rho0 * 0.5*(drho_int(i,k) + drho_int(i,k+1)) / &
1007  (gv%H_to_Z*(h(i,j,k) + h_neglect))
1008  enddo ; enddo
1009  do i=is,ie ; n2_int(i,1) = 0.0 ; n2_int(i,nz+1) = 0.0 ; enddo
1010  do k=2,nz ; do i=is,ie
1011  n2_int(i,k) = g_rho0 * drho_int(i,k) / &
1012  (0.5*gv%H_to_Z*(h(i,j,k-1) + h(i,j,k) + h_neglect))
1013  enddo ; enddo
1014 
1015  ! Find the bottom boundary layer stratification, and use this in the deepest layers.
1016  do i=is,ie
1017  hb(i) = 0.0 ; drho_bot(i) = 0.0 ; h_amp(i) = 0.0
1018  z_from_bot(i) = 0.5*gv%H_to_Z*h(i,j,nz)
1019  do_i(i) = (g%mask2dT(i,j) > 0.5)
1020  enddo
1021  if (cs%use_tidal_mixing) call tidal_mixing_h_amp(h_amp, g, j, cs%tidal_mixing_CSp)
1022 
1023  do k=nz,2,-1
1024  do_any = .false.
1025  do i=is,ie ; if (do_i(i)) then
1026  dz_int = 0.5*gv%H_to_Z*(h(i,j,k) + h(i,j,k-1))
1027  z_from_bot(i) = z_from_bot(i) + dz_int ! middle of the layer above
1028 
1029  hb(i) = hb(i) + dz_int
1030  drho_bot(i) = drho_bot(i) + drho_int(i,k)
1031 
1032  if (z_from_bot(i) > h_amp(i)) then
1033  if (k>2) then
1034  ! Always include at least one full layer.
1035  hb(i) = hb(i) + 0.5*gv%H_to_Z*(h(i,j,k-1) + h(i,j,k-2))
1036  drho_bot(i) = drho_bot(i) + drho_int(i,k-1)
1037  endif
1038  do_i(i) = .false.
1039  else
1040  do_any = .true.
1041  endif
1042  endif ; enddo
1043  if (.not.do_any) exit
1044  enddo
1045 
1046  do i=is,ie
1047  if (hb(i) > 0.0) then
1048  n2_bot(i) = (g_rho0 * drho_bot(i)) / hb(i)
1049  else ; n2_bot(i) = 0.0 ; endif
1050  z_from_bot(i) = 0.5*gv%H_to_Z*h(i,j,nz)
1051  do_i(i) = (g%mask2dT(i,j) > 0.5)
1052  enddo
1053 
1054  do k=nz,2,-1
1055  do_any = .false.
1056  do i=is,ie ; if (do_i(i)) then
1057  dz_int = 0.5*gv%H_to_Z*(h(i,j,k) + h(i,j,k-1))
1058  z_from_bot(i) = z_from_bot(i) + dz_int ! middle of the layer above
1059 
1060  n2_int(i,k) = n2_bot(i)
1061  if (k>2) n2_lay(i,k-1) = n2_bot(i)
1062 
1063  if (z_from_bot(i) > h_amp(i)) then
1064  if (k>2) n2_int(i,k-1) = n2_bot(i)
1065  do_i(i) = .false.
1066  else
1067  do_any = .true.
1068  endif
1069  endif ; enddo
1070  if (.not.do_any) exit
1071  enddo
1072 
1073  if (associated(tv%eqn_of_state)) then
1074  do k=1,nz+1 ; do i=is,ie
1075  drho_int(i,k) = drho_int_unfilt(i,k)
1076  enddo ; enddo
1077  endif
1078 
1079 end subroutine find_n2
1080 
1081 !> This subroutine sets the additional diffusivities of temperature and
1082 !! salinity due to double diffusion, using the same functional form as is
1083 !! used in MOM4.1, and taken from an NCAR technical note (REF?) that updates
1084 !! what was in Large et al. (1994). All the coefficients here should probably
1085 !! be made run-time variables rather than hard-coded constants.
1086 !!
1087 !! \todo Find reference for NCAR tech note above.
1088 subroutine double_diffusion(tv, h, T_f, S_f, j, G, GV, US, CS, Kd_T_dd, Kd_S_dd)
1089  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1090  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1091  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1092  type(thermo_var_ptrs), intent(in) :: tv !< Structure containing pointers to any available
1093  !! thermodynamic fields; absent fields have NULL ptrs.
1094  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1095  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
1096  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1097  intent(in) :: T_f !< layer temperatures with the values in massless layers
1098  !! filled vertically by diffusion [degC].
1099  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1100  intent(in) :: S_f !< Layer salinities with values in massless
1101  !! layers filled vertically by diffusion [ppt].
1102  integer, intent(in) :: j !< Meridional index upon which to work.
1103  type(set_diffusivity_cs), pointer :: CS !< Module control structure.
1104  real, dimension(SZI_(G),SZK_(G)+1), &
1105  intent(out) :: Kd_T_dd !< Interface double diffusion diapycnal
1106  !! diffusivity for temp [Z2 T-1 ~> m2 s-1].
1107  real, dimension(SZI_(G),SZK_(G)+1), &
1108  intent(out) :: Kd_S_dd !< Interface double diffusion diapycnal
1109  !! diffusivity for saln [Z2 T-1 ~> m2 s-1].
1110 
1111  real, dimension(SZI_(G)) :: &
1112  dRho_dT, & ! partial derivatives of density wrt temp [R degC-1 ~> kg m-3 degC-1]
1113  dRho_dS, & ! partial derivatives of density wrt saln [R ppt-1 ~> kg m-3 ppt-1]
1114  pres, & ! pressure at each interface [R L2 T-2 ~> Pa]
1115  Temp_int, & ! temperature at interfaces [degC]
1116  Salin_int ! Salinity at interfaces [ppt]
1117 
1118  real :: alpha_dT ! density difference between layers due to temp diffs [R ~> kg m-3]
1119  real :: beta_dS ! density difference between layers due to saln diffs [R ~> kg m-3]
1120 
1121  real :: Rrho ! vertical density ratio [nondim]
1122  real :: diff_dd ! factor for double-diffusion [nondim]
1123  real :: Kd_dd ! The dominant double diffusive diffusivity [Z2 T-1 ~> m2 s-1]
1124  real :: prandtl ! flux ratio for diffusive convection regime
1125 
1126  real, parameter :: Rrho0 = 1.9 ! limit for double-diffusive density ratio [nondim]
1127 
1128  integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
1129  integer :: i, k, is, ie, nz
1130  is = g%isc ; ie = g%iec ; nz = g%ke
1131 
1132  if (associated(tv%eqn_of_state)) then
1133  do i=is,ie
1134  pres(i) = 0.0 ; kd_t_dd(i,1) = 0.0 ; kd_s_dd(i,1) = 0.0
1135  kd_t_dd(i,nz+1) = 0.0 ; kd_s_dd(i,nz+1) = 0.0
1136  enddo
1137  if (associated(tv%p_surf)) then ; do i=is,ie ; pres(i) = tv%p_surf(i,j) ; enddo ; endif
1138  eosdom(:) = eos_domain(g%HI)
1139  do k=2,nz
1140  do i=is,ie
1141  pres(i) = pres(i) + (gv%g_Earth*gv%H_to_RZ)*h(i,j,k-1)
1142  temp_int(i) = 0.5 * (t_f(i,j,k-1) + t_f(i,j,k))
1143  salin_int(i) = 0.5 * (s_f(i,j,k-1) + s_f(i,j,k))
1144  enddo
1145  call calculate_density_derivs(temp_int, salin_int, pres, drho_dt, drho_ds, &
1146  tv%eqn_of_state, eosdom)
1147 
1148  do i=is,ie
1149  alpha_dt = -1.0*drho_dt(i) * (t_f(i,j,k-1) - t_f(i,j,k))
1150  beta_ds = drho_ds(i) * (s_f(i,j,k-1) - s_f(i,j,k))
1151 
1152  if ((alpha_dt > beta_ds) .and. (beta_ds > 0.0)) then ! salt finger case
1153  rrho = min(alpha_dt / beta_ds, rrho0)
1154  diff_dd = 1.0 - ((rrho-1.0)/(rrho0-1.0))
1155  kd_dd = cs%Max_salt_diff_salt_fingers * diff_dd*diff_dd*diff_dd
1156  kd_t_dd(i,k) = 0.7 * kd_dd
1157  kd_s_dd(i,k) = kd_dd
1158  elseif ((alpha_dt < 0.) .and. (beta_ds < 0.) .and. (alpha_dt > beta_ds)) then ! diffusive convection
1159  rrho = alpha_dt / beta_ds
1160  kd_dd = cs%Kv_molecular * 0.909 * exp(4.6 * exp(-0.54 * (1/rrho - 1)))
1161  prandtl = 0.15*rrho
1162  if (rrho > 0.5) prandtl = (1.85-0.85/rrho)*rrho
1163  kd_t_dd(i,k) = kd_dd
1164  kd_s_dd(i,k) = prandtl * kd_dd
1165  else
1166  kd_t_dd(i,k) = 0.0 ; kd_s_dd(i,k) = 0.0
1167  endif
1168  enddo
1169  enddo
1170  endif
1171 
1172 end subroutine double_diffusion
1173 
1174 !> This routine adds diffusion sustained by flow energy extracted by bottom drag.
1175 subroutine add_drag_diffusivity(h, u, v, tv, fluxes, visc, j, TKE_to_Kd, &
1176  maxTKE, kb, G, GV, US, CS, Kd_lay, Kd_int, Kd_BBL)
1177  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1178  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
1179  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1180  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1181  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1]
1182  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1183  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1]
1184  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1185  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1186  type(thermo_var_ptrs), intent(in) :: tv !< Structure containing pointers to any available
1187  !! thermodynamic fields.
1188  type(forcing), intent(in) :: fluxes !< A structure of thermodynamic surface fluxes
1189  type(vertvisc_type), intent(in) :: visc !< Structure containing vertical viscosities, bottom
1190  !! boundary layer properies, and related fields
1191  integer, intent(in) :: j !< j-index of row to work on
1192  real, dimension(SZI_(G),SZK_(G)), intent(in) :: TKE_to_Kd !< The conversion rate between the TKE
1193  !! TKE dissipated within a layer and the
1194  !! diapycnal diffusivity witin that layer,
1195  !! usually (~Rho_0 / (G_Earth * dRho_lay))
1196  !! [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
1197  real, dimension(SZI_(G),SZK_(G)), intent(in) :: maxTKE !< The energy required to for a layer to entrain
1198  !! to its maximum-realizable thickness [Z3 T-3 ~> m3 s-3]
1199  integer, dimension(SZI_(G)), intent(in) :: kb !< Index of lightest layer denser than the buffer
1200  !! layer, or -1 without a bulk mixed layer
1201  type(set_diffusivity_cs), pointer :: CS !< Diffusivity control structure
1202  real, dimension(SZI_(G),SZK_(G)), intent(inout) :: Kd_lay !< The diapycnal diffusivity in layers,
1203  !! [Z2 T-1 ~> m2 s-1].
1204  real, dimension(SZI_(G),SZK_(G)+1), &
1205  optional, intent(inout) :: Kd_int !< The diapycnal diffusivity at interfaces,
1206  !! [Z2 T-1 ~> m2 s-1].
1207  real, dimension(:,:,:), pointer :: Kd_BBL !< Interface BBL diffusivity [Z2 T-1 ~> m2 s-1].
1208 
1209 ! This routine adds diffusion sustained by flow energy extracted by bottom drag.
1210 
1211  real, dimension(SZK_(G)+1) :: &
1212  Rint ! coordinate density of an interface [R ~> kg m-3]
1213  real, dimension(SZI_(G)) :: &
1214  htot, & ! total thickness above or below a layer, or the
1215  ! integrated thickness in the BBL [Z ~> m].
1216  rho_htot, & ! running integral with depth of density [Z R ~> kg m-2]
1217  gh_sum_top, & ! BBL value of g'h that can be supported by
1218  ! the local ustar, times R0_g [R ~> kg m-2]
1219  rho_top, & ! density at top of the BBL [R ~> kg m-3]
1220  tke, & ! turbulent kinetic energy available to drive
1221  ! bottom-boundary layer mixing in a layer [Z3 T-3 ~> m3 s-3]
1222  i2decay ! inverse of twice the TKE decay scale [Z-1 ~> m-1].
1223 
1224  real :: TKE_to_layer ! TKE used to drive mixing in a layer [Z3 T-3 ~> m3 s-3]
1225  real :: TKE_Ray ! TKE from layer Rayleigh drag used to drive mixing in layer [Z3 T-3 ~> m3 s-3]
1226  real :: TKE_here ! TKE that goes into mixing in this layer [Z3 T-3 ~> m3 s-3]
1227  real :: dRl, dRbot ! temporaries holding density differences [R ~> kg m-3]
1228  real :: cdrag_sqrt ! square root of the drag coefficient [nondim]
1229  real :: ustar_h ! value of ustar at a thickness point [Z T-1 ~> m s-1].
1230  real :: absf ! average absolute Coriolis parameter around a thickness point [T-1 ~> s-1]
1231  real :: R0_g ! Rho0 / G_Earth [R T2 Z-1 m-1 ~> kg s2 m-5]
1232  real :: I_rho0 ! 1 / RHO0 [R-1 ~> m3 kg-1]
1233  real :: delta_Kd ! increment to Kd from the bottom boundary layer mixing [Z2 T-1 ~> m2 s-1].
1234  logical :: Rayleigh_drag ! Set to true if Rayleigh drag velocities
1235  ! defined in visc, on the assumption that this
1236  ! extracted energy also drives diapycnal mixing.
1237 
1238  logical :: domore, do_i(szi_(g))
1239  logical :: do_diag_Kd_BBL
1240 
1241  integer :: i, k, is, ie, nz, i_rem, kb_min
1242  is = g%isc ; ie = g%iec ; nz = g%ke
1243 
1244  do_diag_kd_bbl = associated(kd_bbl)
1245 
1246  if (.not.(cs%bottomdraglaw .and. (cs%BBL_effic>0.0))) return
1247 
1248  cdrag_sqrt = sqrt(cs%cdrag)
1249  tke_ray = 0.0 ; rayleigh_drag = .false.
1250  if (associated(visc%Ray_u) .and. associated(visc%Ray_v)) rayleigh_drag = .true.
1251 
1252  i_rho0 = 1.0 / (gv%Rho0)
1253  r0_g = gv%Rho0 / (us%L_to_Z**2 * gv%g_Earth)
1254 
1255  do k=2,nz ; rint(k) = 0.5*(gv%Rlay(k-1)+gv%Rlay(k)) ; enddo
1256 
1257  kb_min = max(gv%nk_rho_varies+1,2)
1258 
1259  ! The turbulence decay scale is 0.5*ustar/f from K&E & MOM_vertvisc.F90
1260  ! Any turbulence that makes it into the mixed layers is assumed
1261  ! to be relatively small and is discarded.
1262  do i=is,ie
1263  ustar_h = visc%ustar_BBL(i,j)
1264  if (associated(fluxes%ustar_tidal)) &
1265  ustar_h = ustar_h + fluxes%ustar_tidal(i,j)
1266  absf = 0.25 * ((abs(g%CoriolisBu(i-1,j-1)) + abs(g%CoriolisBu(i,j))) + &
1267  (abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j-1))))
1268  if ((ustar_h > 0.0) .and. (absf > 0.5*cs%IMax_decay*ustar_h)) then
1269  i2decay(i) = absf / ustar_h
1270  else
1271  ! The maximum decay scale should be something of order 200 m.
1272  ! If ustar_h = 0, this is land so this value doesn't matter.
1273  i2decay(i) = 0.5*cs%IMax_decay
1274  endif
1275  tke(i) = ((cs%BBL_effic * cdrag_sqrt) * exp(-i2decay(i)*(gv%H_to_Z*h(i,j,nz))) ) * &
1276  visc%TKE_BBL(i,j)
1277 
1278  if (associated(fluxes%TKE_tidal)) &
1279  tke(i) = tke(i) + fluxes%TKE_tidal(i,j) * i_rho0 * &
1280  (cs%BBL_effic * exp(-i2decay(i)*(gv%H_to_Z*h(i,j,nz))))
1281 
1282  ! Distribute the work over a BBL of depth 20^2 ustar^2 / g' following
1283  ! Killworth & Edwards (1999) and Zilitikevich & Mironov (1996).
1284  ! Rho_top is determined by finding the density where
1285  ! integral(bottom, Z) (rho(z') - rho(Z)) dz' = rho_0 400 ustar^2 / g
1286 
1287  gh_sum_top(i) = r0_g * 400.0 * ustar_h**2
1288 
1289  do_i(i) = (g%mask2dT(i,j) > 0.5)
1290  htot(i) = gv%H_to_Z*h(i,j,nz)
1291  rho_htot(i) = gv%Rlay(nz)*(gv%H_to_Z*h(i,j,nz))
1292  rho_top(i) = gv%Rlay(1)
1293  if (cs%bulkmixedlayer .and. do_i(i)) rho_top(i) = gv%Rlay(kb(i)-1)
1294  enddo
1295 
1296  do k=nz-1,2,-1 ; domore = .false.
1297  do i=is,ie ; if (do_i(i)) then
1298  htot(i) = htot(i) + gv%H_to_Z*h(i,j,k)
1299  rho_htot(i) = rho_htot(i) + gv%Rlay(k)*(gv%H_to_Z*h(i,j,k))
1300  if (htot(i)*gv%Rlay(k-1) <= (rho_htot(i) - gh_sum_top(i))) then
1301  ! The top of the mixing is in the interface atop the current layer.
1302  rho_top(i) = (rho_htot(i) - gh_sum_top(i)) / htot(i)
1303  do_i(i) = .false.
1304  elseif (k <= kb(i)) then ; do_i(i) = .false.
1305  else ; domore = .true. ; endif
1306  endif ; enddo
1307  if (.not.domore) exit
1308  enddo ! k-loop
1309 
1310  do i=is,ie ; do_i(i) = (g%mask2dT(i,j) > 0.5) ; enddo
1311  do k=nz-1,kb_min,-1
1312  i_rem = 0
1313  do i=is,ie ; if (do_i(i)) then
1314  if (k<kb(i)) then ; do_i(i) = .false. ; cycle ; endif
1315  i_rem = i_rem + 1 ! Count the i-rows that are still being worked on.
1316  ! Apply vertical decay of the turbulent energy. This energy is
1317  ! simply lost.
1318  tke(i) = tke(i) * exp(-i2decay(i) * (gv%H_to_Z*(h(i,j,k) + h(i,j,k+1))))
1319 
1320 ! if (maxEnt(i,k) <= 0.0) cycle
1321  if (maxtke(i,k) <= 0.0) cycle
1322 
1323  ! This is an analytic integral where diffusity is a quadratic function of
1324  ! rho that goes asymptotically to 0 at Rho_top (vaguely following KPP?).
1325  if (tke(i) > 0.0) then
1326  if (rint(k) <= rho_top(i)) then
1327  tke_to_layer = tke(i)
1328  else
1329  drl = rint(k+1) - rint(k) ; drbot = rint(k+1) - rho_top(i)
1330  tke_to_layer = tke(i) * drl * &
1331  (3.0*drbot*(rint(k) - rho_top(i)) + drl**2) / (drbot**3)
1332  endif
1333  else ; tke_to_layer = 0.0 ; endif
1334 
1335  ! TKE_Ray has been initialized to 0 above.
1336  if (rayleigh_drag) tke_ray = 0.5*cs%BBL_effic * us%L_to_Z**2 * g%IareaT(i,j) * &
1337  ((g%areaCu(i-1,j) * visc%Ray_u(i-1,j,k) * u(i-1,j,k)**2 + &
1338  g%areaCu(i,j) * visc%Ray_u(i,j,k) * u(i,j,k)**2) + &
1339  (g%areaCv(i,j-1) * visc%Ray_v(i,j-1,k) * v(i,j-1,k)**2 + &
1340  g%areaCv(i,j) * visc%Ray_v(i,j,k) * v(i,j,k)**2))
1341 
1342  if (tke_to_layer + tke_ray > 0.0) then
1343  if (cs%BBL_mixing_as_max) then
1344  if (tke_to_layer + tke_ray > maxtke(i,k)) &
1345  tke_to_layer = maxtke(i,k) - tke_ray
1346 
1347  tke(i) = tke(i) - tke_to_layer
1348 
1349  if (kd_lay(i,k) < (tke_to_layer + tke_ray) * tke_to_kd(i,k)) then
1350  delta_kd = (tke_to_layer + tke_ray) * tke_to_kd(i,k) - kd_lay(i,k)
1351  if ((cs%Kd_max >= 0.0) .and. (delta_kd > cs%Kd_max)) then
1352  delta_kd = cs%Kd_max
1353  kd_lay(i,k) = kd_lay(i,k) + delta_kd
1354  else
1355  kd_lay(i,k) = (tke_to_layer + tke_ray) * tke_to_kd(i,k)
1356  endif
1357  if (present(kd_int)) then
1358  kd_int(i,k) = kd_int(i,k) + 0.5 * delta_kd
1359  kd_int(i,k+1) = kd_int(i,k+1) + 0.5 * delta_kd
1360  endif
1361  if (do_diag_kd_bbl) then
1362  kd_bbl(i,j,k) = kd_bbl(i,j,k) + 0.5 * delta_kd
1363  kd_bbl(i,j,k+1) = kd_bbl(i,j,k+1) + 0.5 * delta_kd
1364  endif
1365  endif
1366  else
1367  if (kd_lay(i,k) >= maxtke(i,k) * tke_to_kd(i,k)) then
1368  tke_here = 0.0
1369  tke(i) = tke(i) + tke_ray
1370  elseif (kd_lay(i,k) + (tke_to_layer + tke_ray) * tke_to_kd(i,k) > &
1371  maxtke(i,k) * tke_to_kd(i,k)) then
1372  tke_here = ((tke_to_layer + tke_ray) + kd_lay(i,k) / tke_to_kd(i,k)) - maxtke(i,k)
1373  tke(i) = (tke(i) - tke_here) + tke_ray
1374  else
1375  tke_here = tke_to_layer + tke_ray
1376  tke(i) = tke(i) - tke_to_layer
1377  endif
1378  if (tke(i) < 0.0) tke(i) = 0.0 ! This should be unnecessary?
1379 
1380  if (tke_here > 0.0) then
1381  delta_kd = tke_here * tke_to_kd(i,k)
1382  if (cs%Kd_max >= 0.0) delta_kd = min(delta_kd, cs%Kd_max)
1383  kd_lay(i,k) = kd_lay(i,k) + delta_kd
1384  if (present(kd_int)) then
1385  kd_int(i,k) = kd_int(i,k) + 0.5 * delta_kd
1386  kd_int(i,k+1) = kd_int(i,k+1) + 0.5 * delta_kd
1387  endif
1388  if (do_diag_kd_bbl) then
1389  kd_bbl(i,j,k) = kd_bbl(i,j,k) + 0.5 * delta_kd
1390  kd_bbl(i,j,k+1) = kd_bbl(i,j,k+1) + 0.5 * delta_kd
1391  endif
1392  endif
1393  endif
1394  endif
1395 
1396  ! This may be risky - in the case that there are exactly zero
1397  ! velocities at 4 neighboring points, but nonzero velocities
1398  ! above the iterations would stop too soon. I don't see how this
1399  ! could happen in practice. RWH
1400  if ((tke(i)<= 0.0) .and. (tke_ray == 0.0)) then
1401  do_i(i) = .false. ; i_rem = i_rem - 1
1402  endif
1403 
1404  endif ; enddo
1405  if (i_rem == 0) exit
1406  enddo ! k-loop
1407 
1408 end subroutine add_drag_diffusivity
1409 
1410 !> Calculates a BBL diffusivity use a Prandtl number 1 diffusivity with a law of the
1411 !! wall turbulent viscosity, up to a BBL height where the energy used for mixing has
1412 !! consumed the mechanical TKE input.
1413 subroutine add_lotw_bbl_diffusivity(h, u, v, tv, fluxes, visc, j, N2_int, &
1414  G, GV, US, CS, Kd_BBL, Kd_lay, Kd_int)
1415  type(ocean_grid_type), intent(in) :: G !< Grid structure
1416  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
1417  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1418  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1419  intent(in) :: u !< u component of flow [L T-1 ~> m s-1]
1420  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1421  intent(in) :: v !< v component of flow [L T-1 ~> m s-1]
1422  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1423  intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
1424  type(thermo_var_ptrs), intent(in) :: tv !< Structure containing pointers to any available
1425  !! thermodynamic fields.
1426  type(forcing), intent(in) :: fluxes !< Surface fluxes structure
1427  type(vertvisc_type), intent(in) :: visc !< Structure containing vertical viscosities, bottom
1428  !! boundary layer properies, and related fields.
1429  integer, intent(in) :: j !< j-index of row to work on
1430  real, dimension(SZI_(G),SZK_(G)+1), &
1431  intent(in) :: N2_int !< Square of Brunt-Vaisala at interfaces [T-2 ~> s-2]
1432  type(set_diffusivity_cs), pointer :: CS !< Diffusivity control structure
1433  real, dimension(:,:,:), pointer :: Kd_BBL !< Interface BBL diffusivity [Z2 T-1 ~> m2 s-1]
1434  real, dimension(SZI_(G),SZK_(G)), &
1435  optional, intent(inout) :: Kd_lay !< Layer net diffusivity [Z2 T-1 ~> m2 s-1]
1436  real, dimension(SZI_(G),SZK_(G)+1), &
1437  optional, intent(inout) :: Kd_int !< Interface net diffusivity [Z2 T-1 ~> m2 s-1]
1438 
1439  ! Local variables
1440  real :: TKE_column ! net TKE input into the column [Z3 T-3 ~> m3 s-3]
1441  real :: TKE_remaining ! remaining TKE available for mixing in this layer and above [Z3 T-3 ~> m3 s-3]
1442  real :: TKE_consumed ! TKE used for mixing in this layer [Z3 T-3 ~> m3 s-3]
1443  real :: TKE_Kd_wall ! TKE associated with unlimited law of the wall mixing [Z3 T-3 ~> m3 s-3]
1444  real :: cdrag_sqrt ! square root of the drag coefficient [nondim]
1445  real :: ustar ! value of ustar at a thickness point [Z T-1 ~> m s-1].
1446  real :: ustar2 ! square of ustar, for convenience [Z2 T-2 ~> m2 s-2]
1447  real :: absf ! average absolute value of Coriolis parameter around a thickness point [T-1 ~> s-1]
1448  real :: dh, dhm1 ! thickness of layers k and k-1, respecitvely [Z ~> m].
1449  real :: z_bot ! distance to interface k from bottom [Z ~> m].
1450  real :: D_minus_z ! distance to interface k from surface [Z ~> m].
1451  real :: total_thickness ! total thickness of water column [Z ~> m].
1452  real :: Idecay ! inverse of decay scale used for "Joule heating" loss of TKE with height [Z-1 ~> m-1].
1453  real :: Kd_wall ! Law of the wall diffusivity [Z2 T-1 ~> m2 s-1].
1454  real :: Kd_lower ! diffusivity for lower interface [Z2 T-1 ~> m2 s-1]
1455  real :: ustar_D ! u* x D [Z2 T-1 ~> m2 s-1].
1456  real :: I_Rho0 ! 1 / rho0 [R-1 ~> m3 kg-1]
1457  real :: N2_min ! Minimum value of N2 to use in calculation of TKE_Kd_wall [T-2 ~> s-2]
1458  logical :: Rayleigh_drag ! Set to true if there are Rayleigh drag velocities defined in visc, on
1459  ! the assumption that this extracted energy also drives diapycnal mixing.
1460  integer :: i, k, km1
1461  real, parameter :: von_karm = 0.41 ! Von Karman constant (http://en.wikipedia.org/wiki/Von_Karman_constant)
1462  logical :: do_diag_Kd_BBL
1463 
1464  if (.not.(cs%bottomdraglaw .and. (cs%BBL_effic>0.0))) return
1465  do_diag_kd_bbl = associated(kd_bbl)
1466 
1467  n2_min = 0.
1468  if (cs%LOTW_BBL_use_omega) n2_min = cs%omega**2
1469 
1470  ! Determine whether to add Rayleigh drag contribution to TKE
1471  rayleigh_drag = .false.
1472  if (associated(visc%Ray_u) .and. associated(visc%Ray_v)) rayleigh_drag = .true.
1473  i_rho0 = 1.0 / (gv%Rho0)
1474  cdrag_sqrt = sqrt(cs%cdrag)
1475 
1476  do i=g%isc,g%iec ! Developed in single-column mode
1477 
1478  ! Column-wise parameters.
1479  absf = 0.25 * ((abs(g%CoriolisBu(i-1,j-1)) + abs(g%CoriolisBu(i,j))) + &
1480  (abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j-1)))) ! Non-zero on equator!
1481 
1482  ! u* at the bottom [Z T-1 ~> m s-1].
1483  ustar = visc%ustar_BBL(i,j)
1484  ustar2 = ustar**2
1485  ! In add_drag_diffusivity(), fluxes%ustar_tidal is added in. This might be double counting
1486  ! since ustar_BBL should already include all contributions to u*? -AJA
1487  !### Examine the question of whether there is double counting of fluxes%ustar_tidal.
1488  if (associated(fluxes%ustar_tidal)) ustar = ustar + fluxes%ustar_tidal(i,j)
1489 
1490  ! The maximum decay scale should be something of order 200 m. We use the smaller of u*/f and
1491  ! (IMax_decay)^-1 as the decay scale. If ustar = 0, this is land so this value doesn't matter.
1492  idecay = cs%IMax_decay
1493  if ((ustar > 0.0) .and. (absf > cs%IMax_decay * ustar)) idecay = absf / ustar
1494 
1495  ! Energy input at the bottom [Z3 T-3 ~> m3 s-3].
1496  ! (Note that visc%TKE_BBL is in [Z3 T-3 ~> m3 s-3], set in set_BBL_TKE().)
1497  ! I am still unsure about sqrt(cdrag) in this expressions - AJA
1498  tke_column = cdrag_sqrt * visc%TKE_BBL(i,j)
1499  ! Add in tidal dissipation energy at the bottom [R Z3 T-3 ~> m3 s-3].
1500  ! Note that TKE_tidal is in [R Z3 T-3 ~> W m-2].
1501  if (associated(fluxes%TKE_tidal)) &
1502  tke_column = tke_column + fluxes%TKE_tidal(i,j) * i_rho0
1503  tke_column = cs%BBL_effic * tke_column ! Only use a fraction of the mechanical dissipation for mixing.
1504 
1505  tke_remaining = tke_column
1506  total_thickness = ( sum(h(i,j,:)) + gv%H_subroundoff )* gv%H_to_Z ! Total column thickness [Z ~> m].
1507  ustar_d = ustar * total_thickness
1508  z_bot = 0.
1509  kd_lower = 0. ! Diffusivity on bottom boundary.
1510 
1511  ! Work upwards from the bottom, accumulating work used until it exceeds the available TKE input
1512  ! at the bottom.
1513  do k=g%ke,2,-1
1514  dh = gv%H_to_Z * h(i,j,k) ! Thickness of this level [Z ~> m].
1515  km1 = max(k-1, 1)
1516  dhm1 = gv%H_to_Z * h(i,j,km1) ! Thickness of level above [Z ~> m].
1517 
1518  ! Add in additional energy input from bottom-drag against slopes (sides)
1519  if (rayleigh_drag) tke_remaining = tke_remaining + &
1520  0.5*cs%BBL_effic * us%L_to_Z**2 * g%IareaT(i,j) * &
1521  ((g%areaCu(i-1,j) * visc%Ray_u(i-1,j,k) * u(i-1,j,k)**2 + &
1522  g%areaCu(i,j) * visc%Ray_u(i,j,k) * u(i,j,k)**2) + &
1523  (g%areaCv(i,j-1) * visc%Ray_v(i,j-1,k) * v(i,j-1,k)**2 + &
1524  g%areaCv(i,j) * visc%Ray_v(i,j,k) * v(i,j,k)**2))
1525 
1526  ! Exponentially decay TKE across the thickness of the layer.
1527  ! This is energy loss in addition to work done as mixing, apparently to Joule heating.
1528  tke_remaining = exp(-idecay*dh) * tke_remaining
1529 
1530  z_bot = z_bot + h(i,j,k)*gv%H_to_Z ! Distance between upper interface of layer and the bottom [Z ~> m].
1531  d_minus_z = max(total_thickness - z_bot, 0.) ! Thickness above layer [Z ~> m].
1532 
1533  ! Diffusivity using law of the wall, limited by rotation, at height z [Z2 T-1 ~> m2 s-1].
1534  ! This calculation is at the upper interface of the layer
1535  if ( ustar_d + absf * ( z_bot * d_minus_z ) == 0.) then
1536  kd_wall = 0.
1537  else
1538  kd_wall = ((von_karm * ustar2) * (z_bot * d_minus_z)) &
1539  / (ustar_d + absf * (z_bot * d_minus_z))
1540  endif
1541 
1542  ! TKE associated with Kd_wall [Z3 T-3 ~> m3 s-3].
1543  ! This calculation if for the volume spanning the interface.
1544  tke_kd_wall = kd_wall * 0.5 * (dh + dhm1) * max(n2_int(i,k), n2_min)
1545 
1546  ! Now bound Kd such that the associated TKE is no greater than available TKE for mixing.
1547  if (tke_kd_wall > 0.) then
1548  tke_consumed = min(tke_kd_wall, tke_remaining)
1549  kd_wall = (tke_consumed / tke_kd_wall) * kd_wall ! Scale Kd so that only TKE_consumed is used.
1550  else
1551  ! Either N2=0 or dh = 0.
1552  if (tke_remaining > 0.) then
1553  kd_wall = cs%Kd_max
1554  else
1555  kd_wall = 0.
1556  endif
1557  tke_consumed = 0.
1558  endif
1559 
1560  ! Now use up the appropriate about of TKE associated with the diffusivity chosen
1561  tke_remaining = tke_remaining - tke_consumed ! Note this will be non-negative
1562 
1563  ! Add this BBL diffusivity to the model net diffusivity.
1564  if (present(kd_int)) kd_int(i,k) = kd_int(i,k) + kd_wall
1565  if (present(kd_lay)) kd_lay(i,k) = kd_lay(i,k) + 0.5 * (kd_wall + kd_lower)
1566  kd_lower = kd_wall ! Store for next layer up.
1567  if (do_diag_kd_bbl) kd_bbl(i,j,k) = kd_wall
1568  enddo ! k
1569  enddo ! i
1570 
1571 end subroutine add_lotw_bbl_diffusivity
1572 
1573 !> This routine adds effects of mixed layer radiation to the layer diffusivities.
1574 subroutine add_mlrad_diffusivity(h, fluxes, j, G, GV, US, CS, TKE_to_Kd, Kd_lay, Kd_int)
1575  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1576  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
1577  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1578  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1579  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1580  type(forcing), intent(in) :: fluxes !< Surface fluxes structure
1581  integer, intent(in) :: j !< The j-index to work on
1582  type(set_diffusivity_cs), pointer :: CS !< Diffusivity control structure
1583  real, dimension(SZI_(G),SZK_(G)), intent(in) :: TKE_to_Kd !< The conversion rate between the TKE
1584  !! TKE dissipated within a layer and the
1585  !! diapycnal diffusivity witin that layer,
1586  !! usually (~Rho_0 / (G_Earth * dRho_lay))
1587  !! [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
1588  real, dimension(SZI_(G),SZK_(G)), &
1589  optional, intent(inout) :: Kd_lay !< The diapycnal diffusivity in layers [Z2 T-1 ~> m2 s-1].
1590  real, dimension(SZI_(G),SZK_(G)+1), &
1591  optional, intent(inout) :: Kd_int !< The diapycnal diffusivity at interfaces
1592  !! [Z2 T-1 ~> m2 s-1].
1593 
1594 ! This routine adds effects of mixed layer radiation to the layer diffusivities.
1595 
1596  real, dimension(SZI_(G)) :: h_ml ! Mixed layer thickness [Z ~> m].
1597  real, dimension(SZI_(G)) :: TKE_ml_flux ! Mixed layer TKE flux [Z3 T-3 ~> m3 s-3]
1598  real, dimension(SZI_(G)) :: I_decay ! A decay rate [Z-1 ~> m-1].
1599  real, dimension(SZI_(G)) :: Kd_mlr_ml ! Diffusivities associated with mixed layer radiation [Z2 T-1 ~> m2 s-1].
1600 
1601  real :: f_sq ! The square of the local Coriolis parameter or a related variable [T-2 ~> s-2].
1602  real :: h_ml_sq ! The square of the mixed layer thickness [Z2 ~> m2].
1603  real :: ustar_sq ! ustar squared [Z2 T-2 ~> m2 s-2]
1604  real :: Kd_mlr ! A diffusivity associated with mixed layer turbulence radiation [Z2 T-1 ~> m2 s-1].
1605  real :: C1_6 ! 1/6
1606  real :: Omega2 ! rotation rate squared [T-2 ~> s-2].
1607  real :: z1 ! layer thickness times I_decay [nondim]
1608  real :: dzL ! thickness converted to heights [Z ~> m].
1609  real :: I_decay_len2_TKE ! squared inverse decay lengthscale for
1610  ! TKE, as used in the mixed layer code [Z-2 ~> m-2].
1611  real :: h_neglect ! negligibly small thickness [Z ~> m].
1612 
1613  logical :: do_any, do_i(szi_(g))
1614  integer :: i, k, is, ie, nz, kml
1615  is = g%isc ; ie = g%iec ; nz = g%ke
1616 
1617  omega2 = cs%omega**2
1618  c1_6 = 1.0 / 6.0
1619  kml = gv%nkml
1620  h_neglect = gv%H_subroundoff*gv%H_to_Z
1621 
1622  if (.not.cs%ML_radiation) return
1623 
1624  do i=is,ie ; h_ml(i) = 0.0 ; do_i(i) = (g%mask2dT(i,j) > 0.5) ; enddo
1625  do k=1,kml ; do i=is,ie ; h_ml(i) = h_ml(i) + gv%H_to_Z*h(i,j,k) ; enddo ; enddo
1626 
1627  do i=is,ie ; if (do_i(i)) then
1628  if (cs%ML_omega_frac >= 1.0) then
1629  f_sq = 4.0 * omega2
1630  else
1631  f_sq = 0.25 * ((g%CoriolisBu(i,j)**2 + g%CoriolisBu(i-1,j-1)**2) + &
1632  (g%CoriolisBu(i,j-1)**2 + g%CoriolisBu(i-1,j)**2))
1633  if (cs%ML_omega_frac > 0.0) &
1634  f_sq = cs%ML_omega_frac * 4.0 * omega2 + (1.0 - cs%ML_omega_frac) * f_sq
1635  endif
1636 
1637  ustar_sq = max(fluxes%ustar(i,j), cs%ustar_min)**2
1638 
1639  tke_ml_flux(i) = (cs%mstar * cs%ML_rad_coeff) * (ustar_sq * (fluxes%ustar(i,j)))
1640  i_decay_len2_tke = cs%TKE_decay**2 * (f_sq / ustar_sq)
1641 
1642  if (cs%ML_rad_TKE_decay) &
1643  tke_ml_flux(i) = tke_ml_flux(i) * exp(-h_ml(i) * sqrt(i_decay_len2_tke))
1644 
1645  ! Calculate the inverse decay scale
1646  h_ml_sq = (cs%ML_rad_efold_coeff * (h_ml(i)+h_neglect))**2
1647  i_decay(i) = sqrt((i_decay_len2_tke * h_ml_sq + 1.0) / h_ml_sq)
1648 
1649  ! Average the dissipation layer kml+1, using
1650  ! a more accurate Taylor series approximations for very thin layers.
1651  z1 = (gv%H_to_Z*h(i,j,kml+1)) * i_decay(i)
1652  if (z1 > 1e-5) then
1653  kd_mlr = tke_ml_flux(i) * tke_to_kd(i,kml+1) * (1.0 - exp(-z1))
1654  else
1655  kd_mlr = tke_ml_flux(i) * tke_to_kd(i,kml+1) * (z1 * (1.0 - z1 * (0.5 - c1_6 * z1)))
1656  endif
1657  kd_mlr_ml(i) = min(kd_mlr, cs%ML_rad_kd_max)
1658  tke_ml_flux(i) = tke_ml_flux(i) * exp(-z1)
1659  endif ; enddo
1660 
1661  if (present(kd_lay)) then
1662  do k=1,kml+1 ; do i=is,ie ; if (do_i(i)) then
1663  kd_lay(i,k) = kd_lay(i,k) + kd_mlr_ml(i)
1664  endif ; enddo ; enddo
1665  endif
1666  if (present(kd_int)) then
1667  do k=2,kml+1 ; do i=is,ie ; if (do_i(i)) then
1668  kd_int(i,k) = kd_int(i,k) + kd_mlr_ml(i)
1669  endif ; enddo ; enddo
1670  if (kml<=nz-1) then ; do i=is,ie ; if (do_i(i)) then
1671  kd_int(i,kml+2) = kd_int(i,kml+2) + 0.5 * kd_mlr_ml(i)
1672  endif ; enddo ; endif
1673  endif
1674 
1675  do k=kml+2,nz-1
1676  do_any = .false.
1677  do i=is,ie ; if (do_i(i)) then
1678  dzl = gv%H_to_Z*h(i,j,k) ; z1 = dzl*i_decay(i)
1679  if (cs%ML_Rad_bug) then
1680  ! These expresssions are dimensionally inconsistent. -RWH
1681  ! This is supposed to be the integrated energy deposited in the layer,
1682  ! not the average over the layer as in these expressions.
1683  if (z1 > 1e-5) then
1684  kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * & ! Units of Z2 T-1
1685  us%m_to_Z * ((1.0 - exp(-z1)) / dzl) ! Units of m-1
1686  else
1687  kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * & ! Units of Z2 T-1
1688  us%m_to_Z * (i_decay(i) * (1.0 - z1 * (0.5 - c1_6*z1))) ! Units of m-1
1689  endif
1690  else
1691  if (z1 > 1e-5) then
1692  kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * (1.0 - exp(-z1))
1693  else
1694  kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * (z1 * (1.0 - z1 * (0.5 - c1_6*z1)))
1695  endif
1696  endif
1697  kd_mlr = min(kd_mlr, cs%ML_rad_kd_max)
1698  if (present(kd_lay)) then
1699  kd_lay(i,k) = kd_lay(i,k) + kd_mlr
1700  endif
1701  if (present(kd_int)) then
1702  kd_int(i,k) = kd_int(i,k) + 0.5 * kd_mlr
1703  kd_int(i,k+1) = kd_int(i,k+1) + 0.5 * kd_mlr
1704  endif
1705 
1706  tke_ml_flux(i) = tke_ml_flux(i) * exp(-z1)
1707  if (tke_ml_flux(i) * i_decay(i) < 0.1 * cs%Kd_min * omega2) then
1708  do_i(i) = .false.
1709  else ; do_any = .true. ; endif
1710  endif ; enddo
1711  if (.not.do_any) exit
1712  enddo
1713 
1714 end subroutine add_mlrad_diffusivity
1715 
1716 !> This subroutine calculates several properties related to bottom
1717 !! boundary layer turbulence.
1718 subroutine set_bbl_tke(u, v, h, fluxes, visc, G, GV, US, CS, OBC)
1719  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1720  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
1721  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1722  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1723  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1]
1724  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1725  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1]
1726  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1727  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1728  type(forcing), intent(in) :: fluxes !< A structure of thermodynamic surface fluxes
1729  type(vertvisc_type), intent(in) :: visc !< Structure containing vertical viscosities, bottom
1730  !! boundary layer properies, and related fields.
1731  type(set_diffusivity_cs), pointer :: CS !< Diffusivity control structure
1732  type(ocean_obc_type), optional, pointer :: OBC !< Open boundaries control structure.
1733 
1734  ! This subroutine calculates several properties related to bottom
1735  ! boundary layer turbulence.
1736 
1737  real, dimension(SZI_(G)) :: &
1738  htot ! total thickness above or below a layer, or the
1739  ! integrated thickness in the BBL [Z ~> m].
1740 
1741  real, dimension(SZIB_(G)) :: &
1742  uhtot, & ! running integral of u in the BBL [Z L T-1 ~> m2 s-1]
1743  ustar, & ! bottom boundary layer turbulence speed [Z T-1 ~> m s-1].
1744  u2_bbl ! square of the mean zonal velocity in the BBL [L2 T-2 ~> m2 s-2]
1745 
1746  real :: vhtot(szi_(g)) ! running integral of v in the BBL [Z L T-1 ~> m2 s-1]
1747 
1748  real, dimension(SZI_(G),SZJB_(G)) :: &
1749  vstar, & ! ustar at at v-points [Z T-1 ~> m s-1].
1750  v2_bbl ! square of average meridional velocity in BBL [L2 T-2 ~> m2 s-2]
1751 
1752  real :: cdrag_sqrt ! square root of the drag coefficient [nondim]
1753  real :: hvel ! thickness at velocity points [Z ~> m].
1754 
1755  logical :: domore, do_i(szi_(g))
1756  integer :: i, j, k, is, ie, js, je, nz
1757  integer :: l_seg
1758  logical :: local_open_u_BC, local_open_v_BC
1759  logical :: has_obc
1760 
1761  local_open_u_bc = .false.
1762  local_open_v_bc = .false.
1763  if (present(obc)) then ; if (associated(obc)) then
1764  local_open_u_bc = obc%open_u_BCs_exist_globally
1765  local_open_v_bc = obc%open_v_BCs_exist_globally
1766  endif ; endif
1767 
1768  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1769 
1770  if (.not.associated(cs)) call mom_error(fatal,"set_BBL_TKE: "//&
1771  "Module must be initialized before it is used.")
1772 
1773  if (.not.cs%bottomdraglaw .or. (cs%BBL_effic<=0.0)) then
1774  if (associated(visc%ustar_BBL)) then
1775  do j=js,je ; do i=is,ie ; visc%ustar_BBL(i,j) = 0.0 ; enddo ; enddo
1776  endif
1777  if (associated(visc%TKE_BBL)) then
1778  do j=js,je ; do i=is,ie ; visc%TKE_BBL(i,j) = 0.0 ; enddo ; enddo
1779  endif
1780  return
1781  endif
1782 
1783  cdrag_sqrt = sqrt(cs%cdrag)
1784 
1785  !$OMP parallel default(shared) private(do_i,vhtot,htot,domore,hvel,uhtot,ustar,u2_bbl)
1786  !$OMP do
1787  do j=js-1,je
1788  ! Determine ustar and the square magnitude of the velocity in the
1789  ! bottom boundary layer. Together these give the TKE source and
1790  ! vertical decay scale.
1791  do i=is,ie ; if ((g%mask2dCv(i,j) > 0.5) .and. (cdrag_sqrt*visc%bbl_thick_v(i,j) > 0.0)) then
1792  do_i(i) = .true. ; vhtot(i) = 0.0 ; htot(i) = 0.0
1793  vstar(i,j) = visc%Kv_bbl_v(i,j) / (cdrag_sqrt*visc%bbl_thick_v(i,j))
1794  else
1795  do_i(i) = .false. ; vstar(i,j) = 0.0 ; htot(i) = 0.0
1796  endif ; enddo
1797  do k=nz,1,-1
1798  domore = .false.
1799  do i=is,ie ; if (do_i(i)) then
1800  ! Determine if grid point is an OBC
1801  has_obc = .false.
1802  if (local_open_v_bc) then
1803  l_seg = obc%segnum_v(i,j)
1804  if (l_seg /= obc_none) then
1805  has_obc = obc%segment(l_seg)%open
1806  endif
1807  endif
1808 
1809  ! Compute h based on OBC state
1810  if (has_obc) then
1811  if (obc%segment(l_seg)%direction == obc_direction_n) then
1812  hvel = gv%H_to_Z*h(i,j,k)
1813  else
1814  hvel = gv%H_to_Z*h(i,j+1,k)
1815  endif
1816  else
1817  hvel = 0.5*gv%H_to_Z*(h(i,j,k) + h(i,j+1,k))
1818  endif
1819 
1820  if ((htot(i) + hvel) >= visc%bbl_thick_v(i,j)) then
1821  vhtot(i) = vhtot(i) + (visc%bbl_thick_v(i,j) - htot(i))*v(i,j,k)
1822  htot(i) = visc%bbl_thick_v(i,j)
1823  do_i(i) = .false.
1824  else
1825  vhtot(i) = vhtot(i) + hvel*v(i,j,k)
1826  htot(i) = htot(i) + hvel
1827  domore = .true.
1828  endif
1829  endif ; enddo
1830  if (.not.domore) exit
1831  enddo
1832  do i=is,ie ; if ((g%mask2dCv(i,j) > 0.5) .and. (htot(i) > 0.0)) then
1833  v2_bbl(i,j) = (vhtot(i)*vhtot(i))/(htot(i)*htot(i))
1834  else
1835  v2_bbl(i,j) = 0.0
1836  endif ; enddo
1837  enddo
1838  !$OMP do
1839  do j=js,je
1840  do i=is-1,ie ; if ((g%mask2dCu(i,j) > 0.5) .and. (cdrag_sqrt*visc%bbl_thick_u(i,j) > 0.0)) then
1841  do_i(i) = .true. ; uhtot(i) = 0.0 ; htot(i) = 0.0
1842  ustar(i) = visc%Kv_bbl_u(i,j) / (cdrag_sqrt*visc%bbl_thick_u(i,j))
1843  else
1844  do_i(i) = .false. ; ustar(i) = 0.0 ; htot(i) = 0.0
1845  endif ; enddo
1846  do k=nz,1,-1 ; domore = .false.
1847  do i=is-1,ie ; if (do_i(i)) then
1848  ! Determine if grid point is an OBC
1849  has_obc = .false.
1850  if (local_open_u_bc) then
1851  l_seg = obc%segnum_u(i,j)
1852  if (l_seg /= obc_none) then
1853  has_obc = obc%segment(l_seg)%open
1854  endif
1855  endif
1856 
1857  ! Compute h based on OBC state
1858  if (has_obc) then
1859  if (obc%segment(l_seg)%direction == obc_direction_e) then
1860  hvel = gv%H_to_Z*h(i,j,k)
1861  else ! OBC_DIRECTION_W
1862  hvel = gv%H_to_Z*h(i+1,j,k)
1863  endif
1864  else
1865  hvel = 0.5*gv%H_to_Z*(h(i,j,k) + h(i+1,j,k))
1866  endif
1867 
1868  if ((htot(i) + hvel) >= visc%bbl_thick_u(i,j)) then
1869  uhtot(i) = uhtot(i) + (visc%bbl_thick_u(i,j) - htot(i))*u(i,j,k)
1870  htot(i) = visc%bbl_thick_u(i,j)
1871  do_i(i) = .false.
1872  else
1873  uhtot(i) = uhtot(i) + hvel*u(i,j,k)
1874  htot(i) = htot(i) + hvel
1875  domore = .true.
1876  endif
1877  endif ; enddo
1878  if (.not.domore) exit
1879  enddo
1880  do i=is-1,ie ; if ((g%mask2dCu(i,j) > 0.5) .and. (htot(i) > 0.0)) then
1881  u2_bbl(i) = (uhtot(i)*uhtot(i))/(htot(i)*htot(i))
1882  else
1883  u2_bbl(i) = 0.0
1884  endif ; enddo
1885 
1886  do i=is,ie
1887  visc%ustar_BBL(i,j) = sqrt(0.5*g%IareaT(i,j) * &
1888  ((g%areaCu(i-1,j)*(ustar(i-1)*ustar(i-1)) + &
1889  g%areaCu(i,j)*(ustar(i)*ustar(i))) + &
1890  (g%areaCv(i,j-1)*(vstar(i,j-1)*vstar(i,j-1)) + &
1891  g%areaCv(i,j)*(vstar(i,j)*vstar(i,j))) ) )
1892  visc%TKE_BBL(i,j) = us%L_to_Z**2 * &
1893  (((g%areaCu(i-1,j)*(ustar(i-1)*u2_bbl(i-1)) + &
1894  g%areaCu(i,j) * (ustar(i)*u2_bbl(i))) + &
1895  (g%areaCv(i,j-1)*(vstar(i,j-1)*v2_bbl(i,j-1)) + &
1896  g%areaCv(i,j) * (vstar(i,j)*v2_bbl(i,j))) )*g%IareaT(i,j))
1897  enddo
1898  enddo
1899  !$OMP end parallel
1900 
1901 end subroutine set_bbl_tke
1902 
1903 subroutine set_density_ratios(h, tv, kb, G, GV, US, CS, j, ds_dsp1, rho_0)
1904  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1905  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1906  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1907  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
1908  type(thermo_var_ptrs), intent(in) :: tv !< Structure containing pointers to any
1909  !! available thermodynamic fields; absent
1910  !! fields have NULL ptrs.
1911  integer, dimension(SZI_(G)), intent(in) :: kb !< Index of lightest layer denser than the buffer
1912  !! layer, or -1 without a bulk mixed layer.
1913  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1914  type(set_diffusivity_cs), pointer :: CS !< Control structure returned by previous
1915  !! call to diabatic_entrain_init.
1916  integer, intent(in) :: j !< Meridional index upon which to work.
1917  real, dimension(SZI_(G),SZK_(G)), intent(out) :: ds_dsp1 !< Coordinate variable (sigma-2)
1918  !! difference across an interface divided by
1919  !! the difference across the interface below
1920  !! it [nondim]
1921  real, dimension(SZI_(G),SZK_(G)), &
1922  optional, intent(in) :: rho_0 !< Layer potential densities relative to
1923  !! surface press [R ~> kg m-3].
1924 
1925  ! Local variables
1926  real :: g_R0 ! g_R0 is a rescaled version of g/Rho [L2 Z-1 R-1 T-2 ~> m4 kg-1 s-2]
1927  real :: eps, tmp ! nondimensional temproray variables
1928  real :: a(szk_(g)), a_0(szk_(g)) ! nondimensional temporary variables
1929  real :: p_ref(szi_(g)) ! an array of tv%P_Ref pressures [R L2 T-2 ~> Pa]
1930  real :: Rcv(szi_(g),szk_(g)) ! coordinate density in the mixed and buffer layers [R ~> kg m-3]
1931  real :: I_Drho ! temporary variable [R-1 ~> m3 kg-1]
1932 
1933  integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
1934  integer :: i, k, k3, is, ie, nz, kmb
1935  is = g%isc ; ie = g%iec ; nz = g%ke
1936 
1937  do k=2,nz-1
1938  if (gv%g_prime(k+1) /= 0.0) then
1939  do i=is,ie
1940  ds_dsp1(i,k) = gv%g_prime(k) / gv%g_prime(k+1)
1941  enddo
1942  else
1943  do i=is,ie
1944  ds_dsp1(i,k) = 1.
1945  enddo
1946  endif
1947  enddo
1948 
1949  if (cs%bulkmixedlayer) then
1950  g_r0 = gv%g_Earth / (gv%Rho0)
1951  kmb = gv%nk_rho_varies
1952  eps = 0.1
1953  do i=is,ie ; p_ref(i) = tv%P_Ref ; enddo
1954  eosdom(:) = eos_domain(g%HI)
1955  do k=1,kmb
1956  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_ref, rcv(:,k), tv%eqn_of_state, eosdom)
1957  enddo
1958  do i=is,ie
1959  if (kb(i) <= nz-1) then
1960 ! Set up appropriately limited ratios of the reduced gravities of the
1961 ! interfaces above and below the buffer layer and the next denser layer.
1962  k = kb(i)
1963 
1964  i_drho = g_r0 / gv%g_prime(k+1)
1965  ! The indexing convention for a is appropriate for the interfaces.
1966  do k3=1,kmb
1967  a(k3+1) = (gv%Rlay(k) - rcv(i,k3)) * i_drho
1968  enddo
1969  if ((present(rho_0)) .and. (a(kmb+1) < 2.0*eps*ds_dsp1(i,k))) then
1970 ! If the buffer layer nearly matches the density of the layer below in the
1971 ! coordinate variable (sigma-2), use the sigma-0-based density ratio if it is
1972 ! greater (and stable).
1973  if ((rho_0(i,k) > rho_0(i,kmb)) .and. &
1974  (rho_0(i,k+1) > rho_0(i,k))) then
1975  i_drho = 1.0 / (rho_0(i,k+1)-rho_0(i,k))
1976  a_0(kmb+1) = min((rho_0(i,k)-rho_0(i,kmb)) * i_drho, ds_dsp1(i,k))
1977  if (a_0(kmb+1) > a(kmb+1)) then
1978  do k3=2,kmb
1979  a_0(k3) = a_0(kmb+1) + (rho_0(i,kmb)-rho_0(i,k3-1)) * i_drho
1980  enddo
1981  if (a(kmb+1) <= eps*ds_dsp1(i,k)) then
1982  do k3=2,kmb+1 ; a(k3) = a_0(k3) ; enddo
1983  else
1984 ! Alternative... tmp = 0.5*(1.0 - cos(PI*(a(K2+1)/(eps*ds_dsp1(i,k)) - 1.0)) )
1985  tmp = a(kmb+1)/(eps*ds_dsp1(i,k)) - 1.0
1986  do k3=2,kmb+1 ; a(k3) = tmp*a(k3) + (1.0-tmp)*a_0(k3) ; enddo
1987  endif
1988  endif
1989  endif
1990  endif
1991 
1992  ds_dsp1(i,k) = max(a(kmb+1),1e-5)
1993 
1994  do k3=2,kmb
1995 ! ds_dsp1(i,k3) = MAX(a(k3),1e-5)
1996  ! Deliberately treat convective instabilies of the upper mixed
1997  ! and buffer layers with respect to the deepest buffer layer as
1998  ! though they don't exist. They will be eliminated by the upcoming
1999  ! call to the mixedlayer code anyway.
2000  ! The indexing convention is appropriate for the interfaces.
2001  ds_dsp1(i,k3) = max(a(k3),ds_dsp1(i,k))
2002  enddo
2003  endif ! (kb(i) <= nz-1)
2004  enddo ! I-loop.
2005  endif ! bulkmixedlayer
2006 
2007 end subroutine set_density_ratios
2008 
2009 subroutine set_diffusivity_init(Time, G, GV, US, param_file, diag, CS, int_tide_CSp, halo_TS, &
2010  double_diffuse)
2011  type(time_type), intent(in) :: Time !< The current model time
2012  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
2013  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
2014  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
2015  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
2016  !! parameters.
2017  type(diag_ctrl), target, intent(inout) :: diag !< A structure used to regulate diagnostic output.
2018  type(set_diffusivity_cs), pointer :: CS !< pointer set to point to the module control
2019  !! structure.
2020  type(int_tide_cs), pointer :: int_tide_CSp !< A pointer to the internal tides control
2021  !! structure
2022  integer, optional, intent(out) :: halo_TS !< The halo size of tracer points that must be
2023  !! valid for the calculations in set_diffusivity.
2024  logical, optional, intent(out) :: double_diffuse !< If present, this indicates whether
2025  !! some version of double diffusion is being used.
2026 
2027  ! Local variables
2028  real :: decay_length
2029  logical :: ML_use_omega
2030  logical :: default_2018_answers
2031  ! This include declares and sets the variable "version".
2032 # include "version_variable.h"
2033  character(len=40) :: mdl = "MOM_set_diffusivity" ! This module's name.
2034  real :: omega_frac_dflt
2035  logical :: Bryan_Lewis_diffusivity ! If true, the background diapycnal diffusivity uses
2036  ! the Bryan-Lewis (1979) style tanh profile.
2037  integer :: i, j, is, ie, js, je
2038  integer :: isd, ied, jsd, jed
2039 
2040  if (associated(cs)) then
2041  call mom_error(warning, "diabatic_entrain_init called with an associated "// &
2042  "control structure.")
2043  return
2044  endif
2045  allocate(cs)
2046 
2047  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
2048  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
2049 
2050  cs%diag => diag
2051  if (associated(int_tide_csp)) cs%int_tide_CSp => int_tide_csp
2052 
2053  ! These default values always need to be set.
2054  cs%BBL_mixing_as_max = .true.
2055  cs%cdrag = 0.003 ; cs%BBL_effic = 0.0
2056  cs%bulkmixedlayer = (gv%nkml > 0)
2057 
2058  ! Read all relevant parameters and write them to the model log.
2059  call log_version(param_file, mdl, version, "")
2060 
2061  call get_param(param_file, mdl, "INPUTDIR", cs%inputdir, default=".")
2062  cs%inputdir = slasher(cs%inputdir)
2063  call get_param(param_file, mdl, "FLUX_RI_MAX", cs%FluxRi_max, &
2064  "The flux Richardson number where the stratification is "//&
2065  "large enough that N2 > omega2. The full expression for "//&
2066  "the Flux Richardson number is usually "//&
2067  "FLUX_RI_MAX*N2/(N2+OMEGA2).", units="nondim", default=0.2)
2068  call get_param(param_file, mdl, "OMEGA", cs%omega, &
2069  "The rotation rate of the earth.", units="s-1", default=7.2921e-5, scale=us%T_to_s)
2070 
2071  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
2072  "This sets the default value for the various _2018_ANSWERS parameters.", &
2073  default=.false.)
2074  call get_param(param_file, mdl, "SET_DIFF_2018_ANSWERS", cs%answers_2018, &
2075  "If true, use the order of arithmetic and expressions that recover the "//&
2076  "answers from the end of 2018. Otherwise, use updated and more robust "//&
2077  "forms of the same expressions.", default=default_2018_answers)
2078 
2079  ! CS%use_tidal_mixing is set to True if an internal tidal dissipation scheme is to be used.
2080  cs%use_tidal_mixing = tidal_mixing_init(time, g, gv, us, param_file, diag, &
2081  cs%tidal_mixing_CSp)
2082 
2083  call get_param(param_file, mdl, "ML_RADIATION", cs%ML_radiation, &
2084  "If true, allow a fraction of TKE available from wind "//&
2085  "work to penetrate below the base of the mixed layer "//&
2086  "with a vertical decay scale determined by the minimum "//&
2087  "of: (1) The depth of the mixed layer, (2) an Ekman "//&
2088  "length scale.", default=.false.)
2089  if (cs%ML_radiation) then
2090  ! This give a minimum decay scale that is typically much less than Angstrom.
2091  cs%ustar_min = 2e-4 * cs%omega * (gv%Angstrom_Z + gv%H_subroundoff*gv%H_to_Z)
2092 
2093  call get_param(param_file, mdl, "ML_RAD_EFOLD_COEFF", cs%ML_rad_efold_coeff, &
2094  "A coefficient that is used to scale the penetration "//&
2095  "depth for turbulence below the base of the mixed layer. "//&
2096  "This is only used if ML_RADIATION is true.", units="nondim", default=0.2)
2097  call get_param(param_file, mdl, "ML_RAD_BUG", cs%ML_rad_bug, &
2098  "If true use code with a bug that reduces the energy available "//&
2099  "in the transition layer by a factor of the inverse of the energy "//&
2100  "deposition lenthscale (in m).", default=.false.)
2101  call get_param(param_file, mdl, "ML_RAD_KD_MAX", cs%ML_rad_kd_max, &
2102  "The maximum diapycnal diffusivity due to turbulence "//&
2103  "radiated from the base of the mixed layer. "//&
2104  "This is only used if ML_RADIATION is true.", &
2105  units="m2 s-1", default=1.0e-3, scale=us%m2_s_to_Z2_T)
2106  call get_param(param_file, mdl, "ML_RAD_COEFF", cs%ML_rad_coeff, &
2107  "The coefficient which scales MSTAR*USTAR^3 to obtain "//&
2108  "the energy available for mixing below the base of the "//&
2109  "mixed layer. This is only used if ML_RADIATION is true.", &
2110  units="nondim", default=0.2)
2111  call get_param(param_file, mdl, "ML_RAD_APPLY_TKE_DECAY", cs%ML_rad_TKE_decay, &
2112  "If true, apply the same exponential decay to ML_rad as "//&
2113  "is applied to the other surface sources of TKE in the "//&
2114  "mixed layer code. This is only used if ML_RADIATION is true.", default=.true.)
2115  call get_param(param_file, mdl, "MSTAR", cs%mstar, &
2116  "The ratio of the friction velocity cubed to the TKE "//&
2117  "input to the mixed layer.", units="nondim", default=1.2)
2118  call get_param(param_file, mdl, "TKE_DECAY", cs%TKE_decay, &
2119  "The ratio of the natural Ekman depth to the TKE decay scale.", &
2120  units="nondim", default=2.5)
2121  call get_param(param_file, mdl, "ML_USE_OMEGA", ml_use_omega, &
2122  "If true, use the absolute rotation rate instead of the "//&
2123  "vertical component of rotation when setting the decay "//&
2124  "scale for turbulence.", default=.false., do_not_log=.true.)
2125  omega_frac_dflt = 0.0
2126  if (ml_use_omega) then
2127  call mom_error(warning, "ML_USE_OMEGA is depricated; use ML_OMEGA_FRAC=1.0 instead.")
2128  omega_frac_dflt = 1.0
2129  endif
2130  call get_param(param_file, mdl, "ML_OMEGA_FRAC", cs%ML_omega_frac, &
2131  "When setting the decay scale for turbulence, use this "//&
2132  "fraction of the absolute rotation rate blended with the "//&
2133  "local value of f, as sqrt((1-of)*f^2 + of*4*omega^2).", &
2134  units="nondim", default=omega_frac_dflt)
2135  endif
2136 
2137  call get_param(param_file, mdl, "BOTTOMDRAGLAW", cs%bottomdraglaw, &
2138  "If true, the bottom stress is calculated with a drag "//&
2139  "law of the form c_drag*|u|*u. The velocity magnitude "//&
2140  "may be an assumed value or it may be based on the actual "//&
2141  "velocity in the bottommost HBBL, depending on LINEAR_DRAG.", default=.true.)
2142  if (cs%bottomdraglaw) then
2143  call get_param(param_file, mdl, "CDRAG", cs%cdrag, &
2144  "The drag coefficient relating the magnitude of the "//&
2145  "velocity field to the bottom stress. CDRAG is only used "//&
2146  "if BOTTOMDRAGLAW is true.", units="nondim", default=0.003)
2147  call get_param(param_file, mdl, "BBL_EFFIC", cs%BBL_effic, &
2148  "The efficiency with which the energy extracted by "//&
2149  "bottom drag drives BBL diffusion. This is only "//&
2150  "used if BOTTOMDRAGLAW is true.", units="nondim", default=0.20)
2151  call get_param(param_file, mdl, "BBL_MIXING_MAX_DECAY", decay_length, &
2152  "The maximum decay scale for the BBL diffusion, or 0 to allow the mixing "//&
2153  "to penetrate as far as stratification and rotation permit. The default "//&
2154  "for now is 200 m. This is only used if BOTTOMDRAGLAW is true.", &
2155  units="m", default=200.0, scale=us%m_to_Z)
2156 
2157  cs%IMax_decay = 0.0
2158  if (decay_length > 0.0) cs%IMax_decay = 1.0/decay_length
2159  call get_param(param_file, mdl, "BBL_MIXING_AS_MAX", cs%BBL_mixing_as_max, &
2160  "If true, take the maximum of the diffusivity from the "//&
2161  "BBL mixing and the other diffusivities. Otherwise, "//&
2162  "diffusivity from the BBL_mixing is simply added.", &
2163  default=.true.)
2164  call get_param(param_file, mdl, "USE_LOTW_BBL_DIFFUSIVITY", cs%use_LOTW_BBL_diffusivity, &
2165  "If true, uses a simple, imprecise but non-coordinate dependent, model "//&
2166  "of BBL mixing diffusivity based on Law of the Wall. Otherwise, uses "//&
2167  "the original BBL scheme.", default=.false.)
2168  if (cs%use_LOTW_BBL_diffusivity) then
2169  call get_param(param_file, mdl, "LOTW_BBL_USE_OMEGA", cs%LOTW_BBL_use_omega, &
2170  "If true, use the maximum of Omega and N for the TKE to diffusion "//&
2171  "calculation. Otherwise, N is N.", default=.true.)
2172  endif
2173  else
2174  cs%use_LOTW_BBL_diffusivity = .false. ! This parameterization depends on a u* from viscous BBL
2175  endif
2176  cs%id_Kd_BBL = register_diag_field('ocean_model', 'Kd_BBL', diag%axesTi, time, &
2177  'Bottom Boundary Layer Diffusivity', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
2178  call get_param(param_file, mdl, "SIMPLE_TKE_TO_KD", cs%simple_TKE_to_Kd, &
2179  "If true, uses a simple estimate of Kd/TKE that will "//&
2180  "work for arbitrary vertical coordinates. If false, "//&
2181  "calculates Kd/TKE and bounds based on exact energetics "//&
2182  "for an isopycnal layer-formulation.", default=.false.)
2183 
2184  ! set params related to the background mixing
2185  call bkgnd_mixing_init(time, g, gv, us, param_file, cs%diag, cs%bkgnd_mixing_csp)
2186 
2187  call get_param(param_file, mdl, "KV", cs%Kv, &
2188  "The background kinematic viscosity in the interior. "//&
2189  "The molecular value, ~1e-6 m2 s-1, may be used.", &
2190  units="m2 s-1", scale=us%m2_s_to_Z2_T, fail_if_missing=.true.)
2191 
2192  call get_param(param_file, mdl, "KD", cs%Kd, &
2193  "The background diapycnal diffusivity of density in the "//&
2194  "interior. Zero or the molecular value, ~1e-7 m2 s-1, "//&
2195  "may be used.", default=0.0, units="m2 s-1", scale=us%m2_s_to_Z2_T)
2196  call get_param(param_file, mdl, "KD_MIN", cs%Kd_min, &
2197  "The minimum diapycnal diffusivity.", &
2198  units="m2 s-1", default=0.01*cs%Kd*us%Z2_T_to_m2_s, scale=us%m2_s_to_Z2_T)
2199  call get_param(param_file, mdl, "KD_MAX", cs%Kd_max, &
2200  "The maximum permitted increment for the diapycnal "//&
2201  "diffusivity from TKE-based parameterizations, or a negative "//&
2202  "value for no limit.", units="m2 s-1", default=-1.0, scale=us%m2_s_to_Z2_T)
2203  if (cs%simple_TKE_to_Kd .and. cs%Kd_max<=0.) call mom_error(fatal, &
2204  "set_diffusivity_init: To use SIMPLE_TKE_TO_KD, KD_MAX must be set to >0.")
2205  call get_param(param_file, mdl, "KD_ADD", cs%Kd_add, &
2206  "A uniform diapycnal diffusivity that is added "//&
2207  "everywhere without any filtering or scaling.", &
2208  units="m2 s-1", default=0.0, scale=us%m2_s_to_Z2_T)
2209  if (cs%use_LOTW_BBL_diffusivity .and. cs%Kd_max<=0.) call mom_error(fatal, &
2210  "set_diffusivity_init: KD_MAX must be set (positive) when "// &
2211  "USE_LOTW_BBL_DIFFUSIVITY=True.")
2212  call get_param(param_file, mdl, "KD_SMOOTH", cs%Kd_smooth, &
2213  "A diapycnal diffusivity that is used to interpolate "//&
2214  "more sensible values of T & S into thin layers.", &
2215  units="m2 s-1", default=1.0e-6, scale=us%m2_s_to_Z2_T)
2216 
2217  call get_param(param_file, mdl, "DEBUG", cs%debug, &
2218  "If true, write out verbose debugging data.", &
2219  default=.false., debuggingparam=.true.)
2220 
2221  call get_param(param_file, mdl, "USER_CHANGE_DIFFUSIVITY", cs%user_change_diff, &
2222  "If true, call user-defined code to change the diffusivity.", default=.false.)
2223 
2224  call get_param(param_file, mdl, "DISSIPATION_MIN", cs%dissip_min, &
2225  "The minimum dissipation by which to determine a lower "//&
2226  "bound of Kd (a floor).", &
2227  units="W m-3", default=0.0, scale=us%W_m2_to_RZ3_T3*us%Z_to_m)
2228  call get_param(param_file, mdl, "DISSIPATION_N0", cs%dissip_N0, &
2229  "The intercept when N=0 of the N-dependent expression "//&
2230  "used to set a minimum dissipation by which to determine "//&
2231  "a lower bound of Kd (a floor): A in eps_min = A + B*N.", &
2232  units="W m-3", default=0.0, scale=us%W_m2_to_RZ3_T3*us%Z_to_m)
2233  call get_param(param_file, mdl, "DISSIPATION_N1", cs%dissip_N1, &
2234  "The coefficient multiplying N, following Gargett, used to "//&
2235  "set a minimum dissipation by which to determine a lower "//&
2236  "bound of Kd (a floor): B in eps_min = A + B*N", &
2237  units="J m-3", default=0.0, scale=us%W_m2_to_RZ3_T3*us%Z_to_m*us%s_to_T)
2238  call get_param(param_file, mdl, "DISSIPATION_KD_MIN", cs%dissip_Kd_min, &
2239  "The minimum vertical diffusivity applied as a floor.", &
2240  units="m2 s-1", default=0.0, scale=us%m2_s_to_Z2_T)
2241 
2242  cs%limit_dissipation = (cs%dissip_min>0.) .or. (cs%dissip_N1>0.) .or. &
2243  (cs%dissip_N0>0.) .or. (cs%dissip_Kd_min>0.)
2244  cs%dissip_N2 = 0.0
2245  if (cs%FluxRi_max > 0.0) &
2246  cs%dissip_N2 = cs%dissip_Kd_min * gv%Rho0 / cs%FluxRi_max
2247 
2248  cs%id_Kd_bkgnd = register_diag_field('ocean_model', 'Kd_bkgnd', diag%axesTi, time, &
2249  'Background diffusivity added by MOM_bkgnd_mixing module', 'm2/s', conversion=us%Z2_T_to_m2_s)
2250  cs%id_Kv_bkgnd = register_diag_field('ocean_model', 'Kv_bkgnd', diag%axesTi, time, &
2251  'Background viscosity added by MOM_bkgnd_mixing module', 'm2/s', conversion=us%Z2_T_to_m2_s)
2252 
2253  cs%id_Kd_layer = register_diag_field('ocean_model', 'Kd_layer', diag%axesTL, time, &
2254  'Diapycnal diffusivity of layers (as set)', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
2255 
2256  if (cs%use_tidal_mixing) then
2257  cs%id_Kd_Work = register_diag_field('ocean_model', 'Kd_Work', diag%axesTL, time, &
2258  'Work done by Diapycnal Mixing', 'W m-2', conversion=us%RZ3_T3_to_W_m2)
2259  cs%id_maxTKE = register_diag_field('ocean_model', 'maxTKE', diag%axesTL, time, &
2260  'Maximum layer TKE', 'm3 s-3', conversion=(us%Z_to_m**3*us%s_to_T**3))
2261  cs%id_TKE_to_Kd = register_diag_field('ocean_model', 'TKE_to_Kd', diag%axesTL, time, &
2262  'Convert TKE to Kd', 's2 m', conversion=us%Z2_T_to_m2_s*(us%m_to_Z**3*us%T_to_s**3))
2263  cs%id_N2 = register_diag_field('ocean_model', 'N2', diag%axesTi, time, &
2264  'Buoyancy frequency squared', 's-2', conversion=us%s_to_T**2, cmor_field_name='obvfsq', &
2265  cmor_long_name='Square of seawater buoyancy frequency', &
2266  cmor_standard_name='square_of_brunt_vaisala_frequency_in_sea_water')
2267  endif
2268 
2269  if (cs%user_change_diff) &
2270  cs%id_Kd_user = register_diag_field('ocean_model', 'Kd_user', diag%axesTi, time, &
2271  'User-specified Extra Diffusivity', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
2272 
2273  call get_param(param_file, mdl, "DOUBLE_DIFFUSION", cs%double_diffusion, &
2274  "If true, increase diffusivites for temperature or salinity based on the "//&
2275  "double-diffusive parameterization described in Large et al. (1994).", &
2276  default=.false.)
2277 
2278  if (cs%double_diffusion) then
2279  call get_param(param_file, mdl, "MAX_RRHO_SALT_FINGERS", cs%Max_Rrho_salt_fingers, &
2280  "Maximum density ratio for salt fingering regime.", &
2281  default=2.55, units="nondim")
2282  call get_param(param_file, mdl, "MAX_SALT_DIFF_SALT_FINGERS", cs%Max_salt_diff_salt_fingers, &
2283  "Maximum salt diffusivity for salt fingering regime.", &
2284  default=1.e-4, units="m2 s-1", scale=us%m2_s_to_Z2_T)
2285  call get_param(param_file, mdl, "KV_MOLECULAR", cs%Kv_molecular, &
2286  "Molecular viscosity for calculation of fluxes under double-diffusive "//&
2287  "convection.", default=1.5e-6, units="m2 s-1", scale=us%m2_s_to_Z2_T)
2288  ! The default molecular viscosity follows the CCSM4.0 and MOM4p1 defaults.
2289  endif ! old double-diffusion
2290 
2291  if (cs%user_change_diff) then
2292  call user_change_diff_init(time, g, gv, us, param_file, diag, cs%user_change_diff_CSp)
2293  endif
2294 
2295  call get_param(param_file, mdl, "BRYAN_LEWIS_DIFFUSIVITY", bryan_lewis_diffusivity, &
2296  "If true, use a Bryan & Lewis (JGR 1979) like tanh "//&
2297  "profile of background diapycnal diffusivity with depth. "//&
2298  "This is done via CVMix.", default=.false., do_not_log=.true.)
2299  if (cs%use_tidal_mixing .and. bryan_lewis_diffusivity) &
2300  call mom_error(fatal,"MOM_Set_Diffusivity: "// &
2301  "Bryan-Lewis and internal tidal dissipation are both enabled. Choose one.")
2302 
2303  cs%useKappaShear = kappa_shear_init(time, g, gv, us, param_file, cs%diag, cs%kappaShear_CSp)
2304  cs%Vertex_Shear = kappa_shear_at_vertex(param_file)
2305 
2306  if (cs%useKappaShear) &
2307  id_clock_kappashear = cpu_clock_id('(Ocean kappa_shear)', grain=clock_module)
2308 
2309  ! CVMix shear-driven mixing
2310  cs%use_CVMix_shear = cvmix_shear_init(time, g, gv, us, param_file, cs%diag, cs%CVMix_shear_csp)
2311 
2312  ! CVMix double diffusion mixing
2313  cs%use_CVMix_ddiff = cvmix_ddiff_init(time, g, gv, us, param_file, cs%diag, cs%CVMix_ddiff_csp)
2314  if (cs%use_CVMix_ddiff) &
2315  id_clock_cvmix_ddiff = cpu_clock_id('(Double diffusion via CVMix)', grain=clock_module)
2316 
2317  if (cs%double_diffusion .and. cs%use_CVMix_ddiff) then
2318  call mom_error(fatal, 'set_diffusivity_init: '// &
2319  'Multiple double-diffusion options selected (DOUBLE_DIFFUSION and'//&
2320  'USE_CVMIX_DDIFF), please disable all but one option to proceed.')
2321  endif
2322 
2323  if (cs%double_diffusion .or. cs%use_CVMix_ddiff) then
2324  cs%id_KT_extra = register_diag_field('ocean_model', 'KT_extra', diag%axesTi, time, &
2325  'Double-diffusive diffusivity for temperature', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
2326  cs%id_KS_extra = register_diag_field('ocean_model', 'KS_extra', diag%axesTi, time, &
2327  'Double-diffusive diffusivity for salinity', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
2328  endif
2329  if (cs%use_CVMix_ddiff) then
2330  cs%id_R_rho = register_diag_field('ocean_model', 'R_rho', diag%axesTi, time, &
2331  'Double-diffusion density ratio', 'nondim')
2332  endif
2333 
2334  if (present(halo_ts)) then
2335  halo_ts = 0
2336  if (cs%Vertex_Shear) halo_ts = 1
2337  endif
2338 
2339  if (present(double_diffuse)) then
2340  double_diffuse = (cs%double_diffusion .or. cs%use_CVMix_ddiff)
2341  endif
2342 
2343 end subroutine set_diffusivity_init
2344 
2345 !> Clear pointers and dealocate memory
2346 subroutine set_diffusivity_end(CS)
2347  type(set_diffusivity_cs), pointer :: CS !< Control structure for this module
2348 
2349  if (.not.associated(cs)) return
2350 
2351  call bkgnd_mixing_end(cs%bkgnd_mixing_csp)
2352 
2353  if (cs%use_tidal_mixing) call tidal_mixing_end(cs%tidal_mixing_CSp)
2354 
2355  if (cs%user_change_diff) call user_change_diff_end(cs%user_change_diff_CSp)
2356 
2357  if (cs%use_CVMix_shear) call cvmix_shear_end(cs%CVMix_shear_csp)
2358 
2359  if (cs%use_CVMix_ddiff) call cvmix_ddiff_end(cs%CVMix_ddiff_csp)
2360 
2361  if (associated(cs)) deallocate(cs)
2362 
2363 end subroutine set_diffusivity_end
2364 
2365 end module mom_set_diffusivity
Control structure including parameters for CVMix double diffusion.
This module implements boundary forcing for MOM6.
Control structure including parameters for CVMix interior shear schemes.
Interface to CVMix interior shear schemes.
This control structure has parameters for the MOM_internal_tides module.
Control structure with parameters for the tidal mixing module.
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...
Shear-dependent mixing following Jackson et al. 2008.
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Open boundary segment data structure.
Vertical viscosities, drag coefficients, and related fields.
Interface to background mixing schemes, including the Bryan and Lewis (1979) which is applied via CVM...
Wraps the MPP cpu clock functions.
This control structure contains parameters for MOM_set_diffusivity.
This module contains I/O framework code.
Definition: MOM_io.F90:2
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Make a diagnostic available for averaging or output.
Describes various unit conversion factors.
A module with intrinsic functions that are used by MOM but are not supported by some compilers...
This control structure holds the parameters that regulate shear mixing.
Control structure for user_change_diffusivity.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Subroutines that use the ray-tracing equations to propagate the internal tide energy density...
Routines for error handling and I/O management.
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Definition: MOM_EOS.F90:81
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Control structure including parameters for this module.
Describes the vertical ocean grid, including unit conversion factors.
Does full convective adjustment of unstable regions via a strong diffusivity.
Interface to vertical tidal mixing schemes including CVMix tidal mixing.
Interface to CVMix double diffusion scheme.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Calculate vertical diffusivity from all mixing processes.
Controls where open boundary conditions are applied.
Handy functions for manipulating strings.
Increments the diapycnal diffusivity in a specified band of latitudes and densities.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Read a data field from a file.
Definition: MOM_io.F90:74
This structure has memory for used in calculating diagnostics of diffusivity.
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.
A structure for creating arrays of pointers to 3D arrays.