21 use cvmix_background
, only : cvmix_init_bkgnd, cvmix_coeffs_bkgnd
23 implicit none ;
private 25 #include <MOM_memory.h> 27 public bkgnd_mixing_init
28 public bkgnd_mixing_end
29 public calculate_bkgnd_mixing
40 real :: bryan_lewis_c1
42 real :: bryan_lewis_c2
44 real :: bryan_lewis_c3
46 real :: bryan_lewis_c4
50 real :: bckgrnd_vdc_eq
52 real :: bckgrnd_vdc_psim
54 real :: bckgrnd_vdc_banda
64 real :: kd_tanh_lat_scale
70 logical :: kd_tanh_lat_fn
74 logical :: bryan_lewis_diffusivity
76 logical :: horiz_varying_background
78 logical :: henyey_igw_background
82 logical :: henyey_igw_background_new
96 logical :: bulkmixedlayer
97 logical :: kd_via_kdml_bug
104 character(len=40) :: bkgnd_scheme_str =
"none" 108 character(len=40) :: mdl =
"MOM_bkgnd_mixing" 113 subroutine bkgnd_mixing_init(Time, G, GV, US, param_file, diag, CS)
115 type(time_type),
intent(in) :: Time
120 type(
diag_ctrl),
target,
intent(inout) :: diag
126 real :: prandtl_bkgnd_comp
129 #include "version_variable.h" 131 if (
associated(cs))
then 132 call mom_error(warning,
"bkgnd_mixing_init called with an associated "// &
133 "control structure.")
140 "Adding static vertical background mixing coefficients")
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)
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.)
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)
160 cs%bulkmixedlayer = (gv%nkml > 0)
161 if (cs%bulkmixedlayer)
then 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.")
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.)
179 call get_param(param_file, mdl,
'DEBUG', cs%debug, default=.false., do_not_log=.true.)
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.)
188 if (cs%Bryan_Lewis_diffusivity)
then 189 call check_bkgnd_scheme(cs,
"BRYAN_LEWIS_DIFFUSIVITY")
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.)
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.)
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.)
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.)
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", &
214 if (cs%horiz_varying_background)
then 215 call check_bkgnd_scheme(cs,
"HORIZ_VARYING_BACKGROUND")
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)
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)
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)
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)
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)
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
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")
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")
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")
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))
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)
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.)
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)
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.")
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.", &
305 end subroutine bkgnd_mixing_init
308 subroutine calculate_bkgnd_mixing(h, tv, N2_lay, Kd_lay, Kd_int, Kv_bkgnd, j, G, GV, US, CS)
312 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
314 real,
dimension(SZI_(G),SZK_(GV)),
intent(in) :: N2_lay
316 real,
dimension(SZI_(G),SZK_(GV)),
intent(out) :: Kd_lay
318 real,
dimension(SZI_(G),SZK_(GV)+1),
intent(out) :: Kd_int
320 real,
dimension(SZI_(G),SZK_(GV)+1),
intent(out) :: Kv_bkgnd
322 integer,
intent(in) :: j
328 real,
dimension(SZK_(GV)+1) :: depth_int
329 real,
dimension(SZK_(GV)+1) :: Kd_col
330 real,
dimension(SZK_(GV)+1) :: Kv_col
331 real,
dimension(SZI_(G)) :: Kd_sfc
332 real,
dimension(SZI_(G)) :: depth
342 real :: bckgrnd_vdc_psin
343 real :: bckgrnd_vdc_psis
344 integer :: i, k, is, ie, js, je, nz
346 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
349 deg_to_rad = atan(1.0)/45.0
357 if (cs%Bryan_Lewis_diffusivity)
then 362 depth_int(k) = depth_int(k-1) + gv%H_to_m*h(i,j,k-1)
365 call cvmix_init_bkgnd(max_nlev=nz, &
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)
373 kd_col(:) = 0.0 ; kv_col(:) = 0.0
374 call cvmix_coeffs_bkgnd(mdiff_out=kv_col, tdiff_out=kd_col, nlev=nz, max_nlev=nz)
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)
382 kd_lay(i,k) = kd_lay(i,k) + 0.5 * us%m2_s_to_Z2_T * (kd_col(k) + kd_col(k+1))
386 elseif (cs%horiz_varying_background)
then 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
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
398 kd_int(i,1) = kd_int(i,1) + cs%bckgrnd_vdc1
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
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
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
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
428 do k=1,nz ;
do i=is,ie
429 kd_lay(i,k) = kd_int(i,1)
432 elseif (cs%Henyey_IGW_background_new)
then 433 i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0)
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)
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
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
453 if (cs%Henyey_IGW_background)
then 454 i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0)
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) )
460 elseif (cs%Kd_tanh_lat_fn)
then 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) ))
475 if ((.not.cs%bulkmixedlayer) .and. (cs%Kd /= cs%Kdml))
then 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 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)
488 kd_lay(i,k) = ((kd_sfc(i) - cs%Kdml) * i_hmix) * depth_c + (2.0*cs%Kdml - kd_sfc(i))
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)
494 kd_lay(i,k) = ((kd_sfc(i) - cs%Kdml) * i_hmix) * depth_c + (2.0*cs%Kdml - kd_sfc(i))
498 depth(i) = depth(i) + gv%H_to_Z*h(i,j,k)
502 do k=1,nz ;
do i=is,ie
503 kd_lay(i,k) = kd_sfc(i)
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
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
518 end subroutine calculate_bkgnd_mixing
523 logical function cvmix_bkgnd_is_used(param_file)
525 call get_param(param_file, mdl,
"USE_CVMix_BACKGROUND", cvmix_bkgnd_is_used, &
526 default=.false., do_not_log = .true.)
528 end function cvmix_bkgnd_is_used
532 subroutine check_bkgnd_scheme(CS, str)
534 character(len=*),
intent(in) :: str
537 if (trim(cs%bkgnd_scheme_str)==
"none")
then 538 cs%bkgnd_scheme_str = str
540 call mom_error(fatal,
"bkgnd_mixing_init: Cannot activate both "//trim(str)//
" and "//&
541 trim(cs%bkgnd_scheme_str)//
".")
547 subroutine bkgnd_mixing_end(CS)
551 if (.not.
associated(cs))
return 554 end subroutine bkgnd_mixing_end
This module implements boundary forcing for MOM6.
Ocean grid type. See mom_grid for details.
Calculates density of sea water from T, S and P.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Vertical viscosities, drag coefficients, and related fields.
Interface to background mixing schemes, including the Bryan and Lewis (1979) which is applied via CVM...
The MOM6 facility to parse input files for runtime 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...
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
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...
Provides subroutines for quantities specific to the equation of state.
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.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.