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
16 use mom_debugging,
only : hchksum, uvchksum, bchksum, hchksum_pair
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
43 implicit none ;
private
45 #include <MOM_memory.h>
47 public set_diffusivity
49 public set_diffusivity_init
50 public set_diffusivity_end
61 logical :: bulkmixedlayer
67 logical :: bottomdraglaw
69 logical :: bbl_mixing_as_max
73 logical :: use_lotw_bbl_diffusivity
74 logical :: lotw_bbl_use_omega
91 logical :: limit_dissipation
100 logical :: ml_radiation
113 real :: ml_rad_kd_max
115 real :: ml_rad_efold_coeff
119 logical :: ml_rad_bug
122 logical :: ml_rad_tke_decay
131 logical :: ml_use_omega
134 real :: ml_omega_frac
137 logical :: user_change_diff
138 logical :: usekappashear
140 logical :: vertex_shear
142 logical :: use_cvmix_shear
144 logical :: double_diffusion
145 logical :: use_cvmix_ddiff
146 logical :: use_tidal_mixing
147 logical :: simple_tke_to_kd
149 real :: max_rrho_salt_fingers
150 real :: max_salt_diff_salt_fingers
153 logical :: answers_2018
157 character(len=200) :: inputdir
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
177 real,
pointer,
dimension(:,:,:) :: &
183 kd_bkgnd => null(), &
184 kv_bkgnd => null(), &
185 kt_extra => null(), &
186 ks_extra => null(), &
188 real,
pointer,
dimension(:,:,:) :: tke_to_kd => null()
196 integer :: id_clock_kappashear, id_clock_cvmix_ddiff
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)
216 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
218 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
220 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
222 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
224 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
228 type(
forcing),
intent(in) :: fluxes
229 type(optics_type),
pointer :: optics
233 real,
intent(in) :: dt
235 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
236 optional,
intent(out) :: kd_lay
237 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
238 optional,
intent(out) :: kd_int
239 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
240 optional,
intent(out) :: kd_extra_t
243 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
244 optional,
intent(out) :: kd_extra_s
249 real,
dimension(SZI_(G)) :: &
254 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
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]
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]
278 logical :: tke_to_kd_used
279 integer :: kb(szi_(g))
281 logical :: showcalltree
283 integer :: i, j, k, is, ie, js, je, nz, isd, ied, jsd, jed
285 real :: kappa_dt_fill
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")
292 if (.not.
associated(cs))
call mom_error(fatal,
"set_diffusivity: "//&
293 "Module must be initialized before it is used.")
295 if (cs%answers_2018)
then
297 kappa_dt_fill = us%m_to_Z**2 * 1.e-3 * 7200.
299 kappa_dt_fill = cs%Kd_smooth * dt
301 omega2 = cs%omega * cs%omega
303 use_eos =
associated(tv%eqn_of_state)
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.")
310 tke_to_kd_used = (cs%use_tidal_mixing .or. cs%ML_radiation .or. &
311 (cs%bottomdraglaw .and. .not.cs%use_LOTW_BBL_diffusivity))
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
322 if (cs%id_N2 > 0)
then
323 allocate(dd%N2_3d(isd:ied,jsd:jed,nz+1)) ; dd%N2_3d(:,:,:) = 0.0
325 if (cs%id_Kd_user > 0)
then
326 allocate(dd%Kd_user(isd:ied,jsd:jed,nz+1)) ; dd%Kd_user(:,:,:) = 0.0
328 if (cs%id_Kd_work > 0)
then
329 allocate(dd%Kd_work(isd:ied,jsd:jed,nz)) ; dd%Kd_work(:,:,:) = 0.0
331 if (cs%id_maxTKE > 0)
then
332 allocate(dd%maxTKE(isd:ied,jsd:jed,nz)) ; dd%maxTKE(:,:,:) = 0.0
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
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
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
343 if (cs%id_R_rho > 0)
then
344 allocate(dd%drho_rat(isd:ied,jsd:jed,nz+1)) ; dd%drho_rat(:,:,:) = 0.0
346 if (cs%id_Kd_BBL > 0)
then
347 allocate(dd%Kd_BBL(isd:ied,jsd:jed,nz+1)) ; dd%Kd_BBL(:,:,:) = 0.0
350 if (cs%id_Kd_bkgnd > 0)
then
351 allocate(dd%Kd_bkgnd(isd:ied,jsd:jed,nz+1)) ; dd%Kd_bkgnd(:,:,:) = 0.
353 if (cs%id_Kv_bkgnd > 0)
then
354 allocate(dd%Kv_bkgnd(isd:ied,jsd:jed,nz+1)) ; dd%Kv_bkgnd(:,:,:) = 0.
358 call setup_tidal_diagnostics(g, cs%tidal_mixing_CSp)
360 if (cs%useKappaShear)
then
362 call hchksum_pair(
"before calc_KS [uv]_h", u_h, v_h, g%HI, scale=us%L_T_to_m_s)
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)
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
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)
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)
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)
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
391 call calculate_cvmix_shear(u_h, v_h, h, tv, visc%Kd_shear, visc%Kv_shear, g, gv, us, cs%CVMix_shear_CSp)
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)
396 elseif (
associated(visc%Kv_shear))
then
397 visc%Kv_shear(:,:,:) = 0.0
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)
407 call vert_fill_ts(h, tv%T, tv%S, kappa_dt_fill, t_f, s_f, g, gv, larger_h_denom=.true.)
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)
424 call find_n2(h, tv, t_f, s_f, fluxes, j, g, gv, us, cs, drho_int, n2_lay, n2_int, n2_bot)
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
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)
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
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)
448 do k=2,nz ;
do i=is,ie
449 if (ks_extra(i,k) > kt_extra(i,k))
then
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
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
460 kd_extra_t(i,j,k) = 0.0
461 kd_extra_s(i,j,k) = 0.0
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
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
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)
481 call compute_ddiff_coeffs(h, tv, g, gv, us, j, kd_extra_t, kd_extra_s, cs%CVMix_ddiff_csp)
483 call cpu_clock_end(id_clock_cvmix_ddiff)
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
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))
504 kd_int_2d(i,1) = visc%Kd_shear(i,j,1)
505 kd_int_2d(i,nz+1) = 0.0
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))
512 if (
present(kd_int))
then
514 kd_int_2d(i,1) = kd_lay_2d(i,1) ; kd_int_2d(i,nz+1) = 0.0
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))
522 if (
present(kd_int))
then
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)
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)
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)
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)
543 if (cs%limit_dissipation)
then
549 do k=2,nz ;
do i=is,ie
550 dissip = max( cs%dissip_min, &
551 cs%dissip_N0 + cs%dissip_N1 * sqrt(n2_int(i,k)), &
552 cs%dissip_N2 * n2_int(i,k))
553 kd_int_2d(i,k) = max(kd_int_2d(i,k) , &
554 dissip * (cs%FluxRi_max / (gv%Rho0 * (n2_int(i,k) + omega2))))
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
564 do k=1,nz+1 ;
do i=is,ie
565 kd_int(i,j,k) = kd_int_2d(i,k)
571 if (cs%ML_radiation) &
572 call add_mlrad_diffusivity(h, fluxes, j, g, gv, us, cs, tke_to_kd, kd_lay_2d)
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)
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)
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)
592 if (cs%limit_dissipation)
then
598 do k=2,nz-1 ;
do i=is,ie
599 dissip = max( cs%dissip_min, &
600 cs%dissip_N0 + cs%dissip_N1 * sqrt(n2_lay(i,k)), &
601 cs%dissip_N2 * n2_lay(i,k))
602 kd_lay_2d(i,k) = max(kd_lay_2d(i,k) , &
603 dissip * (cs%FluxRi_max / (gv%Rho0 * (n2_lay(i,k) + omega2))))
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) * &
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
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
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)
633 if (
present(kd_lay))
call hchksum(kd_lay,
"Kd_lay", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
635 if (cs%useKappaShear)
call hchksum(visc%Kd_shear,
"Turbulent Kd", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
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)
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, &
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, &
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)
661 if (
present(kd_lay) .and. (cs%id_Kd_layer > 0))
call post_data(cs%id_Kd_layer, kd_lay, cs%diag)
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)
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)
674 if (cs%id_Kd_user > 0)
call post_data(cs%id_Kd_user, dd%Kd_user, cs%diag)
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)
685 if (cs%id_Kd_BBL > 0)
call post_data(cs%id_Kd_BBL, dd%Kd_BBL, cs%diag)
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)
697 if (showcalltree)
call calltree_leave(
"set_diffusivity()")
699 end subroutine set_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)
707 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
711 real,
dimension(SZI_(G),SZK_(G)+1),
intent(in) :: dRho_int
713 real,
dimension(SZI_(G),SZK_(G)),
intent(in) :: N2_lay
715 integer,
intent(in) :: j
716 real,
intent(in) :: dt
718 real,
dimension(SZI_(G),SZK_(G)),
intent(out) :: TKE_to_Kd
723 real,
dimension(SZI_(G),SZK_(G)),
intent(out) :: maxTKE
725 integer,
dimension(SZI_(G)),
intent(out) :: kb
728 real,
dimension(SZI_(G),SZK_(G)) :: &
740 real,
dimension(SZI_(G)) :: &
759 logical :: do_i(SZI_(G))
761 integer,
dimension(2) :: EOSdom
762 integer :: i, k, is, ie, nz, i_rem, kmb, kb_min
763 is = g%isc ; ie = g%iec ; nz = g%ke
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
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)
781 tke_to_kd(i,k) = 1.0 / hn2po2
782 else; tke_to_kd(i,k) = 0.;
endif
785 maxtke(i,k) = hn2po2 * cs%Kd_max
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)
797 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_0, rho_0(:,k), tv%eqn_of_state, eosdom)
799 call calculate_density(tv%T(:,j,kmb), tv%S(:,j,kmb), p_ref, rcv_kmb, tv%eqn_of_state, eosdom)
805 do k=kmb+1,nz-1 ;
if (rcv_kmb(i) <= gv%Rlay(k))
exit ;
enddo
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
816 call set_density_ratios(h, tv, kb, g, gv, us, cs, j, ds_dsp1, rho_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)
825 do k=2,nz-1 ;
do i=is,ie
826 dsp1_ds(i,k) = 1.0 / ds_dsp1(i,k)
828 do i=is,ie ; dsp1_ds(i,nz) = 0.0 ;
enddo
830 if (cs%bulkmixedlayer)
then
831 kmb = gv%nk_rho_varies
833 htot(i) = gv%H_to_Z*h(i,j,kmb)
836 mfkb(i) = ds_dsp1(i,kb(i)) * (gv%H_to_Z*(h(i,j,kmb) - gv%Angstrom_H))
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))
844 maxent(i,1) = 0.0 ; htot(i) = gv%H_to_Z*(h(i,j,1) - gv%Angstrom_H)
847 do k=kb_min,nz-1 ;
do i=is,ie
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))
854 maxent(i,k) = ds_dsp1(i,k)*(maxent(i,k-1) + htot(i))
856 htot(i) = htot(i) + gv%H_to_Z*(h(i,j,k) - gv%Angstrom_H)
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)
866 do i=is,ie ;
if (do_i(i))
then
867 if (k<kb(i))
then ; do_i(i) = .false. ; cycle ;
endif
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)
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
880 do k=2,kmb ;
do i=is,ie
882 tke_to_kd(i,k) = 1.0 / ((n2_lay(i,k) + omega2) * &
883 (gv%H_to_Z*(h(i,j,k) + h_neglect)))
885 do k=kmb+1,kb_min-1 ;
do i=is,ie
888 maxtke(i,k) = 0.0 ; tke_to_kd(i,k) = 0.0
890 do k=kb_min,nz-1 ;
do i=is,ie
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))
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))
913 end subroutine find_tke_to_kd
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)
921 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
925 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
928 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
931 type(
forcing),
intent(in) :: fluxes
932 integer,
intent(in) :: j
934 real,
dimension(SZI_(G),SZK_(G)+1), &
935 intent(out) :: dRho_int
937 real,
dimension(SZI_(G),SZK_(G)+1), &
938 intent(out) :: N2_int
939 real,
dimension(SZI_(G),SZK_(G)), &
940 intent(out) :: N2_lay
941 real,
dimension(SZI_(G)),
intent(out) :: N2_bot
943 real,
dimension(SZI_(G),SZK_(G)+1) :: &
948 real,
dimension(SZI_(G)) :: &
963 logical :: do_i(SZI_(G)), do_any
964 integer,
dimension(2) :: EOSdom
965 integer :: i, k, is, ie, nz
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
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
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
980 do i=is,ie ; pres(i) = 0.0 ;
enddo
982 eosdom(:) = eos_domain(g%HI)
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))
990 tv%eqn_of_state, eosdom)
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)
999 do k=2,nz ;
do i=is,ie
1000 drho_int(i,k) = gv%Rlay(k) - gv%Rlay(k-1)
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))
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))
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)
1021 if (cs%use_tidal_mixing)
call tidal_mixing_h_amp(h_amp, g, j, cs%tidal_mixing_CSp)
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
1029 hb(i) = hb(i) + dz_int
1030 drho_bot(i) = drho_bot(i) + drho_int(i,k)
1032 if (z_from_bot(i) > h_amp(i))
then
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)
1043 if (.not.do_any)
exit
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)
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
1060 n2_int(i,k) = n2_bot(i)
1061 if (k>2) n2_lay(i,k-1) = n2_bot(i)
1063 if (z_from_bot(i) > h_amp(i))
then
1064 if (k>2) n2_int(i,k-1) = n2_bot(i)
1070 if (.not.do_any)
exit
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)
1079 end subroutine find_n2
1088 subroutine double_diffusion(tv, h, T_f, S_f, j, G, GV, US, CS, Kd_T_dd, Kd_S_dd)
1094 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1096 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1099 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1102 integer,
intent(in) :: j
1104 real,
dimension(SZI_(G),SZK_(G)+1), &
1105 intent(out) :: Kd_T_dd
1107 real,
dimension(SZI_(G),SZK_(G)+1), &
1108 intent(out) :: Kd_S_dd
1111 real,
dimension(SZI_(G)) :: &
1126 real,
parameter :: Rrho0 = 1.9
1128 integer,
dimension(2) :: EOSdom
1129 integer :: i, k, is, ie, nz
1130 is = g%isc ; ie = g%iec ; nz = g%ke
1132 if (
associated(tv%eqn_of_state))
then
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
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)
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))
1146 tv%eqn_of_state, eosdom)
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))
1152 if ((alpha_dt > beta_ds) .and. (beta_ds > 0.0))
then
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
1159 rrho = alpha_dt / beta_ds
1160 kd_dd = cs%Kv_molecular * 0.909 * exp(4.6 * exp(-0.54 * (1/rrho - 1)))
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
1166 kd_t_dd(i,k) = 0.0 ; kd_s_dd(i,k) = 0.0
1172 end subroutine double_diffusion
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)
1180 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1182 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1184 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1188 type(
forcing),
intent(in) :: fluxes
1191 integer,
intent(in) :: j
1192 real,
dimension(SZI_(G),SZK_(G)),
intent(in) :: TKE_to_Kd
1197 real,
dimension(SZI_(G),SZK_(G)),
intent(in) :: maxTKE
1199 integer,
dimension(SZI_(G)),
intent(in) :: kb
1202 real,
dimension(SZI_(G),SZK_(G)),
intent(inout) :: Kd_lay
1204 real,
dimension(SZI_(G),SZK_(G)+1), &
1205 optional,
intent(inout) :: Kd_int
1207 real,
dimension(:,:,:),
pointer :: Kd_BBL
1211 real,
dimension(SZK_(G)+1) :: &
1213 real,
dimension(SZI_(G)) :: &
1224 real :: TKE_to_layer
1234 logical :: Rayleigh_drag
1238 logical :: domore, do_i(SZI_(G))
1239 logical :: do_diag_Kd_BBL
1241 integer :: i, k, is, ie, nz, i_rem, kb_min
1242 is = g%isc ; ie = g%iec ; nz = g%ke
1244 do_diag_kd_bbl =
associated(kd_bbl)
1246 if (.not.(cs%bottomdraglaw .and. (cs%BBL_effic>0.0)))
return
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.
1252 i_rho0 = 1.0 / (gv%Rho0)
1253 r0_g = gv%Rho0 / (us%L_to_Z**2 * gv%g_Earth)
1255 do k=2,nz ; rint(k) = 0.5*(gv%Rlay(k-1)+gv%Rlay(k)) ;
enddo
1257 kb_min = max(gv%nk_rho_varies+1,2)
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
1273 i2decay(i) = 0.5*cs%IMax_decay
1275 tke(i) = ((cs%BBL_effic * cdrag_sqrt) * exp(-i2decay(i)*(gv%H_to_Z*h(i,j,nz))) ) * &
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))))
1287 gh_sum_top(i) = r0_g * 400.0 * ustar_h**2
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)
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
1302 rho_top(i) = (rho_htot(i) - gh_sum_top(i)) / htot(i)
1304 elseif (k <= kb(i))
then ; do_i(i) = .false.
1305 else ; domore = .true. ;
endif
1307 if (.not.domore)
exit
1310 do i=is,ie ; do_i(i) = (g%mask2dT(i,j) > 0.5) ;
enddo
1313 do i=is,ie ;
if (do_i(i))
then
1314 if (k<kb(i))
then ; do_i(i) = .false. ; cycle ;
endif
1318 tke(i) = tke(i) * exp(-i2decay(i) * (gv%H_to_Z*(h(i,j,k) + h(i,j,k+1))))
1321 if (maxtke(i,k) <= 0.0) cycle
1325 if (tke(i) > 0.0)
then
1326 if (rint(k) <= rho_top(i))
then
1327 tke_to_layer = tke(i)
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)
1333 else ; tke_to_layer = 0.0 ;
endif
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))
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
1347 tke(i) = tke(i) - tke_to_layer
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
1355 kd_lay(i,k) = (tke_to_layer + tke_ray) * tke_to_kd(i,k)
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
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
1367 if (kd_lay(i,k) >= maxtke(i,k) * tke_to_kd(i,k))
then
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
1375 tke_here = tke_to_layer + tke_ray
1376 tke(i) = tke(i) - tke_to_layer
1378 if (tke(i) < 0.0) tke(i) = 0.0
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
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
1400 if ((tke(i)<= 0.0) .and. (tke_ray == 0.0))
then
1401 do_i(i) = .false. ; i_rem = i_rem - 1
1405 if (i_rem == 0)
exit
1408 end subroutine add_drag_diffusivity
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)
1418 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1420 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1422 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1426 type(
forcing),
intent(in) :: fluxes
1429 integer,
intent(in) :: j
1430 real,
dimension(SZI_(G),SZK_(G)+1), &
1431 intent(in) :: N2_int
1433 real,
dimension(:,:,:),
pointer :: Kd_BBL
1434 real,
dimension(SZI_(G),SZK_(G)), &
1435 optional,
intent(inout) :: Kd_lay
1436 real,
dimension(SZI_(G),SZK_(G)+1), &
1437 optional,
intent(inout) :: Kd_int
1441 real :: TKE_remaining
1442 real :: TKE_consumed
1451 real :: total_thickness
1458 logical :: Rayleigh_drag
1460 integer :: i, k, km1
1461 real,
parameter :: von_karm = 0.41
1462 logical :: do_diag_Kd_BBL
1464 if (.not.(cs%bottomdraglaw .and. (cs%BBL_effic>0.0)))
return
1465 do_diag_kd_bbl =
associated(kd_bbl)
1468 if (cs%LOTW_BBL_use_omega) n2_min = cs%omega**2
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)
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))))
1483 ustar = visc%ustar_BBL(i,j)
1488 if (
associated(fluxes%ustar_tidal)) ustar = ustar + fluxes%ustar_tidal(i,j)
1492 idecay = cs%IMax_decay
1493 if ((ustar > 0.0) .and. (absf > cs%IMax_decay * ustar)) idecay = absf / ustar
1498 tke_column = cdrag_sqrt * visc%TKE_BBL(i,j)
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
1505 tke_remaining = tke_column
1506 total_thickness = ( sum(h(i,j,:)) + gv%H_subroundoff )* gv%H_to_Z
1507 ustar_d = ustar * total_thickness
1514 dh = gv%H_to_Z * h(i,j,k)
1516 dhm1 = gv%H_to_Z * h(i,j,km1)
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))
1528 tke_remaining = exp(-idecay*dh) * tke_remaining
1530 z_bot = z_bot + h(i,j,k)*gv%H_to_Z
1531 d_minus_z = max(total_thickness - z_bot, 0.)
1535 if ( ustar_d + absf * ( z_bot * d_minus_z ) == 0.)
then
1538 kd_wall = ((von_karm * ustar2) * (z_bot * d_minus_z)) &
1539 / (ustar_d + absf * (z_bot * d_minus_z))
1544 tke_kd_wall = kd_wall * 0.5 * (dh + dhm1) * max(n2_int(i,k), n2_min)
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
1552 if (tke_remaining > 0.)
then
1561 tke_remaining = tke_remaining - tke_consumed
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)
1567 if (do_diag_kd_bbl) kd_bbl(i,j,k) = kd_wall
1571 end subroutine add_lotw_bbl_diffusivity
1574 subroutine add_mlrad_diffusivity(h, fluxes, j, G, GV, US, CS, TKE_to_Kd, Kd_lay, Kd_int)
1578 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1580 type(
forcing),
intent(in) :: fluxes
1581 integer,
intent(in) :: j
1583 real,
dimension(SZI_(G),SZK_(G)),
intent(in) :: TKE_to_Kd
1588 real,
dimension(SZI_(G),SZK_(G)), &
1589 optional,
intent(inout) :: Kd_lay
1590 real,
dimension(SZI_(G),SZK_(G)+1), &
1591 optional,
intent(inout) :: Kd_int
1596 real,
dimension(SZI_(G)) :: h_ml
1597 real,
dimension(SZI_(G)) :: TKE_ml_flux
1598 real,
dimension(SZI_(G)) :: I_decay
1599 real,
dimension(SZI_(G)) :: Kd_mlr_ml
1609 real :: I_decay_len2_TKE
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
1617 omega2 = cs%omega**2
1620 h_neglect = gv%H_subroundoff*gv%H_to_Z
1622 if (.not.cs%ML_radiation)
return
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
1627 do i=is,ie ;
if (do_i(i))
then
1628 if (cs%ML_omega_frac >= 1.0)
then
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
1637 ustar_sq = max(fluxes%ustar(i,j), cs%ustar_min)**2
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)
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))
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)
1651 z1 = (gv%H_to_Z*h(i,j,kml+1)) * i_decay(i)
1653 kd_mlr = tke_ml_flux(i) * tke_to_kd(i,kml+1) * (1.0 - exp(-z1))
1655 kd_mlr = tke_ml_flux(i) * tke_to_kd(i,kml+1) * (z1 * (1.0 - z1 * (0.5 - c1_6 * z1)))
1657 kd_mlr_ml(i) = min(kd_mlr, cs%ML_rad_kd_max)
1658 tke_ml_flux(i) = tke_ml_flux(i) * exp(-z1)
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
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
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
1684 kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * &
1685 us%m_to_Z * ((1.0 - exp(-z1)) / dzl)
1687 kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * &
1688 us%m_to_Z * (i_decay(i) * (1.0 - z1 * (0.5 - c1_6*z1)))
1692 kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * (1.0 - exp(-z1))
1694 kd_mlr = (tke_ml_flux(i) * tke_to_kd(i,k)) * (z1 * (1.0 - z1 * (0.5 - c1_6*z1)))
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
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
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
1709 else ; do_any = .true. ;
endif
1711 if (.not.do_any)
exit
1714 end subroutine add_mlrad_diffusivity
1718 subroutine set_bbl_tke(u, v, h, fluxes, visc, G, GV, US, CS, OBC)
1722 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1724 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1726 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1728 type(
forcing),
intent(in) :: fluxes
1737 real,
dimension(SZI_(G)) :: &
1741 real,
dimension(SZIB_(G)) :: &
1746 real :: vhtot(szi_(g))
1748 real,
dimension(SZI_(G),SZJB_(G)) :: &
1755 logical :: domore, do_i(szi_(g))
1756 integer :: i, j, k, is, ie, js, je, nz
1758 logical :: local_open_u_bc, local_open_v_bc
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
1768 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1770 if (.not.
associated(cs))
call mom_error(fatal,
"set_BBL_TKE: "//&
1771 "Module must be initialized before it is used.")
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
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
1783 cdrag_sqrt = sqrt(cs%cdrag)
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))
1795 do_i(i) = .false. ; vstar(i,j) = 0.0 ; htot(i) = 0.0
1799 do i=is,ie ;
if (do_i(i))
then
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
1811 if (obc%segment(l_seg)%direction == obc_direction_n)
then
1812 hvel = gv%H_to_Z*h(i,j,k)
1814 hvel = gv%H_to_Z*h(i,j+1,k)
1817 hvel = 0.5*gv%H_to_Z*(h(i,j,k) + h(i,j+1,k))
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)
1825 vhtot(i) = vhtot(i) + hvel*v(i,j,k)
1826 htot(i) = htot(i) + hvel
1830 if (.not.domore)
exit
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))
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))
1844 do_i(i) = .false. ; ustar(i) = 0.0 ; htot(i) = 0.0
1846 do k=nz,1,-1 ; domore = .false.
1847 do i=is-1,ie ;
if (do_i(i))
then
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
1859 if (obc%segment(l_seg)%direction == obc_direction_e)
then
1860 hvel = gv%H_to_Z*h(i,j,k)
1862 hvel = gv%H_to_Z*h(i+1,j,k)
1865 hvel = 0.5*gv%H_to_Z*(h(i,j,k) + h(i+1,j,k))
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)
1873 uhtot(i) = uhtot(i) + hvel*u(i,j,k)
1874 htot(i) = htot(i) + hvel
1878 if (.not.domore)
exit
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))
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))
1901 end subroutine set_bbl_tke
1903 subroutine set_density_ratios(h, tv, kb, G, GV, US, CS, j, ds_dsp1, rho_0)
1906 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1911 integer,
dimension(SZI_(G)),
intent(in) :: kb
1916 integer,
intent(in) :: j
1917 real,
dimension(SZI_(G),SZK_(G)),
intent(out) :: ds_dsp1
1921 real,
dimension(SZI_(G),SZK_(G)), &
1922 optional,
intent(in) :: rho_0
1928 real :: a(SZK_(G)), a_0(SZK_(G))
1929 real :: p_ref(SZI_(G))
1930 real :: Rcv(SZI_(G),SZK_(G))
1933 integer,
dimension(2) :: EOSdom
1934 integer :: i, k, k3, is, ie, nz, kmb
1935 is = g%isc ; ie = g%iec ; nz = g%ke
1938 if (gv%g_prime(k+1) /= 0.0)
then
1940 ds_dsp1(i,k) = gv%g_prime(k) / gv%g_prime(k+1)
1949 if (cs%bulkmixedlayer)
then
1950 g_r0 = gv%g_Earth / (gv%Rho0)
1951 kmb = gv%nk_rho_varies
1953 do i=is,ie ; p_ref(i) = tv%P_Ref ;
enddo
1954 eosdom(:) = eos_domain(g%HI)
1956 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_ref, rcv(:,k), tv%eqn_of_state, eosdom)
1959 if (kb(i) <= nz-1)
then
1964 i_drho = g_r0 / gv%g_prime(k+1)
1967 a(k3+1) = (gv%Rlay(k) - rcv(i,k3)) * i_drho
1969 if ((
present(rho_0)) .and. (a(kmb+1) < 2.0*eps*ds_dsp1(i,k)))
then
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
1979 a_0(k3) = a_0(kmb+1) + (rho_0(i,kmb)-rho_0(i,k3-1)) * i_drho
1981 if (a(kmb+1) <= eps*ds_dsp1(i,k))
then
1982 do k3=2,kmb+1 ; a(k3) = a_0(k3) ;
enddo
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
1992 ds_dsp1(i,k) = max(a(kmb+1),1e-5)
2001 ds_dsp1(i,k3) = max(a(k3),ds_dsp1(i,k))
2007 end subroutine set_density_ratios
2009 subroutine set_diffusivity_init(Time, G, GV, US, param_file, diag, CS, int_tide_CSp, halo_TS, &
2011 type(time_type),
intent(in) :: time
2017 type(
diag_ctrl),
target,
intent(inout) :: diag
2022 integer,
optional,
intent(out) :: halo_ts
2024 logical,
optional,
intent(out) :: double_diffuse
2028 real :: decay_length
2029 logical :: ml_use_omega
2030 logical :: default_2018_answers
2032 # include "version_variable.h"
2033 character(len=40) :: mdl =
"MOM_set_diffusivity"
2034 real :: omega_frac_dflt
2035 logical :: bryan_lewis_diffusivity
2037 integer :: i, j, is, ie, js, je
2038 integer :: isd, ied, jsd, jed
2040 if (
associated(cs))
then
2041 call mom_error(warning,
"diabatic_entrain_init called with an associated "// &
2042 "control structure.")
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
2051 if (
associated(int_tide_csp)) cs%int_tide_CSp => int_tide_csp
2054 cs%BBL_mixing_as_max = .true.
2055 cs%cdrag = 0.003 ; cs%BBL_effic = 0.0
2056 cs%bulkmixedlayer = (gv%nkml > 0)
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)
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.", &
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)
2080 cs%use_tidal_mixing = tidal_mixing_init(time, g, gv, us, param_file, diag, &
2081 cs%tidal_mixing_CSp)
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
2091 cs%ustar_min = 2e-4 * cs%omega * (gv%Angstrom_Z + gv%H_subroundoff*gv%H_to_Z)
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
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)
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)
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.", &
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.)
2174 cs%use_LOTW_BBL_diffusivity = .false.
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.)
2185 call bkgnd_mixing_init(time, g, gv, us, param_file, cs%diag, cs%bkgnd_mixing_csp)
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.)
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)
2217 call get_param(param_file, mdl,
"DEBUG", cs%debug, &
2218 "If true, write out verbose debugging data.", &
2219 default=.false., debuggingparam=.true.)
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.)
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)
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.)
2245 if (cs%FluxRi_max > 0.0) &
2246 cs%dissip_N2 = cs%dissip_Kd_min * gv%Rho0 / cs%FluxRi_max
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)
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)
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')
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)
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).", &
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)
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)
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.")
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)
2306 if (cs%useKappaShear) &
2307 id_clock_kappashear = cpu_clock_id(
'(Ocean kappa_shear)', grain=clock_module)
2310 cs%use_CVMix_shear = cvmix_shear_init(time, g, gv, us, param_file, cs%diag, cs%CVMix_shear_csp)
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)
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.')
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)
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')
2334 if (
present(halo_ts))
then
2336 if (cs%Vertex_Shear) halo_ts = 1
2339 if (
present(double_diffuse))
then
2340 double_diffuse = (cs%double_diffusion .or. cs%use_CVMix_ddiff)
2343 end subroutine set_diffusivity_init
2346 subroutine set_diffusivity_end(CS)
2349 if (.not.
associated(cs))
return
2351 call bkgnd_mixing_end(cs%bkgnd_mixing_csp)
2353 if (cs%use_tidal_mixing)
call tidal_mixing_end(cs%tidal_mixing_CSp)
2355 if (cs%user_change_diff)
call user_change_diff_end(cs%user_change_diff_CSp)
2357 if (cs%use_CVMix_shear)
call cvmix_shear_end(cs%CVMix_shear_csp)
2359 if (cs%use_CVMix_ddiff)
call cvmix_ddiff_end(cs%CVMix_ddiff_csp)
2361 if (
associated(cs))
deallocate(cs)
2363 end subroutine set_diffusivity_end