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 sfc_bkgnd_mixing (G, US, CS)
 Get surface vertical background diffusivities/viscosities. More...
 
subroutine, public calculate_bkgnd_mixing (h, tv, N2_lay, Kd_lay, Kv, 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. More...
 

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 576 of file MOM_bkgnd_mixing.F90.

576  type(bkgnd_mixing_cs), pointer :: CS !< Control structure for this module that
577  !! will be deallocated in this subroutine
578 
579  if (.not. associated(cs)) return
580 
581  deallocate(cs%kd_bkgnd)
582  deallocate(cs%kv_bkgnd)
583  deallocate(cs)
584 

Referenced by mom_set_diffusivity::set_diffusivity_end().

Here is the caller graph for this function:

◆ 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 119 of file MOM_bkgnd_mixing.F90.

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

References check_bkgnd_scheme(), mdl, and mom_error_handler::mom_error().

Here is the call graph for this function:

◆ calculate_bkgnd_mixing()

subroutine, public mom_bkgnd_mixing::calculate_bkgnd_mixing ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension( g %isd: g %ied, g %ke), intent(in)  N2_lay,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  Kd_lay,
real, dimension(:,:,:), pointer  Kv,
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]usA dimensional unit scaling type
[in]hLayer thickness [H ~> m or kg m-2].
[in]tvThermodynamics structure.
[in]n2_laysquared buoyancy frequency associated with layers [T-2 ~> s-2]
[in,out]kd_layDiapycnal diffusivity of each layer [Z2 T-1 ~> m2 s-1].
kvThe "slow" vertical viscosity at each interface (not layer!) [Z2 T-1 ~> m2 s-1]
[in]jMeridional grid index
csThe control structure returned by a previous call to bkgnd_mixing_init.

Definition at line 376 of file MOM_bkgnd_mixing.F90.

376 
377  type(ocean_grid_type), intent(in) :: G !< Grid structure.
378  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
379  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
380  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
381  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure.
382  real, dimension(SZI_(G),SZK_(G)), intent(in) :: N2_lay !< squared buoyancy frequency associated
383  !! with layers [T-2 ~> s-2]
384  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: Kd_lay !< Diapycnal diffusivity of each layer
385  !! [Z2 T-1 ~> m2 s-1].
386  real, dimension(:,:,:), pointer :: Kv !< The "slow" vertical viscosity at each interface
387  !! (not layer!) [Z2 T-1 ~> m2 s-1]
388  integer, intent(in) :: j !< Meridional grid index
389  type(bkgnd_mixing_cs), pointer :: CS !< The control structure returned by
390  !! a previous call to bkgnd_mixing_init.
391 
392  ! local variables
393  real, dimension(SZK_(G)+1) :: depth_int !< distance from surface of the interfaces [m]
394  real, dimension(SZK_(G)+1) :: Kd_col !< Diffusivities at the interfaces [m2 s-1]
395  real, dimension(SZK_(G)+1) :: Kv_col !< Viscosities at the interfaces [m2 s-1]
396  real, dimension(SZI_(G)) :: depth !< distance from surface of an interface [Z ~> m]
397  real :: depth_c !< depth of the center of a layer [Z ~> m]
398  real :: I_Hmix !< inverse of fixed mixed layer thickness [Z-1 ~> m-1]
399  real :: I_2Omega !< 1/(2 Omega) [T ~> s]
400  real :: N_2Omega ! The ratio of the stratification to the Earth's rotation rate [nondim]
401  real :: N02_N2 ! The ratio a reference stratification to the actual stratification [nondim]
402  real :: I_x30 !< 2/acos(2) = 1/(sin(30 deg) * acosh(1/sin(30 deg)))
403  real :: deg_to_rad !< factor converting degrees to radians, pi/180.
404  real :: abs_sin !< absolute value of sine of latitude [nondim]
405  real :: epsilon ! The minimum value of the sine of latitude [nondim]
406  real :: bckgrnd_vdc_psin !< PSI diffusivity in northern hemisphere [Z2 T-1 ~> m2 s-1]
407  real :: bckgrnd_vdc_psis !< PSI diffusivity in southern hemisphere [Z2 T-1 ~> m2 s-1]
408  integer :: i, k, is, ie, js, je, nz
409 
410  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
411 
412  ! set some parameters
413  deg_to_rad = atan(1.0)/45.0 ! = PI/180
414  epsilon = 1.e-10
415 
416  ! Set up the background diffusivity.
417  if (cs%Bryan_Lewis_diffusivity) then
418 
419  do i=is,ie
420  depth_int(1) = 0.0
421  do k=2,nz+1
422  depth_int(k) = depth_int(k-1) + gv%H_to_m*h(i,j,k-1)
423  enddo
424 
425  call cvmix_init_bkgnd(max_nlev=nz, &
426  zw = depth_int(:), & !< interface depths relative to the surface in m, must be positive.
427  bl1 = cs%Bryan_Lewis_c1, &
428  bl2 = cs%Bryan_Lewis_c2, &
429  bl3 = cs%Bryan_Lewis_c3, &
430  bl4 = cs%Bryan_Lewis_c4, &
431  prandtl = cs%prandtl_bkgnd)
432 
433  kd_col(:) = 0.0 ; kv_col(:) = 0.0 ! Is this line necessary?
434  call cvmix_coeffs_bkgnd(mdiff_out=kv_col, tdiff_out=kd_col, nlev=nz, max_nlev=nz)
435 
436  ! Update Kd and Kv.
437  do k=1,nz+1
438  cs%Kv_bkgnd(i,j,k) = us%m2_s_to_Z2_T*kv_col(k)
439  cs%Kd_bkgnd(i,j,k) = us%m2_s_to_Z2_T*kd_col(k)
440  enddo
441  do k=1,nz
442  kd_lay(i,j,k) = kd_lay(i,j,k) + 0.5 * us%m2_s_to_Z2_T * (kd_col(k) + kd_col(k+1))
443  enddo
444  enddo ! i loop
445 
446  elseif ((.not. cs%Bryan_Lewis_diffusivity) .and. (.not.cs%bulkmixedlayer) .and. &
447  (.not. cs%horiz_varying_background) .and. (cs%Kd /= cs%Kdml)) then
448  i_hmix = 1.0 / cs%Hmix
449  do i=is,ie ; depth(i) = 0.0 ; enddo
450  do k=1,nz ; do i=is,ie
451  depth_c = depth(i) + 0.5*gv%H_to_Z*h(i,j,k)
452  if (depth_c <= cs%Hmix) then ; cs%Kd_bkgnd(i,j,k) = cs%Kdml
453  elseif (depth_c >= 2.0*cs%Hmix) then ; cs%Kd_bkgnd(i,j,k) = cs%Kd_sfc(i,j)
454  else
455  kd_lay(i,j,k) = ((cs%Kd_sfc(i,j) - cs%Kdml) * i_hmix) * depth_c + &
456  (2.0*cs%Kdml - cs%Kd_sfc(i,j))
457  endif
458 
459  depth(i) = depth(i) + gv%H_to_Z*h(i,j,k)
460  enddo ; enddo
461 
462  elseif (cs%horiz_varying_background) then
463  !### Note that there are lots of hard-coded parameters (mostly latitudes and longitudes) here.
464  do i=is,ie
465  bckgrnd_vdc_psis = cs%bckgrnd_vdc_psim * exp(-(0.4*(g%geoLatT(i,j)+28.9))**2)
466  bckgrnd_vdc_psin = cs%bckgrnd_vdc_psim * exp(-(0.4*(g%geoLatT(i,j)-28.9))**2)
467  !### Add parentheses.
468  cs%Kd_bkgnd(i,j,:) = cs%bckgrnd_vdc_eq + bckgrnd_vdc_psin + bckgrnd_vdc_psis
469 
470  if (g%geoLatT(i,j) < -10.0) then
471  cs%Kd_bkgnd(i,j,:) = cs%Kd_bkgnd(i,j,:) + cs%bckgrnd_vdc1
472  elseif (g%geoLatT(i,j) <= 10.0) then
473  cs%Kd_bkgnd(i,j,:) = cs%Kd_bkgnd(i,j,:) + cs%bckgrnd_vdc1 * (g%geoLatT(i,j)/10.0)**2
474  else
475  cs%Kd_bkgnd(i,j,:) = cs%Kd_bkgnd(i,j,:) + cs%bckgrnd_vdc1
476  endif
477 
478  ! North Banda Sea
479  if ( (g%geoLatT(i,j) < -1.0) .and. (g%geoLatT(i,j) > -4.0) .and. &
480  ( mod(g%geoLonT(i,j)+360.0,360.0) > 103.0) .and. &
481  ( mod(g%geoLonT(i,j)+360.0,360.0) < 134.0) ) then
482  cs%Kd_bkgnd(i,j,:) = cs%bckgrnd_vdc_Banda
483  endif
484 
485  ! Middle Banda Sea
486  if ( (g%geoLatT(i,j) <= -4.0) .and. (g%geoLatT(i,j) > -7.0) .and. &
487  ( mod(g%geoLonT(i,j)+360.0,360.0) > 106.0) .and. &
488  ( mod(g%geoLonT(i,j)+360.0,360.0) < 140.0) ) then
489  cs%Kd_bkgnd(i,j,:) = cs%bckgrnd_vdc_Banda
490  endif
491 
492  ! South Banda Sea
493  if ( (g%geoLatT(i,j) <= -7.0) .and. (g%geoLatT(i,j) > -8.3) .and. &
494  ( mod(g%geoLonT(i,j)+360.0,360.0) > 111.0) .and. &
495  ( mod(g%geoLonT(i,j)+360.0,360.0) < 142.0) ) then
496  cs%Kd_bkgnd(i,j,:) = cs%bckgrnd_vdc_Banda
497  endif
498 
499  ! Compute kv_bkgnd
500  cs%kv_bkgnd(i,j,:) = cs%Kd_bkgnd(i,j,:) * cs%prandtl_bkgnd
501 
502  ! Update Kd (uniform profile; no interpolation needed)
503  kd_lay(i,j,:) = cs%Kd_bkgnd(i,j,1)
504 
505  enddo
506 
507  elseif (cs%Henyey_IGW_background_new) then
508  i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0) ! This is evaluated at 30 deg.
509  i_2omega = 0.5 / cs%omega
510  do k=1,nz ; do i=is,ie
511  abs_sin = max(epsilon, abs(sin(g%geoLatT(i,j)*deg_to_rad)))
512  n_2omega = max(abs_sin, sqrt(n2_lay(i,k))*i_2omega)
513  n02_n2 = (cs%N0_2Omega/n_2omega)**2
514  kd_lay(i,j,k) = max(cs%Kd_min, cs%Kd_sfc(i,j) * &
515  ((abs_sin * invcosh(n_2omega/abs_sin)) * i_x30)*n02_n2)
516  enddo ; enddo
517 
518  else
519  do k=1,nz ; do i=is,ie
520  kd_lay(i,j,k) = cs%Kd_sfc(i,j)
521  enddo ; enddo
522  endif
523 
524  ! Update CS%kd_bkgnd and CS%kv_bkgnd for diagnostic purposes
525  if (.not. (cs%Bryan_Lewis_diffusivity .or. cs%horiz_varying_background)) then
526  do i=is,ie
527  cs%kd_bkgnd(i,j,1) = 0.0; cs%kv_bkgnd(i,j,1) = 0.0
528  cs%kd_bkgnd(i,j,nz+1) = 0.0; cs%kv_bkgnd(i,j,nz+1) = 0.0
529  do k=2,nz
530  cs%Kd_bkgnd(i,j,k) = 0.5*(kd_lay(i,j,k-1) + kd_lay(i,j,k))
531  cs%Kv_bkgnd(i,j,k) = cs%Kd_bkgnd(i,j,k) * cs%prandtl_bkgnd
532  enddo
533  enddo
534  endif
535 
536  ! Update Kv
537  if (associated(kv)) then
538  do k=1,nz+1 ; do i=is,ie
539  kv(i,j,k) = kv(i,j,k) + cs%Kv_bkgnd(i,j,k)
540  enddo ; enddo
541  endif
542 
543  ! TODO: In both CS%Bryan_Lewis_diffusivity and CS%horiz_varying_background, KV and KD at surface
544  ! and bottom interfaces are set to be nonzero. Make sure this is not problematic.
545 

References mom_intrinsic_functions::invcosh().

Here is the call graph for this function:

◆ 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 561 of file MOM_bkgnd_mixing.F90.

561  type(bkgnd_mixing_cs), pointer :: CS !< Control structure
562  character(len=*), intent(in) :: str !< Background scheme identifier deducted from MOM_input
563  !! parameters
564 
565  if (trim(cs%bkgnd_scheme_str)=="none") then
566  cs%bkgnd_scheme_str = str
567  else
568  call mom_error(fatal, "set_diffusivity_init: Cannot activate both "//trim(str)//" and "//&
569  trim(cs%bkgnd_scheme_str)//".")
570  endif
571 

References mom_error_handler::mom_error().

Referenced by bkgnd_mixing_init().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 552 of file MOM_bkgnd_mixing.F90.

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

References mdl.

◆ sfc_bkgnd_mixing()

subroutine, public mom_bkgnd_mixing::sfc_bkgnd_mixing ( type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(bkgnd_mixing_cs), intent(inout), pointer  CS 
)

Get surface vertical background diffusivities/viscosities.

Parameters
[in]gGrid structure.
[in]usA dimensional unit scaling type
[in,out]csThe control structure returned by a previous call to bkgnd_mixing_init.

Definition at line 322 of file MOM_bkgnd_mixing.F90.

322 
323  type(ocean_grid_type), intent(in) :: G !< Grid structure.
324  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
325  type(bkgnd_mixing_cs), pointer, intent(inout) :: CS !< The control structure returned by
326  !! a previous call to bkgnd_mixing_init.
327  ! local variables
328  real :: I_x30 !< 2/acos(2) = 1/(sin(30 deg) * acosh(1/sin(30 deg)))
329  real :: deg_to_rad !< factor converting degrees to radians, pi/180.
330  real :: abs_sin !< absolute value of sine of latitude [nondim]
331  real :: epsilon
332  integer :: i, j, k, is, ie, js, je
333 
334  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
335 
336  ! set some parameters
337  deg_to_rad = atan(1.0)/45.0 ! = PI/180
338  epsilon = 1.e-10
339 
340 
341  if (.not. (cs%Bryan_Lewis_diffusivity .or. cs%horiz_varying_background)) then
342 !$OMP parallel do default(none) shared(is,ie,js,je,CS)
343  do j=js,je ; do i=is,ie
344  cs%Kd_sfc(i,j) = cs%Kd
345  enddo ; enddo
346  endif
347 
348  if (cs%Henyey_IGW_background) then
349  i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0) ! This is evaluated at 30 deg.
350 !$OMP parallel do default(none) &
351 !$OMP shared(is,ie,js,je,CS,G,deg_to_rad,epsilon,I_x30) &
352 !$OMP private(abs_sin)
353  do j=js,je ; do i=is,ie
354  abs_sin = abs(sin(g%geoLatT(i,j)*deg_to_rad))
355  cs%Kd_sfc(i,j) = max(cs%Kd_min, cs%Kd_sfc(i,j) * &
356  ((abs_sin * invcosh(cs%N0_2Omega/max(epsilon,abs_sin))) * i_x30) )
357  enddo ; enddo
358  elseif (cs%Kd_tanh_lat_fn) then
359 !$OMP parallel do default(none) shared(is,ie,js,je,CS,G)
360  do j=js,je ; do i=is,ie
361  ! The transition latitude and latitude range are hard-scaled here, since
362  ! this is not really intended for wide-spread use, but rather for
363  ! comparison with CM2M / CM2.1 settings.
364  cs%Kd_sfc(i,j) = max(cs%Kd_min, cs%Kd_sfc(i,j) * (1.0 + &
365  cs%Kd_tanh_lat_scale * 0.5*tanh((abs(g%geoLatT(i,j)) - 35.0)/5.0) ))
366  enddo ; enddo
367  endif
368 
369  if (cs%debug) call hchksum(cs%Kd_sfc,"After sfc_bkgnd_mixing: Kd_sfc",g%HI,haloshift=0, scale=us%Z2_T_to_m2_s)
370 

References mom_intrinsic_functions::invcosh().

Referenced by mom_set_diffusivity::set_diffusivity().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ mdl

character(len=40) mom_bkgnd_mixing::mdl = "MOM_bkgnd_mixing"

This module's name.

Definition at line 113 of file MOM_bkgnd_mixing.F90.

113 character(len=40) :: mdl = "MOM_bkgnd_mixing" !< This module's name.

Referenced by bkgnd_mixing_init(), and cvmix_bkgnd_is_used().