MOM6
mom_bkgnd_mixing Module Reference

Detailed Description

Interface to background mixing schemes, including the Bryan and Lewis (1979) which is applied via CVMix.

Data Types

type  bkgnd_mixing_cs
 Control structure including parameters for this module. More...
 

Functions/Subroutines

subroutine, public bkgnd_mixing_init (Time, G, GV, US, param_file, diag, CS)
 Initialize the background mixing routine. More...
 
subroutine, public calculate_bkgnd_mixing (h, tv, N2_lay, Kd_lay, Kd_int, Kv_bkgnd, j, G, GV, US, CS)
 Calculates the vertical background diffusivities/viscosities. More...
 
logical function cvmix_bkgnd_is_used (param_file)
 Reads the parameter "USE_CVMix_BACKGROUND" and returns state. This function allows other modules to know whether this parameterization will be used without needing to duplicate the log entry. More...
 
subroutine check_bkgnd_scheme (CS, str)
 Sets CSbkgnd_scheme_str to check whether multiple background diffusivity schemes are activated. The string is also for error/log messages. More...
 
subroutine, public bkgnd_mixing_end (CS)
 Clear pointers and dealocate memory. More...
 

Variables

character(len=40) mdl = "MOM_bkgnd_mixing"
 This module's name.
 

Function/Subroutine Documentation

◆ bkgnd_mixing_end()

subroutine, public mom_bkgnd_mixing::bkgnd_mixing_end ( type(bkgnd_mixing_cs), pointer  CS)

Clear pointers and dealocate memory.

Parameters
csControl structure for this module that will be deallocated in this subroutine

Definition at line 547 of file MOM_bkgnd_mixing.F90.

548  type(bkgnd_mixing_cs), pointer :: CS !< Control structure for this module that
549  !! will be deallocated in this subroutine
550 
551  if (.not. associated(cs)) return
552  deallocate(cs)
553 

◆ bkgnd_mixing_init()

subroutine, public mom_bkgnd_mixing::bkgnd_mixing_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(bkgnd_mixing_cs), pointer  CS 
)

Initialize the background mixing routine.

Parameters
[in]timeThe current time.
[in]gGrid structure.
[in]gvVertical grid structure.
[in]usA dimensional unit scaling type
[in]param_fileRun-time parameter file handle
[in,out]diagDiagnostics control structure.
csThis module's control structure.

Definition at line 113 of file MOM_bkgnd_mixing.F90.

114 
115  type(time_type), intent(in) :: Time !< The current time.
116  type(ocean_grid_type), intent(in) :: G !< Grid structure.
117  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
118  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
119  type(param_file_type), intent(in) :: param_file !< Run-time parameter file handle
120  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
121  type(bkgnd_mixing_cs), pointer :: CS !< This module's control structure.
122 
123  ! Local variables
124  real :: Kv ! The interior vertical viscosity [Z2 T-1 ~> m2 s-1] - read to set Prandtl
125  ! number unless it is provided as a parameter
126  real :: prandtl_bkgnd_comp ! Kv/CS%Kd. Gets compared with user-specified prandtl_bkgnd.
127 
128 ! This include declares and sets the variable "version".
129 #include "version_variable.h"
130 
131  if (associated(cs)) then
132  call mom_error(warning, "bkgnd_mixing_init called with an associated "// &
133  "control structure.")
134  return
135  endif
136  allocate(cs)
137 
138  ! Read parameters
139  call log_version(param_file, mdl, version, &
140  "Adding static vertical background mixing coefficients")
141 
142  call get_param(param_file, mdl, "KD", cs%Kd, &
143  "The background diapycnal diffusivity of density in the "//&
144  "interior. Zero or the molecular value, ~1e-7 m2 s-1, "//&
145  "may be used.", default=0.0, units="m2 s-1", scale=us%m2_s_to_Z2_T)
146 
147  call get_param(param_file, mdl, "KV", kv, &
148  "The background kinematic viscosity in the interior. "//&
149  "The molecular value, ~1e-6 m2 s-1, may be used.", &
150  units="m2 s-1", scale=us%m2_s_to_Z2_T, fail_if_missing=.true.)
151 
152  call get_param(param_file, mdl, "KD_MIN", cs%Kd_min, &
153  "The minimum diapycnal diffusivity.", &
154  units="m2 s-1", default=0.01*cs%Kd*us%Z2_T_to_m2_s, scale=us%m2_s_to_Z2_T)
155 
156  ! The following is needed to set one of the choices of vertical background mixing
157 
158  ! BULKMIXEDLAYER is not always defined (e.g., CM2G63L), so the following line by passes
159  ! the need to include BULKMIXEDLAYER in MOM_input
160  cs%bulkmixedlayer = (gv%nkml > 0)
161  if (cs%bulkmixedlayer) then
162  ! Check that Kdml is not set when using bulk mixed layer
163  call get_param(param_file, mdl, "KDML", cs%Kdml, default=-1.)
164  if (cs%Kdml>0.) call mom_error(fatal, &
165  "bkgnd_mixing_init: KDML cannot be set when using bulk mixed layer.")
166  cs%Kdml = cs%Kd ! This is not used with a bulk mixed layer, but also cannot be a NaN.
167  else
168  call get_param(param_file, mdl, "KDML", cs%Kdml, &
169  "If BULKMIXEDLAYER is false, KDML is the elevated "//&
170  "diapycnal diffusivity in the topmost HMIX of fluid. "//&
171  "KDML is only used if BULKMIXEDLAYER is false.", &
172  units="m2 s-1", default=cs%Kd*us%Z2_T_to_m2_s, scale=us%m2_s_to_Z2_T)
173  call get_param(param_file, mdl, "HMIX_FIXED", cs%Hmix, &
174  "The prescribed depth over which the near-surface "//&
175  "viscosity and diffusivity are elevated when the bulk "//&
176  "mixed layer is not used.", units="m", scale=us%m_to_Z, fail_if_missing=.true.)
177  endif
178 
179  call get_param(param_file, mdl, 'DEBUG', cs%debug, default=.false., do_not_log=.true.)
180 
181 ! call openParameterBlock(param_file,'MOM_BACKGROUND_MIXING')
182 
183  call get_param(param_file, mdl, "BRYAN_LEWIS_DIFFUSIVITY", cs%Bryan_Lewis_diffusivity, &
184  "If true, use a Bryan & Lewis (JGR 1979) like tanh "//&
185  "profile of background diapycnal diffusivity with depth. "//&
186  "This is done via CVMix.", default=.false.)
187 
188  if (cs%Bryan_Lewis_diffusivity) then
189  call check_bkgnd_scheme(cs, "BRYAN_LEWIS_DIFFUSIVITY")
190 
191  call get_param(param_file, mdl, "BRYAN_LEWIS_C1", cs%Bryan_Lewis_c1, &
192  "The vertical diffusivity values for Bryan-Lewis profile at |z|=D.", &
193  units="m2 s-1", fail_if_missing=.true.)
194 
195  call get_param(param_file, mdl, "BRYAN_LEWIS_C2", cs%Bryan_Lewis_c2, &
196  "The amplitude of variation in diffusivity for the Bryan-Lewis profile", &
197  units="m2 s-1", fail_if_missing=.true.)
198 
199  call get_param(param_file, mdl, "BRYAN_LEWIS_C3", cs%Bryan_Lewis_c3, &
200  "The inverse length scale for transition region in the Bryan-Lewis profile", &
201  units="m-1", fail_if_missing=.true.)
202 
203  call get_param(param_file, mdl, "BRYAN_LEWIS_C4", cs%Bryan_Lewis_c4, &
204  "The depth where diffusivity is BRYAN_LEWIS_C1 in the Bryan-Lewis profile",&
205  units="m", fail_if_missing=.true.)
206 
207  endif ! CS%Bryan_Lewis_diffusivity
208 
209  call get_param(param_file, mdl, "HORIZ_VARYING_BACKGROUND", cs%horiz_varying_background, &
210  "If true, apply vertically uniform, latitude-dependent background "//&
211  "diffusivity, as described in Danabasoglu et al., 2012", &
212  default=.false.)
213 
214  if (cs%horiz_varying_background) then
215  call check_bkgnd_scheme(cs, "HORIZ_VARYING_BACKGROUND")
216 
217  call get_param(param_file, mdl, "BCKGRND_VDC1", cs%bckgrnd_vdc1, &
218  "Background diffusivity (Ledwell) when HORIZ_VARYING_BACKGROUND=True", &
219  units="m2 s-1",default = 0.16e-04, scale=us%m2_s_to_Z2_T)
220 
221  call get_param(param_file, mdl, "BCKGRND_VDC_EQ", cs%bckgrnd_vdc_eq, &
222  "Equatorial diffusivity (Gregg) when HORIZ_VARYING_BACKGROUND=True", &
223  units="m2 s-1",default = 0.01e-04, scale=us%m2_s_to_Z2_T)
224 
225  call get_param(param_file, mdl, "BCKGRND_VDC_PSIM", cs%bckgrnd_vdc_psim, &
226  "Max. PSI induced diffusivity (MacKinnon) when HORIZ_VARYING_BACKGROUND=True", &
227  units="m2 s-1",default = 0.13e-4, scale=us%m2_s_to_Z2_T)
228 
229  call get_param(param_file, mdl, "BCKGRND_VDC_BAN", cs%bckgrnd_vdc_Banda, &
230  "Banda Sea diffusivity (Gordon) when HORIZ_VARYING_BACKGROUND=True", &
231  units="m2 s-1",default = 1.0e-4, scale=us%m2_s_to_Z2_T)
232  endif
233 
234  call get_param(param_file, mdl, "PRANDTL_BKGND", cs%prandtl_bkgnd, &
235  "Turbulent Prandtl number used to convert vertical "//&
236  "background diffusivities into viscosities.", &
237  units="nondim", default=1.0)
238 
239  if (cs%Bryan_Lewis_diffusivity .or. cs%horiz_varying_background) then
240  prandtl_bkgnd_comp = cs%prandtl_bkgnd
241  if (cs%Kd /= 0.0) prandtl_bkgnd_comp = kv / cs%Kd
242 
243  if ( abs(cs%prandtl_bkgnd - prandtl_bkgnd_comp)>1.e-14) then
244  call mom_error(fatal, "bkgnd_mixing_init: The provided KD, KV and PRANDTL_BKGND values "//&
245  "are incompatible. The following must hold: KD*PRANDTL_BKGND==KV")
246  endif
247  endif
248 
249  call get_param(param_file, mdl, "HENYEY_IGW_BACKGROUND", cs%Henyey_IGW_background, &
250  "If true, use a latitude-dependent scaling for the near "//&
251  "surface background diffusivity, as described in "//&
252  "Harrison & Hallberg, JPO 2008.", default=.false.)
253  if (cs%Henyey_IGW_background) call check_bkgnd_scheme(cs, "HENYEY_IGW_BACKGROUND")
254 
255 
256  call get_param(param_file, mdl, "HENYEY_IGW_BACKGROUND_NEW", cs%Henyey_IGW_background_new, &
257  "If true, use a better latitude-dependent scaling for the "//&
258  "background diffusivity, as described in "//&
259  "Harrison & Hallberg, JPO 2008.", default=.false.)
260  if (cs%Henyey_IGW_background_new) call check_bkgnd_scheme(cs, "HENYEY_IGW_BACKGROUND_NEW")
261 
262  if (cs%Kd>0.0 .and. (trim(cs%bkgnd_scheme_str)=="BRYAN_LEWIS_DIFFUSIVITY" .or.&
263  trim(cs%bkgnd_scheme_str)=="HORIZ_VARYING_BACKGROUND" )) then
264  call mom_error(warning, "bkgnd_mixing_init: a nonzero constant background "//&
265  "diffusivity (KD) is specified along with "//trim(cs%bkgnd_scheme_str))
266  endif
267 
268  if (cs%Henyey_IGW_background) then
269  call get_param(param_file, mdl, "HENYEY_N0_2OMEGA", cs%N0_2Omega, &
270  "The ratio of the typical Buoyancy frequency to twice "//&
271  "the Earth's rotation period, used with the Henyey "//&
272  "scaling from the mixing.", units="nondim", default=20.0)
273  call get_param(param_file, mdl, "OMEGA", cs%omega, &
274  "The rotation rate of the earth.", units="s-1", &
275  default=7.2921e-5, scale=us%T_to_s)
276  endif
277 
278  call get_param(param_file, mdl, "KD_TANH_LAT_FN", cs%Kd_tanh_lat_fn, &
279  "If true, use a tanh dependence of Kd_sfc on latitude, "//&
280  "like CM2.1/CM2M. There is no physical justification "//&
281  "for this form, and it can not be used with "//&
282  "HENYEY_IGW_BACKGROUND.", default=.false.)
283 
284  if (cs%Kd_tanh_lat_fn) &
285  call get_param(param_file, mdl, "KD_TANH_LAT_SCALE", cs%Kd_tanh_lat_scale, &
286  "A nondimensional scaling for the range ofdiffusivities "//&
287  "with KD_TANH_LAT_FN. Valid values are in the range of "//&
288  "-2 to 2; 0.4 reproduces CM2M.", units="nondim", default=0.0)
289 
290  if (cs%Henyey_IGW_background .and. cs%Kd_tanh_lat_fn) call mom_error(fatal, &
291  "MOM_bkgnd_mixing: KD_TANH_LAT_FN can not be used with HENYEY_IGW_BACKGROUND.")
292 
293  cs%Kd_via_Kdml_bug = .false.
294  if ((cs%Kd /= cs%Kdml) .and. .not.(cs%Kd_tanh_lat_fn .or. cs%bulkmixedlayer .or. &
295  cs%Henyey_IGW_background .or. cs%Henyey_IGW_background_new .or. &
296  cs%horiz_varying_background .or. cs%Bryan_Lewis_diffusivity)) then
297  call get_param(param_file, mdl, "KD_BACKGROUND_VIA_KDML_BUG", cs%Kd_via_Kdml_bug, &
298  "If true and KDML /= KD and several other conditions apply, the background "//&
299  "diffusivity is set incorrectly using a bug that was introduced in March, 2018.", &
300  default=.true.) ! The default should be changed to false and this parameter obsoleted.
301  endif
302 
303 ! call closeParameterBlock(param_file)
304 

◆ calculate_bkgnd_mixing()

subroutine, public mom_bkgnd_mixing::calculate_bkgnd_mixing ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension( g %isd: g %ied, gv %ke), intent(in)  N2_lay,
real, dimension( g %isd: g %ied, gv %ke), intent(out)  Kd_lay,
real, dimension( g %isd: g %ied, gv %ke+1), intent(out)  Kd_int,
real, dimension( g %isd: g %ied, gv %ke+1), intent(out)  Kv_bkgnd,
integer, intent(in)  j,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(bkgnd_mixing_cs), pointer  CS 
)

Calculates the vertical background diffusivities/viscosities.

Parameters
[in]gGrid structure.
[in]gvVertical grid structure.
[in]hLayer thickness [H ~> m or kg m-2].
[in]tvThermodynamics structure.
[in]n2_laysquared buoyancy frequency associated with layers [T-2 ~> s-2]
[out]kd_layThe background diapycnal diffusivity of each layer [Z2 T-1 ~> m2 s-1].
[out]kd_intThe background diapycnal diffusivity of each interface [Z2 T-1 ~> m2 s-1].
[out]kv_bkgndThe background vertical viscosity at each interface [Z2 T-1 ~> m2 s-1]
[in]jMeridional grid index
[in]usA dimensional unit scaling type
csThe control structure returned by a previous call to bkgnd_mixing_init.

Definition at line 308 of file MOM_bkgnd_mixing.F90.

309 
310  type(ocean_grid_type), intent(in) :: G !< Grid structure.
311  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
312  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
313  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure.
314  real, dimension(SZI_(G),SZK_(GV)), intent(in) :: N2_lay !< squared buoyancy frequency associated
315  !! with layers [T-2 ~> s-2]
316  real, dimension(SZI_(G),SZK_(GV)), intent(out) :: Kd_lay !< The background diapycnal diffusivity
317  !! of each layer [Z2 T-1 ~> m2 s-1].
318  real, dimension(SZI_(G),SZK_(GV)+1), intent(out) :: Kd_int !< The background diapycnal diffusivity
319  !! of each interface [Z2 T-1 ~> m2 s-1].
320  real, dimension(SZI_(G),SZK_(GV)+1), intent(out) :: Kv_bkgnd !< The background vertical viscosity at
321  !! each interface [Z2 T-1 ~> m2 s-1]
322  integer, intent(in) :: j !< Meridional grid index
323  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
324  type(bkgnd_mixing_cs), pointer :: CS !< The control structure returned by
325  !! a previous call to bkgnd_mixing_init.
326 
327  ! local variables
328  real, dimension(SZK_(GV)+1) :: depth_int !< Distance from surface of the interfaces [m]
329  real, dimension(SZK_(GV)+1) :: Kd_col !< Diffusivities at the interfaces [m2 s-1]
330  real, dimension(SZK_(GV)+1) :: Kv_col !< Viscosities at the interfaces [m2 s-1]
331  real, dimension(SZI_(G)) :: Kd_sfc !< Surface value of the diffusivity [Z2 T-1 ~> m2 s-1]
332  real, dimension(SZI_(G)) :: depth !< Distance from surface of an interface [Z ~> m]
333  real :: depth_c !< depth of the center of a layer [Z ~> m]
334  real :: I_Hmix !< inverse of fixed mixed layer thickness [Z-1 ~> m-1]
335  real :: I_2Omega !< 1/(2 Omega) [T ~> s]
336  real :: N_2Omega ! The ratio of the stratification to the Earth's rotation rate [nondim]
337  real :: N02_N2 ! The ratio a reference stratification to the actual stratification [nondim]
338  real :: I_x30 !< 2/acos(2) = 1/(sin(30 deg) * acosh(1/sin(30 deg)))
339  real :: deg_to_rad !< factor converting degrees to radians, pi/180.
340  real :: abs_sinlat !< absolute value of sine of latitude [nondim]
341  real :: min_sinlat ! The minimum value of the sine of latitude [nondim]
342  real :: bckgrnd_vdc_psin !< PSI diffusivity in northern hemisphere [Z2 T-1 ~> m2 s-1]
343  real :: bckgrnd_vdc_psis !< PSI diffusivity in southern hemisphere [Z2 T-1 ~> m2 s-1]
344  integer :: i, k, is, ie, js, je, nz
345 
346  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
347 
348  ! set some parameters
349  deg_to_rad = atan(1.0)/45.0 ! = PI/180
350  min_sinlat = 1.e-10
351 
352  ! Start with a constant value that may be replaced below.
353  kd_lay(:,:) = cs%Kd
354  kv_bkgnd(:,:) = 0.0
355 
356  ! Set up the background diffusivity.
357  if (cs%Bryan_Lewis_diffusivity) then
358 
359  do i=is,ie
360  depth_int(1) = 0.0
361  do k=2,nz+1
362  depth_int(k) = depth_int(k-1) + gv%H_to_m*h(i,j,k-1)
363  enddo
364 
365  call cvmix_init_bkgnd(max_nlev=nz, &
366  zw = depth_int(:), & !< interface depths relative to the surface in m, must be positive.
367  bl1 = cs%Bryan_Lewis_c1, &
368  bl2 = cs%Bryan_Lewis_c2, &
369  bl3 = cs%Bryan_Lewis_c3, &
370  bl4 = cs%Bryan_Lewis_c4, &
371  prandtl = cs%prandtl_bkgnd)
372 
373  kd_col(:) = 0.0 ; kv_col(:) = 0.0 ! Is this line necessary?
374  call cvmix_coeffs_bkgnd(mdiff_out=kv_col, tdiff_out=kd_col, nlev=nz, max_nlev=nz)
375 
376  ! Update Kd and Kv.
377  do k=1,nz+1
378  kv_bkgnd(i,k) = us%m2_s_to_Z2_T*kv_col(k)
379  kd_int(i,k) = us%m2_s_to_Z2_T*kd_col(k)
380  enddo
381  do k=1,nz
382  kd_lay(i,k) = kd_lay(i,k) + 0.5 * us%m2_s_to_Z2_T * (kd_col(k) + kd_col(k+1))
383  enddo
384  enddo ! i loop
385 
386  elseif (cs%horiz_varying_background) then
387  !### Note that there are lots of hard-coded parameters (mostly latitudes and longitudes) here.
388  do i=is,ie
389  bckgrnd_vdc_psis = cs%bckgrnd_vdc_psim * exp(-(0.4*(g%geoLatT(i,j)+28.9))**2)
390  bckgrnd_vdc_psin = cs%bckgrnd_vdc_psim * exp(-(0.4*(g%geoLatT(i,j)-28.9))**2)
391  kd_int(i,1) = (cs%bckgrnd_vdc_eq + bckgrnd_vdc_psin) + bckgrnd_vdc_psis
392 
393  if (g%geoLatT(i,j) < -10.0) then
394  kd_int(i,1) = kd_int(i,1) + cs%bckgrnd_vdc1
395  elseif (g%geoLatT(i,j) <= 10.0) then
396  kd_int(i,1) = kd_int(i,1) + cs%bckgrnd_vdc1 * (g%geoLatT(i,j)/10.0)**2
397  else
398  kd_int(i,1) = kd_int(i,1) + cs%bckgrnd_vdc1
399  endif
400 
401  ! North Banda Sea
402  if ( (g%geoLatT(i,j) < -1.0) .and. (g%geoLatT(i,j) > -4.0) .and. &
403  ( mod(g%geoLonT(i,j)+360.0,360.0) > 103.0) .and. &
404  ( mod(g%geoLonT(i,j)+360.0,360.0) < 134.0) ) then
405  kd_int(i,1) = cs%bckgrnd_vdc_Banda
406  endif
407 
408  ! Middle Banda Sea
409  if ( (g%geoLatT(i,j) <= -4.0) .and. (g%geoLatT(i,j) > -7.0) .and. &
410  ( mod(g%geoLonT(i,j)+360.0,360.0) > 106.0) .and. &
411  ( mod(g%geoLonT(i,j)+360.0,360.0) < 140.0) ) then
412  kd_int(i,1) = cs%bckgrnd_vdc_Banda
413  endif
414 
415  ! South Banda Sea
416  if ( (g%geoLatT(i,j) <= -7.0) .and. (g%geoLatT(i,j) > -8.3) .and. &
417  ( mod(g%geoLonT(i,j)+360.0,360.0) > 111.0) .and. &
418  ( mod(g%geoLonT(i,j)+360.0,360.0) < 142.0) ) then
419  kd_int(i,1) = cs%bckgrnd_vdc_Banda
420  endif
421 
422  enddo
423  ! Update interior values of Kd and Kv (uniform profile; no interpolation needed)
424  do k=1,nz+1 ; do i=is,ie
425  kd_int(i,k) = kd_int(i,1)
426  kv_bkgnd(i,k) = kd_int(i,1) * cs%prandtl_bkgnd
427  enddo ; enddo
428  do k=1,nz ; do i=is,ie
429  kd_lay(i,k) = kd_int(i,1)
430  enddo ; enddo
431 
432  elseif (cs%Henyey_IGW_background_new) then
433  i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0) ! This is evaluated at 30 deg.
434  i_2omega = 0.5 / cs%omega
435  do k=1,nz ; do i=is,ie
436  abs_sinlat = max(min_sinlat, abs(sin(g%geoLatT(i,j)*deg_to_rad)))
437  n_2omega = max(abs_sinlat, sqrt(n2_lay(i,k))*i_2omega)
438  n02_n2 = (cs%N0_2Omega/n_2omega)**2
439  kd_lay(i,k) = max(cs%Kd_min, cs%Kd * &
440  ((abs_sinlat * invcosh(n_2omega/abs_sinlat)) * i_x30)*n02_n2)
441  enddo ; enddo
442  ! Update Kd_int and Kv_bkgnd, based on Kd_lay. These might be just used for diagnostic purposes.
443  do i=is,ie
444  kd_int(i,1) = 0.0; kv_bkgnd(i,1) = 0.0
445  kd_int(i,nz+1) = 0.0; kv_bkgnd(i,nz+1) = 0.0
446  enddo
447  do k=2,nz ; do i=is,ie
448  kd_int(i,k) = 0.5*(kd_lay(i,k-1) + kd_lay(i,k))
449  kv_bkgnd(i,k) = kd_int(i,k) * cs%prandtl_bkgnd
450  enddo ; enddo
451  else
452  ! Set a potentially spatially varying surface value of diffusivity.
453  if (cs%Henyey_IGW_background) then
454  i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0) ! This is evaluated at 30 deg.
455  do i=is,ie
456  abs_sinlat = abs(sin(g%geoLatT(i,j)*deg_to_rad))
457  kd_sfc(i) = max(cs%Kd_min, cs%Kd * &
458  ((abs_sinlat * invcosh(cs%N0_2Omega / max(min_sinlat, abs_sinlat))) * i_x30) )
459  enddo
460  elseif (cs%Kd_tanh_lat_fn) then
461  do i=is,ie
462  ! The transition latitude and latitude range are hard-scaled here, since
463  ! this is not really intended for wide-spread use, but rather for
464  ! comparison with CM2M / CM2.1 settings.
465  kd_sfc(i) = max(cs%Kd_min, cs%Kd * (1.0 + &
466  cs%Kd_tanh_lat_scale * 0.5*tanh((abs(g%geoLatT(i,j)) - 35.0)/5.0) ))
467  enddo
468  else ! Use a spatially constant surface value.
469  do i=is,ie
470  kd_sfc(i) = cs%Kd
471  enddo
472  endif
473 
474  ! Now set background diffusivies based on these surface values, possibly with vertical structure.
475  if ((.not.cs%bulkmixedlayer) .and. (cs%Kd /= cs%Kdml)) then
476  ! This is a crude way to put in a diffusive boundary layer without an explicit boundary
477  ! layer turbulence scheme. It should not be used for any realistic ocean models.
478  i_hmix = 1.0 / (cs%Hmix + gv%H_subroundoff*gv%H_to_Z)
479  do i=is,ie ; depth(i) = 0.0 ; enddo
480  do k=1,nz ; do i=is,ie
481  depth_c = depth(i) + 0.5*gv%H_to_Z*h(i,j,k)
482  if (cs%Kd_via_Kdml_bug) then
483  ! These two lines should update Kd_lay, not Kd_int. They were correctly working on the
484  ! same variables until MOM6 commit 7a818716 (PR#750), which was added on March 26, 2018.
485  if (depth_c <= cs%Hmix) then ; kd_int(i,k) = cs%Kdml
486  elseif (depth_c >= 2.0*cs%Hmix) then ; kd_int(i,k) = kd_sfc(i)
487  else
488  kd_lay(i,k) = ((kd_sfc(i) - cs%Kdml) * i_hmix) * depth_c + (2.0*cs%Kdml - kd_sfc(i))
489  endif
490  else
491  if (depth_c <= cs%Hmix) then ; kd_lay(i,k) = cs%Kdml
492  elseif (depth_c >= 2.0*cs%Hmix) then ; kd_lay(i,k) = kd_sfc(i)
493  else
494  kd_lay(i,k) = ((kd_sfc(i) - cs%Kdml) * i_hmix) * depth_c + (2.0*cs%Kdml - kd_sfc(i))
495  endif
496  endif
497 
498  depth(i) = depth(i) + gv%H_to_Z*h(i,j,k)
499  enddo ; enddo
500 
501  else ! There is no vertical structure to the background diffusivity.
502  do k=1,nz ; do i=is,ie
503  kd_lay(i,k) = kd_sfc(i)
504  enddo ; enddo
505  endif
506 
507  ! Update Kd_int and Kv_bkgnd, based on Kd_lay. These might be just used for diagnostic purposes.
508  do i=is,ie
509  kd_int(i,1) = 0.0; kv_bkgnd(i,1) = 0.0
510  kd_int(i,nz+1) = 0.0; kv_bkgnd(i,nz+1) = 0.0
511  enddo
512  do k=2,nz ; do i=is,ie
513  kd_int(i,k) = 0.5*(kd_lay(i,k-1) + kd_lay(i,k))
514  kv_bkgnd(i,k) = kd_int(i,k) * cs%prandtl_bkgnd
515  enddo ; enddo
516  endif
517 

◆ check_bkgnd_scheme()

subroutine mom_bkgnd_mixing::check_bkgnd_scheme ( type(bkgnd_mixing_cs), pointer  CS,
character(len=*), intent(in)  str 
)
private

Sets CSbkgnd_scheme_str to check whether multiple background diffusivity schemes are activated. The string is also for error/log messages.

Parameters
csControl structure
[in]strBackground scheme identifier deducted from MOM_input parameters

Definition at line 532 of file MOM_bkgnd_mixing.F90.

533  type(bkgnd_mixing_cs), pointer :: CS !< Control structure
534  character(len=*), intent(in) :: str !< Background scheme identifier deducted from MOM_input
535  !! parameters
536 
537  if (trim(cs%bkgnd_scheme_str)=="none") then
538  cs%bkgnd_scheme_str = str
539  else
540  call mom_error(fatal, "bkgnd_mixing_init: Cannot activate both "//trim(str)//" and "//&
541  trim(cs%bkgnd_scheme_str)//".")
542  endif
543 

◆ cvmix_bkgnd_is_used()

logical function mom_bkgnd_mixing::cvmix_bkgnd_is_used ( type(param_file_type), intent(in)  param_file)
private

Reads the parameter "USE_CVMix_BACKGROUND" and returns state. This function allows other modules to know whether this parameterization will be used without needing to duplicate the log entry.

Parameters
[in]param_fileA structure to parse for run-time parameters

Definition at line 523 of file MOM_bkgnd_mixing.F90.

524  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
525  call get_param(param_file, mdl, "USE_CVMix_BACKGROUND", cvmix_bkgnd_is_used, &
526  default=.false., do_not_log = .true.)
527