8 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end, clock_routine
23 use mom_restart,
only : register_restart_field_as_obsolete
32 implicit none ;
private
34 #include <MOM_memory.h>
36 public set_viscous_bbl, set_viscous_ml, set_visc_init, set_visc_end
37 public set_visc_register_restarts
60 real :: htbl_shelf_min
63 logical :: bottomdraglaw
69 logical :: bbl_use_eos
71 logical :: linear_drag
73 logical :: channel_drag
76 logical :: correct_bbl_bounds
79 logical :: dynamic_viscous_ml
92 logical :: answers_2018
96 logical :: bbl_use_tidal_bg
98 character(len=200) :: inputdir
103 real,
allocatable,
dimension(:,:) :: tideamp
105 real,
allocatable,
dimension(:,:) :: bbl_u
106 real,
allocatable,
dimension(:,:) :: bbl_v
108 integer :: id_bbl_thick_u = -1, id_kv_bbl_u = -1, id_bbl_u = -1
109 integer :: id_bbl_thick_v = -1, id_kv_bbl_v = -1, id_bbl_v = -1
110 integer :: id_ray_u = -1, id_ray_v = -1
111 integer :: id_nkml_visc_u = -1, id_nkml_visc_v = -1
192 subroutine set_viscous_bbl(u, v, h, tv, visc, G, GV, US, CS, symmetrize)
196 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
198 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
200 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
209 logical,
optional,
intent(in) :: symmetrize
214 real,
dimension(SZIB_(G)) :: &
230 real,
dimension(SZIB_(G),SZJ_(G)) :: &
234 real,
dimension(SZI_(G),SZJB_(G)) :: &
238 real,
dimension(SZIB_(G),SZK_(G)) :: &
283 real :: v_at_u, u_at_v
287 real,
dimension(SZI_(G),SZJ_(G),max(GV%nk_rho_varies,1)) :: &
289 real :: p_ref(szi_(g))
300 real :: apb_4a, ax2_3apb
301 real :: a2x48_apb3, iapb, ibma_2
341 real :: bbl_visc_frac
343 real,
parameter :: c1_3 = 1.0/3.0, c1_6 = 1.0/6.0, c1_12 = 1.0/12.0
346 real :: tmp_val_m1_to_p1
349 logical :: use_l0, do_one_l_iter
350 logical :: use_bbl_eos, do_i(szib_(g))
351 integer,
dimension(2) :: eosdom
352 integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz, m, n, k2, nkmb, nkml
353 integer :: itt, maxitt=20
356 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
357 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
358 nkmb = gv%nk_rho_varies ; nkml = gv%nkml
359 h_neglect = gv%H_subroundoff
360 rho0x400_g = 400.0*(gv%Rho0 / (us%L_to_Z**2 * gv%g_Earth)) * gv%Z_to_H
361 vol_quit = 0.9*gv%Angstrom_H + h_neglect
362 c2pi_3 = 8.0*atan(1.0)/3.0
364 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_set_viscosity(BBL): "//&
365 "Module must be initialized before it is used.")
366 if (.not.cs%bottomdraglaw)
return
368 if (
present(symmetrize))
then ;
if (symmetrize)
then
369 jsq = js-1 ; isq = is-1
373 call uvchksum(
"Start set_viscous_BBL [uv]", u, v, g%HI, haloshift=1, scale=us%L_T_to_m_s)
374 call hchksum(h,
"Start set_viscous_BBL h", g%HI, haloshift=1, scale=gv%H_to_m)
375 if (
associated(tv%T))
call hchksum(tv%T,
"Start set_viscous_BBL T", g%HI, haloshift=1)
376 if (
associated(tv%S))
call hchksum(tv%S,
"Start set_viscous_BBL S", g%HI, haloshift=1)
379 use_bbl_eos =
associated(tv%eqn_of_state) .and. cs%BBL_use_EOS
382 u_bg_sq = cs%drag_bg_vel * cs%drag_bg_vel
383 cdrag_sqrt = sqrt(cs%cdrag)
384 cdrag_sqrt_z = us%L_to_Z * sqrt(cs%cdrag)
390 if ((nkml>0) .and. .not.use_bbl_eos)
then
391 eosdom(1) = isq - (g%isd-1) ; eosdom(2) = g%iec+1 - (g%isd-1)
392 do i=isq,ieq+1 ; p_ref(i) = tv%P_Ref ;
enddo
394 do k=1,nkmb ;
do j=jsq,jeq+1
395 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_ref, rml(:,j,k), tv%eqn_of_state, &
401 do j=js-1,je ;
do i=is-1,ie+1
402 d_v(i,j) = 0.5*(g%bathyT(i,j) + g%bathyT(i,j+1))
403 mask_v(i,j) = g%mask2dCv(i,j)
406 do j=js-1,je+1 ;
do i=is-1,ie
407 d_u(i,j) = 0.5*(g%bathyT(i,j) + g%bathyT(i+1,j))
408 mask_u(i,j) = g%mask2dCu(i,j)
411 if (
associated(obc))
then ;
do n=1,obc%number_of_segments
412 if (.not. obc%segment(n)%on_pe) cycle
414 i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
415 if (obc%segment(n)%is_N_or_S .and. (j >= js-1) .and. (j <= je))
then
416 do i = max(is-1,obc%segment(n)%HI%isd), min(ie+1,obc%segment(n)%HI%ied)
417 if (obc%segment(n)%direction == obc_direction_n) d_v(i,j) = g%bathyT(i,j)
418 if (obc%segment(n)%direction == obc_direction_s) d_v(i,j) = g%bathyT(i,j+1)
420 elseif (obc%segment(n)%is_E_or_W .and. (i >= is-1) .and. (i <= ie))
then
421 do j = max(js-1,obc%segment(n)%HI%jsd), min(je+1,obc%segment(n)%HI%jed)
422 if (obc%segment(n)%direction == obc_direction_e) d_u(i,j) = g%bathyT(i,j)
423 if (obc%segment(n)%direction == obc_direction_w) d_u(i,j) = g%bathyT(i+1,j)
427 if (
associated(obc))
then ;
do n=1,obc%number_of_segments
430 if (.not. obc%segment(n)%on_pe) cycle
432 i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
433 if (obc%segment(n)%is_N_or_S .and. (j >= js-1) .and. (j <= je))
then
434 do i = max(is-1,obc%segment(n)%HI%IsdB), min(ie,obc%segment(n)%HI%IedB)
435 if (obc%segment(n)%direction == obc_direction_n)
then
436 d_u(i,j+1) = d_u(i,j) ; mask_u(i,j+1) = 0.0
437 elseif (obc%segment(n)%direction == obc_direction_s)
then
438 d_u(i,j) = d_u(i,j+1) ; mask_u(i,j) = 0.0
441 elseif (obc%segment(n)%is_E_or_W .and. (i >= is-1) .and. (i <= ie))
then
442 do j = max(js-1,obc%segment(n)%HI%JsdB), min(je,obc%segment(n)%HI%JedB)
443 if (obc%segment(n)%direction == obc_direction_e)
then
444 d_v(i+1,j) = d_v(i,j) ; mask_v(i+1,j) = 0.0
445 elseif (obc%segment(n)%direction == obc_direction_w)
then
446 d_v(i,j) = d_v(i+1,j) ; mask_v(i,j) = 0.0
452 if (.not.use_bbl_eos) rml_vel(:,:) = 0.0
459 do j=jsq,jeq ;
do m=1,2
467 if (g%mask2dCu(i,j) > 0) do_i(i) = .true.
471 is = g%isc ; ie = g%iec
474 if (g%mask2dCv(i,j) > 0) do_i(i) = .true.
481 do k=1,nz ;
do i=is,ie
483 if (u(i,j,k) *(h(i+1,j,k) - h(i,j,k)) >= 0)
then
485 h_at_vel(i,k) = 2.0*h(i,j,k)*h(i+1,j,k) / &
486 (h(i,j,k) + h(i+1,j,k) + h_neglect)
489 h_at_vel(i,k) = 0.5 * (h(i,j,k) + h(i+1,j,k))
492 h_vel(i,k) = 0.5 * (h(i,j,k) + h(i+1,j,k))
494 if (use_bbl_eos)
then ;
do k=1,nz ;
do i=is,ie
496 t_vel(i,k) = 0.5 * (tv%T(i,j,k) + tv%T(i+1,j,k))
497 s_vel(i,k) = 0.5 * (tv%S(i,j,k) + tv%S(i+1,j,k))
498 enddo ;
enddo ;
else ;
do k=1,nkmb ;
do i=is,ie
499 rml_vel(i,k) = 0.5 * (rml(i,j,k) + rml(i+1,j,k))
500 enddo ;
enddo ;
endif
502 do k=1,nz ;
do i=is,ie
504 if (v(i,j,k) * (h(i,j+1,k) - h(i,j,k)) >= 0)
then
506 h_at_vel(i,k) = 2.0*h(i,j,k)*h(i,j+1,k) / &
507 (h(i,j,k) + h(i,j+1,k) + h_neglect)
510 h_at_vel(i,k) = 0.5 * (h(i,j,k) + h(i,j+1,k))
513 h_vel(i,k) = 0.5 * (h(i,j,k) + h(i,j+1,k))
515 if (use_bbl_eos)
then ;
do k=1,nz ;
do i=is,ie
516 t_vel(i,k) = 0.5 * (tv%T(i,j,k) + tv%T(i,j+1,k))
517 s_vel(i,k) = 0.5 * (tv%S(i,j,k) + tv%S(i,j+1,k))
518 enddo ;
enddo ;
else ;
do k=1,nkmb ;
do i=is,ie
519 rml_vel(i,k) = 0.5 * (rml(i,j,k) + rml(i,j+1,k))
520 enddo ;
enddo ;
endif
523 if (
associated(obc))
then ;
if (obc%number_of_segments > 0)
then
526 do i=is,ie ;
if (do_i(i) .and. (obc%segnum_u(i,j) /= obc_none))
then
527 if (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_e)
then
529 h_at_vel(i,k) = h(i,j,k) ; h_vel(i,k) = h(i,j,k)
531 if (use_bbl_eos)
then
533 t_vel(i,k) = tv%T(i,j,k) ; s_vel(i,k) = tv%S(i,j,k)
537 rml_vel(i,k) = rml(i,j,k)
540 elseif (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_w)
then
542 h_at_vel(i,k) = h(i+1,j,k) ; h_vel(i,k) = h(i+1,j,k)
544 if (use_bbl_eos)
then
546 t_vel(i,k) = tv%T(i+1,j,k) ; s_vel(i,k) = tv%S(i+1,j,k)
550 rml_vel(i,k) = rml(i+1,j,k)
556 do i=is,ie ;
if (do_i(i) .and. (obc%segnum_v(i,j) /= obc_none))
then
557 if (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_n)
then
559 h_at_vel(i,k) = h(i,j,k) ; h_vel(i,k) = h(i,j,k)
561 if (use_bbl_eos)
then
563 t_vel(i,k) = tv%T(i,j,k) ; s_vel(i,k) = tv%S(i,j,k)
567 rml_vel(i,k) = rml(i,j,k)
570 elseif (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_s)
then
572 h_at_vel(i,k) = h(i,j+1,k) ; h_vel(i,k) = h(i,j+1,k)
574 if (use_bbl_eos)
then
576 t_vel(i,k) = tv%T(i,j+1,k) ; s_vel(i,k) = tv%S(i,j+1,k)
580 rml_vel(i,k) = rml(i,j+1,k)
588 if (use_bbl_eos .or. .not.cs%linear_drag)
then
592 do i=is,ie ;
if (do_i(i))
then
593 htot_vel = 0.0 ; hwtot = 0.0 ; hutot = 0.0
594 thtot = 0.0 ; shtot = 0.0
597 if (htot_vel>=cs%Hbbl)
exit
599 hweight = min(cs%Hbbl - htot_vel, h_at_vel(i,k))
600 if (hweight < 1.5*gv%Angstrom_H + h_neglect) cycle
602 htot_vel = htot_vel + h_at_vel(i,k)
603 hwtot = hwtot + hweight
605 if ((.not.cs%linear_drag) .and. (hweight >= 0.0))
then ;
if (m==1)
then
606 v_at_u = set_v_at_u(v, h, g, i, j, k, mask_v, obc)
607 if (cs%BBL_use_tidal_bg)
then
608 u_bg_sq = 0.5*( g%mask2dT(i,j)*(cs%tideamp(i,j)*cs%tideamp(i,j))+ &
609 g%mask2dT(i+1,j)*(cs%tideamp(i+1,j)*cs%tideamp(i+1,j)) )
611 hutot = hutot + hweight * sqrt(u(i,j,k)*u(i,j,k) + &
612 v_at_u*v_at_u + u_bg_sq)
614 u_at_v = set_u_at_v(u, h, g, i, j, k, mask_u, obc)
615 if (cs%BBL_use_tidal_bg)
then
616 u_bg_sq = 0.5*( g%mask2dT(i,j)*(cs%tideamp(i,j)*cs%tideamp(i,j))+ &
617 g%mask2dT(i,j+1)*(cs%tideamp(i,j+1)*cs%tideamp(i,j+1)) )
619 hutot = hutot + hweight * sqrt(v(i,j,k)*v(i,j,k) + &
620 u_at_v*u_at_v + u_bg_sq)
623 if (use_bbl_eos .and. (hweight >= 0.0))
then
624 thtot = thtot + hweight * t_vel(i,k)
625 shtot = shtot + hweight * s_vel(i,k)
630 if (.not.cs%linear_drag .and. (hwtot > 0.0))
then
631 ustar(i) = cdrag_sqrt_z*hutot/hwtot
633 ustar(i) = cdrag_sqrt_z*cs%drag_bg_vel
636 if (use_bbl_eos)
then ;
if (hwtot > 0.0)
then
637 t_eos(i) = thtot/hwtot ; s_eos(i) = shtot/hwtot
639 t_eos(i) = 0.0 ; s_eos(i) = 0.0
643 if (cs%id_bbl_u>0 .and. m==1)
then
644 if (hwtot > 0.0) cs%bbl_u(i,j) = hutot/hwtot
645 elseif (cs%id_bbl_v>0 .and. m==2)
then
646 if (hwtot > 0.0) cs%bbl_v(i,j) = hutot/hwtot
651 do i=is,ie ; ustar(i) = cdrag_sqrt_z*cs%drag_bg_vel ;
enddo
654 if (use_bbl_eos)
then
655 if (
associated(tv%p_surf))
then
656 if (m==1)
then ;
do i=is,ie ; press(i) = 0.5*(tv%p_surf(i,j) + tv%p_surf(i+1,j)) ;
enddo
657 else ;
do i=is,ie ; press(i) = 0.5*(tv%p_surf(i,j) + tv%p_surf(i,j+1)) ;
enddo ;
endif
659 do i=is,ie ; press(i) = 0.0 ;
enddo
661 do i=is,ie ;
if (.not.do_i(i))
then ; t_eos(i) = 0.0 ; s_eos(i) = 0.0 ;
endif ;
enddo
662 do k=1,nz ;
do i=is,ie
663 press(i) = press(i) + (gv%H_to_RZ*gv%g_Earth) * h_vel(i,k)
666 (/is-g%IsdB+1,ie-g%IsdB+1/) )
676 do i=is,ie ;
if (do_i(i))
then
678 ustarsq = rho0x400_g * ustar(i)**2
690 if (use_bbl_eos)
then
691 thtot = 0.0 ; shtot = 0.0 ; oldfn = 0.0
693 if (h_at_vel(i,k) <= 0.0) cycle
696 oldfn = dr_dt(i)*(thtot - t_vel(i,k)*htot) + &
697 dr_ds(i)*(shtot - s_vel(i,k)*htot)
698 if (oldfn >= ustarsq)
exit
701 dfn = (dr_dt(i)*(t_vel(i,k) - t_vel(i,k-1)) + &
702 dr_ds(i)*(s_vel(i,k) - s_vel(i,k-1))) * &
703 (h_at_vel(i,k) + htot)
705 if ((oldfn + dfn) <= ustarsq)
then
710 dh = h_at_vel(i,k) * sqrt((ustarsq-oldfn) / (dfn))
715 thtot = thtot + t_vel(i,k)*dh ; shtot = shtot + s_vel(i,k)*dh
717 if ((oldfn < ustarsq) .and. h_at_vel(i,1) > 0.0)
then
719 if (dr_dt(i) * (thtot - t_vel(i,1)*htot) + &
720 dr_ds(i) * (shtot - s_vel(i,1)*htot) < ustarsq) &
721 htot = htot + h_at_vel(i,1)
726 oldfn = rhtot - gv%Rlay(k)*htot
727 dfn = (gv%Rlay(k) - gv%Rlay(k-1))*(h_at_vel(i,k)+htot)
729 if (oldfn >= ustarsq)
then
731 elseif ((oldfn + dfn) <= ustarsq)
then
734 dh = h_at_vel(i,k) * sqrt((ustarsq-oldfn) / (dfn))
738 rhtot = rhtot + gv%Rlay(k)*dh
742 oldfn = rhtot - rml_vel(i,k)*htot
743 dfn = (rml_vel(i,k) - rml_vel(i,k-1)) * (h_at_vel(i,k)+htot)
745 if (oldfn >= ustarsq)
then
747 elseif ((oldfn + dfn) <= ustarsq)
then
750 dh = h_at_vel(i,k) * sqrt((ustarsq-oldfn) / (dfn))
754 rhtot = rhtot + rml_vel(i,k)*dh
756 if (rhtot - rml_vel(i,1)*htot < ustarsq) htot = htot + h_at_vel(i,1)
758 if (rhtot - gv%Rlay(1)*htot < ustarsq) htot = htot + h_at_vel(i,1)
763 if (m==1)
then ; c2f = g%CoriolisBu(i,j-1) + g%CoriolisBu(i,j)
764 else ; c2f = g%CoriolisBu(i-1,j) + g%CoriolisBu(i,j) ;
endif
777 if (cs%cdrag * u_bg_sq <= 0.0)
then
780 usth = ustar(i)*gv%Z_to_H ; root = sqrt(0.25*usth**2 + (htot*c2f)**2)
781 if (htot*usth <= (cs%BBL_thick_min+h_neglect) * (0.5*usth + root))
then
782 bbl_thick = cs%BBL_thick_min
787 bbl_thick = (htot * usth) / (0.5*usth + root)
793 bbl_thick = htot / (0.5 + sqrt(0.25 + htot*htot*c2f*c2f/ &
794 ((ustar(i)*ustar(i)) * (gv%Z_to_H**2)) ) )
796 if (bbl_thick < cs%BBL_thick_min) bbl_thick = cs%BBL_thick_min
806 if ((bbl_thick > 0.5*cs%Hbbl) .and. (cs%RiNo_mix)) bbl_thick = 0.5*cs%Hbbl
808 if (cs%Channel_drag)
then
814 tmp = g%mask2dCu(i,j+1) * d_u(i,j+1)
815 dp = 2.0 * d_vel * tmp / (d_vel + tmp)
816 tmp = g%mask2dCu(i,j-1) * d_u(i,j-1)
817 dm = 2.0 * d_vel * tmp / (d_vel + tmp)
820 tmp = g%mask2dCv(i+1,j) * d_v(i+1,j)
821 dp = 2.0 * d_vel * tmp / (d_vel + tmp)
822 tmp = g%mask2dCv(i-1,j) * d_v(i-1,j)
823 dm = 2.0 * d_vel * tmp / (d_vel + tmp)
825 if (dm > dp)
then ; tmp = dp ; dp = dm ; dm = tmp ;
endif
828 dp = gv%Z_to_H*dp ; dm = gv%Z_to_H*dm ; d_vel = gv%Z_to_H*d_vel
830 a_3 = (dp + dm - 2.0*d_vel) ; a = 3.0*a_3 ; a_12 = 0.25*a_3
834 if (abs(a) < 1e-2*(slope + cs%BBL_thick_min)) a = 0.0
842 vol_open = d_vel - dm ; vol_2_reg = vol_open
845 vol_open = 0.25*slope*tmp + c1_12*a
846 vol_2_reg = 0.5*tmp**2 * (a - c1_3*slope)
849 c24_a = 24.0/a ; iapb = 1.0/(a+slope)
850 apb_4a = (slope+a)/(4.0*a) ; a2x48_apb3 = (48.0*(a*a))*(iapb**3)
851 ax2_3apb = 2.0*c1_3*a*iapb
852 elseif (a == 0.0)
then
854 if (slope > 0) iapb = 1.0/slope
856 vol_open = d_vel - dm
857 if (slope >= -a)
then
858 iapb = 1.0e30 ;
if (slope+a /= 0.0) iapb = 1.0/(a+slope)
859 vol_direct = 0.0 ; l_direct = 0.0 ; c24_a = 0.0
861 c24_a = 24.0/a ; iapb = 1.0/(a+slope)
862 l_direct = 1.0 + slope/a
863 vol_direct = -c1_6*a*l_direct**3
865 ibma_2 = 2.0 / (slope - a)
868 l(nz+1) = 0.0 ; vol = 0.0 ; vol_err = 0.0 ; bbl_visc_frac = 0.0
873 vol = vol + h_vel(i,k)
874 h_vel_pos = h_vel(i,k) + h_neglect
876 if (vol >= vol_open)
then ; l(k) = 1.0
878 l(k) = sqrt(2.0*vol*iapb)
882 if (vol < vol_2_reg)
then
885 if (a2x48_apb3*vol < 1e-8)
then
887 l0 = sqrt(2.0*vol*iapb) ; l(k) = l0*(1.0 + ax2_3apb*l0)
889 l(k) = apb_4a * (1.0 - &
890 2.0 * cos(c1_3*acos(a2x48_apb3*vol - 1.0) - c2pi_3))
898 tmp_val_m1_to_p1 = 1.0 - c24_a*(vol_open - vol)
899 tmp_val_m1_to_p1 = max(-1., min(1., tmp_val_m1_to_p1))
900 l(k) = 0.5 - cos(c1_3*acos(tmp_val_m1_to_p1) - c2pi_3)
905 if (vol <= vol_direct)
then
907 l(k) = (-0.25*c24_a*vol)**c1_3
915 if (vol_below + vol_err <= vol_direct)
then
916 l0 = l_direct ; vol_0 = vol_direct
918 l0 = l(k+1) ; vol_0 = vol_below + vol_err
924 dv_dl2 = 0.5*(slope+a) - a*l0 ; dvol = (vol-vol_0)
928 do_one_l_iter = .false.
929 if (cs%answers_2018)
then
930 curv_tol = gv%Angstrom_H*dv_dl2**2 &
931 * (0.25 * dv_dl2 * gv%Angstrom_H - a * l0 * dvol)
932 do_one_l_iter = (a * a * dvol**3) < curv_tol
936 use_l0 = (dvol <= 0.)
938 vol_tol = max(0.5 * gv%Angstrom_H + gv%H_subroundoff, 1e-14 * vol)
939 vol_quit = max(0.9 * gv%Angstrom_H + gv%H_subroundoff, 1e-14 * vol)
941 curv_tol = vol_tol * dv_dl2**2 &
942 * (dv_dl2 * vol_tol - 2.0 * a * l0 * dvol)
943 do_one_l_iter = (a * a * dvol**3) < curv_tol
948 vol_err = 0.5*(l(k)*l(k))*(slope + a_3*(3.0-4.0*l(k))) - vol
949 elseif (do_one_l_iter)
then
952 l(k) = sqrt(l0*l0 + dvol / dv_dl2)
953 vol_err = 0.5*(l(k)*l(k))*(slope + a_3*(3.0-4.0*l(k))) - vol
955 if (dv_dl2*(1.0-l0*l0) < dvol + &
956 dv_dl2 * (vol_open - vol)*ibma_2)
then
957 l_max = sqrt(1.0 - (vol_open - vol)*ibma_2)
959 l_max = sqrt(l0*l0 + dvol / dv_dl2)
961 l_min = sqrt(l0*l0 + dvol / (0.5*(slope+a) - a*l_max))
963 vol_err_min = 0.5*(l_min**2)*(slope + a_3*(3.0-4.0*l_min)) - vol
964 vol_err_max = 0.5*(l_max**2)*(slope + a_3*(3.0-4.0*l_max)) - vol
966 if (abs(vol_err_min) <= vol_quit)
then
967 l(k) = l_min ; vol_err = vol_err_min
969 l(k) = sqrt((l_min**2*vol_err_max - l_max**2*vol_err_min) / &
970 (vol_err_max - vol_err_min))
972 vol_err = 0.5*(l(k)*l(k))*(slope + a_3*(3.0-4.0*l(k))) - vol
973 if (abs(vol_err) <= vol_quit)
exit
976 l(k) = l(k) - vol_err / (l(k)* (slope + a - 2.0*a*l(k)))
987 if (l(k) > l(k+1))
then
988 if (vol_below < bbl_thick)
then
989 bbl_frac = (1.0-vol_below/bbl_thick)**2
990 bbl_visc_frac = bbl_visc_frac + bbl_frac*(l(k) - l(k+1))
995 if (m==1)
then ; cell_width = g%dy_Cu(i,j)
996 else ; cell_width = g%dx_Cv(i,j) ;
endif
997 gam = 1.0 - l(k+1)/l(k)
998 rayleigh = us%L_to_Z * cs%cdrag * (l(k)-l(k+1)) * (1.0-bbl_frac) * &
999 (12.0*cs%c_Smag*h_vel_pos) / (12.0*cs%c_Smag*h_vel_pos + &
1000 us%L_to_Z*gv%Z_to_H * cs%cdrag * gam*(1.0-gam)*(1.0-1.5*gam) * l(k)**2 * cell_width)
1006 if (rayleigh > 0.0)
then
1007 v_at_u = set_v_at_u(v, h, g, i, j, k, mask_v, obc)
1008 visc%Ray_u(i,j,k) = rayleigh*sqrt(u(i,j,k)*u(i,j,k) + &
1009 v_at_u*v_at_u + u_bg_sq)
1010 else ; visc%Ray_u(i,j,k) = 0.0 ;
endif
1012 if (rayleigh > 0.0)
then
1013 u_at_v = set_u_at_v(u, h, g, i, j, k, mask_u, obc)
1014 visc%Ray_v(i,j,k) = rayleigh*sqrt(v(i,j,k)*v(i,j,k) + &
1015 u_at_v*u_at_v + u_bg_sq)
1016 else ; visc%Ray_v(i,j,k) = 0.0 ;
endif
1024 bbl_thick_z = bbl_thick * gv%H_to_Z
1025 if (cs%correct_BBL_bounds .and. &
1026 cdrag_sqrt*ustar(i)*bbl_thick_z*bbl_visc_frac <= cs%Kv_BBL_min)
then
1030 kv_bbl = cs%Kv_BBL_min
1031 if (cdrag_sqrt*ustar(i)*bbl_visc_frac*g%Rad_Earth*us%m_to_Z > kv_bbl)
then
1032 bbl_thick_z = kv_bbl / ( cdrag_sqrt*ustar(i)*bbl_visc_frac )
1034 bbl_thick_z = g%Rad_Earth * us%m_to_Z
1037 kv_bbl = cdrag_sqrt*ustar(i)*bbl_thick_z*bbl_visc_frac
1056 bbl_thick_z = bbl_thick * gv%H_to_Z
1057 if (cs%correct_BBL_bounds .and. &
1058 cdrag_sqrt*ustar(i)*bbl_thick_z <= cs%Kv_BBL_min)
then
1062 kv_bbl = cs%Kv_BBL_min
1063 if (cdrag_sqrt*ustar(i)*g%Rad_Earth*us%m_to_Z > kv_bbl)
then
1064 bbl_thick_z = kv_bbl / ( cdrag_sqrt*ustar(i) )
1066 bbl_thick_z = g%Rad_Earth * us%m_to_Z
1069 kv_bbl = cdrag_sqrt*ustar(i)*bbl_thick_z
1072 kv_bbl = max(cs%Kv_BBL_min, kv_bbl)
1074 visc%Kv_bbl_u(i,j) = kv_bbl
1075 visc%bbl_thick_u(i,j) = bbl_thick_z
1077 visc%Kv_bbl_v(i,j) = kv_bbl
1078 visc%bbl_thick_v(i,j) = bbl_thick_z
1084 if (cs%id_bbl_thick_u > 0) &
1085 call post_data(cs%id_bbl_thick_u, visc%bbl_thick_u, cs%diag)
1086 if (cs%id_kv_bbl_u > 0) &
1087 call post_data(cs%id_kv_bbl_u, visc%kv_bbl_u, cs%diag)
1088 if (cs%id_bbl_u > 0) &
1089 call post_data(cs%id_bbl_u, cs%bbl_u, cs%diag)
1090 if (cs%id_bbl_thick_v > 0) &
1091 call post_data(cs%id_bbl_thick_v, visc%bbl_thick_v, cs%diag)
1092 if (cs%id_kv_bbl_v > 0) &
1093 call post_data(cs%id_kv_bbl_v, visc%kv_bbl_v, cs%diag)
1094 if (cs%id_bbl_v > 0) &
1095 call post_data(cs%id_bbl_v, cs%bbl_v, cs%diag)
1096 if (cs%id_Ray_u > 0) &
1097 call post_data(cs%id_Ray_u, visc%Ray_u, cs%diag)
1098 if (cs%id_Ray_v > 0) &
1099 call post_data(cs%id_Ray_v, visc%Ray_v, cs%diag)
1102 if (
associated(visc%Ray_u) .and.
associated(visc%Ray_v)) &
1103 call uvchksum(
"Ray [uv]", visc%Ray_u, visc%Ray_v, g%HI, haloshift=0, scale=us%Z_to_m*us%s_to_T)
1104 if (
associated(visc%kv_bbl_u) .and.
associated(visc%kv_bbl_v)) &
1105 call uvchksum(
"kv_bbl_[uv]", visc%kv_bbl_u, visc%kv_bbl_v, g%HI, &
1106 haloshift=0, scale=us%Z2_T_to_m2_s, scalar_pair=.true.)
1107 if (
associated(visc%bbl_thick_u) .and.
associated(visc%bbl_thick_v)) &
1108 call uvchksum(
"bbl_thick_[uv]", visc%bbl_thick_u, visc%bbl_thick_v, &
1109 g%HI, haloshift=0, scale=us%Z_to_m, scalar_pair=.true.)
1112 end subroutine set_viscous_bbl
1115 function set_v_at_u(v, h, G, i, j, k, mask2dCv, OBC)
1117 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1119 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1121 integer,
intent(in) :: i
1122 integer,
intent(in) :: j
1123 integer,
intent(in) :: k
1124 real,
dimension(SZI_(G),SZJB_(G)),&
1125 intent(in) :: mask2dcv
1131 real :: hwt(0:1,-1:0)
1133 integer :: i0, j0, i1, j1
1135 do j0 = -1,0 ;
do i0 = 0,1 ; i1 = i+i0 ; j1 = j+j0
1136 hwt(i0,j0) = (h(i1,j1,k) + h(i1,j1+1,k)) * mask2dcv(i1,j1)
1139 if (
associated(obc))
then ;
if (obc%number_of_segments > 0)
then
1140 do j0 = -1,0 ;
do i0 = 0,1 ;
if ((obc%segnum_v(i+i0,j+j0) /= obc_none))
then
1141 i1 = i+i0 ; j1 = j+j0
1142 if (obc%segment(obc%segnum_v(i1,j1))%direction == obc_direction_n)
then
1143 hwt(i0,j0) = 2.0 * h(i1,j1,k) * mask2dcv(i1,j1)
1144 elseif (obc%segment(obc%segnum_v(i1,j1))%direction == obc_direction_s)
then
1145 hwt(i0,j0) = 2.0 * h(i1,j1+1,k) * mask2dcv(i1,j1)
1147 endif ;
enddo ;
enddo
1150 hwt_tot = (hwt(0,-1) + hwt(1,0)) + (hwt(1,-1) + hwt(0,0))
1152 if (hwt_tot > 0.0) set_v_at_u = &
1153 ((hwt(0,0) * v(i,j,k) + hwt(1,-1) * v(i+1,j-1,k)) + &
1154 (hwt(1,0) * v(i+1,j,k) + hwt(0,-1) * v(i,j-1,k))) / hwt_tot
1156 end function set_v_at_u
1159 function set_u_at_v(u, h, G, i, j, k, mask2dCu, OBC)
1161 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1163 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1165 integer,
intent(in) :: i
1166 integer,
intent(in) :: j
1167 integer,
intent(in) :: k
1168 real,
dimension(SZIB_(G),SZJ_(G)), &
1169 intent(in) :: mask2dcu
1175 real :: hwt(-1:0,0:1)
1177 integer :: i0, j0, i1, j1
1179 do j0 = 0,1 ;
do i0 = -1,0 ; i1 = i+i0 ; j1 = j+j0
1180 hwt(i0,j0) = (h(i1,j1,k) + h(i1+1,j1,k)) * mask2dcu(i1,j1)
1183 if (
associated(obc))
then ;
if (obc%number_of_segments > 0)
then
1184 do j0 = 0,1 ;
do i0 = -1,0 ;
if ((obc%segnum_u(i+i0,j+j0) /= obc_none))
then
1185 i1 = i+i0 ; j1 = j+j0
1186 if (obc%segment(obc%segnum_u(i1,j1))%direction == obc_direction_e)
then
1187 hwt(i0,j0) = 2.0 * h(i1,j1,k) * mask2dcu(i1,j1)
1188 elseif (obc%segment(obc%segnum_u(i1,j1))%direction == obc_direction_w)
then
1189 hwt(i0,j0) = 2.0 * h(i1+1,j1,k) * mask2dcu(i1,j1)
1191 endif ;
enddo ;
enddo
1194 hwt_tot = (hwt(-1,0) + hwt(0,1)) + (hwt(0,0) + hwt(-1,1))
1196 if (hwt_tot > 0.0) set_u_at_v = &
1197 ((hwt(0,0) * u(i,j,k) + hwt(-1,1) * u(i-1,j+1,k)) + &
1198 (hwt(-1,0) * u(i-1,j,k) + hwt(0,1) * u(i,j+1,k))) / hwt_tot
1200 end function set_u_at_v
1207 subroutine set_viscous_ml(u, v, h, tv, forces, visc, dt, G, GV, US, CS, symmetrize)
1211 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1213 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1215 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1223 real,
intent(in) :: dt
1226 logical,
optional,
intent(in) :: symmetrize
1231 real,
dimension(SZIB_(G)) :: &
1252 real,
dimension(SZIB_(G),SZJ_(G)) :: &
1255 real,
dimension(SZI_(G),SZJB_(G)) :: &
1258 real :: h_at_vel(szib_(g),szk_(g))
1261 integer :: k_massive(szib_(g))
1298 real :: cdrag_sqrt_z
1322 logical :: use_eos, do_any, do_any_shelf, do_i(szib_(g))
1323 integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz, k2, nkmb, nkml, n
1326 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1327 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1328 nkmb = gv%nk_rho_varies ; nkml = gv%nkml
1330 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_set_viscosity(visc_ML): "//&
1331 "Module must be initialized before it is used.")
1332 if (.not.(cs%dynamic_viscous_ML .or.
associated(forces%frac_shelf_u) .or. &
1333 associated(forces%frac_shelf_v)) )
return
1335 if (
present(symmetrize))
then ;
if (symmetrize)
then
1336 jsq = js-1 ; isq = is-1
1339 rho0x400_g = 400.0*(gv%Rho0/(us%L_to_Z**2 * gv%g_Earth)) * gv%Z_to_H
1340 u_bg_sq = cs%drag_bg_vel * cs%drag_bg_vel
1341 cdrag_sqrt = sqrt(cs%cdrag)
1342 cdrag_sqrt_z = us%L_to_Z * sqrt(cs%cdrag)
1345 use_eos =
associated(tv%eqn_of_state)
1346 dt_rho0 = dt / gv%H_to_RZ
1347 h_neglect = gv%H_subroundoff
1348 h_tiny = 2.0*gv%Angstrom_H + h_neglect
1349 g_h_rho0 = (gv%g_Earth*gv%H_to_Z) / (gv%Rho0)
1351 if (
associated(forces%frac_shelf_u) .neqv.
associated(forces%frac_shelf_v)) &
1352 call mom_error(fatal,
"set_viscous_ML: one of forces%frac_shelf_u and "//&
1353 "forces%frac_shelf_v is associated, but the other is not.")
1355 if (
associated(forces%frac_shelf_u))
then
1358 call safe_alloc_ptr(visc%tauy_shelf, g%isd, g%ied, g%JsdB, g%JedB)
1359 call safe_alloc_ptr(visc%tbl_thick_shelf_u, g%IsdB, g%IedB, g%jsd, g%jed)
1360 call safe_alloc_ptr(visc%tbl_thick_shelf_v, g%isd, g%ied, g%JsdB, g%JedB)
1361 call safe_alloc_ptr(visc%kv_tbl_shelf_u, g%IsdB, g%IedB, g%jsd, g%jed)
1362 call safe_alloc_ptr(visc%kv_tbl_shelf_v, g%isd, g%ied, g%JsdB, g%JedB)
1363 call safe_alloc_ptr(visc%taux_shelf, g%IsdB, g%IedB, g%jsd, g%jed)
1364 call safe_alloc_ptr(visc%tauy_shelf, g%isd, g%ied, g%JsdB, g%JedB)
1371 do j=js-1,je ;
do i=is-1,ie+1
1372 mask_v(i,j) = g%mask2dCv(i,j)
1375 do j=js-1,je+1 ;
do i=is-1,ie
1376 mask_u(i,j) = g%mask2dCu(i,j)
1379 if (
associated(obc))
then ;
do n=1,obc%number_of_segments
1382 if (.not. obc%segment(n)%on_pe) cycle
1384 i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
1385 if (obc%segment(n)%is_N_or_S .and. (j >= js-1) .and. (j <= je))
then
1386 do i = max(is-1,obc%segment(n)%HI%IsdB), min(ie,obc%segment(n)%HI%IedB)
1387 if (obc%segment(n)%direction == obc_direction_n) mask_u(i,j+1) = 0.0
1388 if (obc%segment(n)%direction == obc_direction_s) mask_u(i,j) = 0.0
1390 elseif (obc%segment(n)%is_E_or_W .and. (i >= is-1) .and. (i <= je))
then
1391 do j = max(js-1,obc%segment(n)%HI%JsdB), min(je,obc%segment(n)%HI%JedB)
1392 if (obc%segment(n)%direction == obc_direction_e) mask_v(i+1,j) = 0.0
1393 if (obc%segment(n)%direction == obc_direction_w) mask_v(i,j) = 0.0
1402 if (cs%dynamic_viscous_ML)
then
1406 if (g%mask2dCu(i,j) < 0.5)
then
1407 do_i(i) = .false. ; visc%nkml_visc_u(i,j) = nkml
1409 do_i(i) = .true. ; do_any = .true.
1411 thtot(i) = 0.0 ; shtot(i) = 0.0 ; rhtot(i) = 0.0
1412 uhtot(i) = dt_rho0 * forces%taux(i,j)
1413 vhtot(i) = 0.25 * dt_rho0 * ((forces%tauy(i,j) + forces%tauy(i+1,j-1)) + &
1414 (forces%tauy(i,j-1) + forces%tauy(i+1,j)))
1416 if (cs%omega_frac >= 1.0)
then ; absf = 2.0*cs%omega ;
else
1417 absf = 0.5*(abs(g%CoriolisBu(i,j)) + abs(g%CoriolisBu(i,j-1)))
1418 if (cs%omega_frac > 0.0) &
1419 absf = sqrt(cs%omega_frac*4.0*cs%omega**2 + (1.0-cs%omega_frac)*absf**2)
1421 u_star = max(cs%ustar_min, 0.5 * (forces%ustar(i,j) + forces%ustar(i+1,j)))
1422 idecay_len_tke(i) = ((absf / u_star) * cs%TKE_decay) * gv%H_to_Z
1426 if (do_any)
then ;
do k=1,nz
1429 if (use_eos .and. (k==nkml+1))
then
1432 press(i) = (gv%H_to_RZ*gv%g_Earth) * htot(i)
1433 if (
associated(tv%p_surf)) press(i) = press(i) + 0.5*(tv%p_surf(i,j)+tv%p_surf(i+1,j))
1435 i_2hlay = 1.0 / (h(i,j,k2) + h(i+1,j,k2) + h_neglect)
1436 t_eos(i) = (h(i,j,k2)*tv%T(i,j,k2) + h(i+1,j,k2)*tv%T(i+1,j,k2)) * i_2hlay
1437 s_eos(i) = (h(i,j,k2)*tv%S(i,j,k2) + h(i+1,j,k2)*tv%S(i+1,j,k2)) * i_2hlay
1440 (/isq-g%IsdB+1,ieq-g%IsdB+1/) )
1443 do i=isq,ieq ;
if (do_i(i))
then
1445 hlay = 0.5*(h(i,j,k) + h(i+1,j,k))
1446 if (hlay > h_tiny)
then
1447 i_2hlay = 1.0 / (h(i,j,k) + h(i+1,j,k))
1448 v_at_u = 0.5 * (h(i,j,k) * (v(i,j,k) + v(i,j-1,k)) + &
1449 h(i+1,j,k) * (v(i+1,j,k) + v(i+1,j-1,k))) * i_2hlay
1450 uh2 = ((uhtot(i) - htot(i)*u(i,j,k))**2 + (vhtot(i) - htot(i)*v_at_u)**2)
1453 t_lay = (h(i,j,k)*tv%T(i,j,k) + h(i+1,j,k)*tv%T(i+1,j,k)) * i_2hlay
1454 s_lay = (h(i,j,k)*tv%S(i,j,k) + h(i+1,j,k)*tv%S(i+1,j,k)) * i_2hlay
1455 ghprime = g_h_rho0 * (dr_dt(i) * (t_lay*htot(i) - thtot(i)) + &
1456 dr_ds(i) * (s_lay*htot(i) - shtot(i)))
1458 ghprime = g_h_rho0 * (gv%Rlay(k)*htot(i) - rhtot(i))
1461 if (ghprime > 0.0)
then
1462 ribulk = cs%bulk_Ri_ML * exp(-htot(i) * idecay_len_tke(i))
1463 if (ribulk * uh2 <= (htot(i)**2) * ghprime)
then
1464 visc%nkml_visc_u(i,j) = real(k_massive(i))
1466 elseif (ribulk * uh2 <= (htot(i) + hlay)**2 * ghprime)
then
1467 visc%nkml_visc_u(i,j) = real(k-1) + &
1468 ( sqrt(ribulk * uh2 / ghprime) - htot(i) ) / hlay
1475 if (do_i(i)) do_any = .true.
1478 if (.not.do_any)
exit
1481 do i=isq,ieq ;
if (do_i(i))
then
1482 htot(i) = htot(i) + 0.5 * (h(i,j,k) + h(i+1,j,k))
1483 uhtot(i) = uhtot(i) + 0.5 * (h(i,j,k) + h(i+1,j,k)) * u(i,j,k)
1484 vhtot(i) = vhtot(i) + 0.25 * (h(i,j,k) * (v(i,j,k) + v(i,j-1,k)) + &
1485 h(i+1,j,k) * (v(i+1,j,k) + v(i+1,j-1,k)))
1487 thtot(i) = thtot(i) + 0.5 * (h(i,j,k)*tv%T(i,j,k) + h(i+1,j,k)*tv%T(i+1,j,k))
1488 shtot(i) = shtot(i) + 0.5 * (h(i,j,k)*tv%S(i,j,k) + h(i+1,j,k)*tv%S(i+1,j,k))
1490 rhtot(i) = rhtot(i) + 0.5 * (h(i,j,k) + h(i+1,j,k)) * gv%Rlay(k)
1495 if (do_any)
then ;
do i=isq,ieq ;
if (do_i(i))
then
1496 visc%nkml_visc_u(i,j) = k_massive(i)
1497 endif ;
enddo ;
endif
1500 do_any_shelf = .false.
1501 if (
associated(forces%frac_shelf_u))
then
1503 if (forces%frac_shelf_u(i,j)*g%mask2dCu(i,j) == 0.0)
then
1505 visc%tbl_thick_shelf_u(i,j) = 0.0 ; visc%kv_tbl_shelf_u(i,j) = 0.0
1507 do_i(i) = .true. ; do_any_shelf = .true.
1512 if (do_any_shelf)
then
1513 do k=1,nz ;
do i=isq,ieq ;
if (do_i(i))
then
1514 if (u(i,j,k) *(h(i+1,j,k) - h(i,j,k)) >= 0)
then
1515 h_at_vel(i,k) = 2.0*h(i,j,k)*h(i+1,j,k) / &
1516 (h(i,j,k) + h(i+1,j,k) + h_neglect)
1518 h_at_vel(i,k) = 0.5 * (h(i,j,k) + h(i+1,j,k))
1521 h_at_vel(i,k) = 0.0 ; ustar(i) = 0.0
1522 endif ;
enddo ;
enddo
1524 do i=isq,ieq ;
if (do_i(i))
then
1525 htot_vel = 0.0 ; hwtot = 0.0 ; hutot = 0.0
1526 thtot(i) = 0.0 ; shtot(i) = 0.0
1527 if (use_eos .or. .not.cs%linear_drag)
then ;
do k=1,nz
1528 if (htot_vel>=cs%Htbl_shelf)
exit
1529 hweight = min(cs%Htbl_shelf - htot_vel, h_at_vel(i,k))
1530 if (hweight <= 1.5*gv%Angstrom_H + h_neglect) cycle
1532 htot_vel = htot_vel + h_at_vel(i,k)
1533 hwtot = hwtot + hweight
1535 if (.not.cs%linear_drag)
then
1536 v_at_u = set_v_at_u(v, h, g, i, j, k, mask_v, obc)
1537 hutot = hutot + hweight * sqrt(u(i,j,k)**2 + &
1538 v_at_u**2 + u_bg_sq)
1541 thtot(i) = thtot(i) + hweight * 0.5 * (tv%T(i,j,k) + tv%T(i+1,j,k))
1542 shtot(i) = shtot(i) + hweight * 0.5 * (tv%S(i,j,k) + tv%S(i+1,j,k))
1546 if ((.not.cs%linear_drag) .and. (hwtot > 0.0))
then
1547 ustar(i) = cdrag_sqrt_z * hutot/hwtot
1549 ustar(i) = cdrag_sqrt_z * cs%drag_bg_vel
1552 if (use_eos)
then ;
if (hwtot > 0.0)
then
1553 t_eos(i) = thtot(i)/hwtot ; s_eos(i) = shtot(i)/hwtot
1555 t_eos(i) = 0.0 ; s_eos(i) = 0.0
1561 tv%eqn_of_state, (/isq-g%IsdB+1,ieq-g%IsdB+1/) )
1564 do i=isq,ieq ;
if (do_i(i))
then
1567 ustarsq = rho0x400_g * ustar(i)**2
1570 thtot(i) = 0.0 ; shtot(i) = 0.0
1572 if (h_at_vel(i,k) <= 0.0) cycle
1573 t_lay = 0.5 * (tv%T(i,j,k) + tv%T(i+1,j,k))
1574 s_lay = 0.5 * (tv%S(i,j,k) + tv%S(i+1,j,k))
1575 oldfn = dr_dt(i)*(t_lay*htot(i) - thtot(i)) + dr_ds(i)*(s_lay*htot(i) - shtot(i))
1576 if (oldfn >= ustarsq)
exit
1578 dfn = (dr_dt(i)*(0.5*(tv%T(i,j,k+1)+tv%T(i+1,j,k+1)) - t_lay) + &
1579 dr_ds(i)*(0.5*(tv%S(i,j,k+1)+tv%S(i+1,j,k+1)) - s_lay)) * &
1580 (h_at_vel(i,k)+htot(i))
1581 if ((oldfn + dfn) <= ustarsq)
then
1584 dh = h_at_vel(i,k) * sqrt((ustarsq-oldfn) / (dfn))
1587 htot(i) = htot(i) + dh
1588 thtot(i) = thtot(i) + t_lay*dh ; shtot(i) = shtot(i) + s_lay*dh
1590 if ((oldfn < ustarsq) .and. (h_at_vel(i,nz) > 0.0))
then
1591 t_lay = 0.5*(tv%T(i,j,nz) + tv%T(i+1,j,nz))
1592 s_lay = 0.5*(tv%S(i,j,nz) + tv%S(i+1,j,nz))
1593 if (dr_dt(i)*(t_lay*htot(i) - thtot(i)) + &
1594 dr_ds(i)*(s_lay*htot(i) - shtot(i)) < ustarsq) &
1595 htot(i) = htot(i) + h_at_vel(i,nz)
1600 rlay = gv%Rlay(k) ; rlb = gv%Rlay(k+1)
1602 oldfn = rlay*htot(i) - rhtot(i)
1603 if (oldfn >= ustarsq)
exit
1605 dfn = (rlb - rlay)*(h_at_vel(i,k)+htot(i))
1606 if ((oldfn + dfn) <= ustarsq)
then
1609 dh = h_at_vel(i,k) * sqrt((ustarsq-oldfn) / (dfn))
1612 htot(i) = htot(i) + dh
1613 rhtot(i) = rhtot(i) + rlay*dh
1615 if (gv%Rlay(nz)*htot(i) - rhtot(i) < ustarsq) &
1616 htot(i) = htot(i) + h_at_vel(i,nz)
1623 ustar1 = ustar(i)*gv%Z_to_H
1624 h2f2 = (htot(i)*(g%CoriolisBu(i,j-1)+g%CoriolisBu(i,j)) + h_neglect*cs%omega)**2
1625 tbl_thick_z = gv%H_to_Z * max(cs%Htbl_shelf_min, &
1626 ( htot(i)*ustar1 ) / ( 0.5*ustar1 + sqrt((0.5*ustar1)**2 + h2f2 ) ) )
1627 visc%tbl_thick_shelf_u(i,j) = tbl_thick_z
1628 visc%Kv_tbl_shelf_u(i,j) = max(cs%Kv_TBL_min, cdrag_sqrt*ustar(i)*tbl_thick_z)
1638 if (cs%dynamic_viscous_ML)
then
1642 if (g%mask2dCv(i,j) < 0.5)
then
1643 do_i(i) = .false. ; visc%nkml_visc_v(i,j) = nkml
1645 do_i(i) = .true. ; do_any = .true.
1647 thtot(i) = 0.0 ; shtot(i) = 0.0 ; rhtot(i) = 0.0
1648 vhtot(i) = dt_rho0 * forces%tauy(i,j)
1649 uhtot(i) = 0.25 * dt_rho0 * ((forces%taux(i,j) + forces%taux(i-1,j+1)) + &
1650 (forces%taux(i-1,j) + forces%taux(i,j+1)))
1652 if (cs%omega_frac >= 1.0)
then ; absf = 2.0*cs%omega ;
else
1653 absf = 0.5*(abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j)))
1654 if (cs%omega_frac > 0.0) &
1655 absf = sqrt(cs%omega_frac*4.0*cs%omega**2 + (1.0-cs%omega_frac)*absf**2)
1658 u_star = max(cs%ustar_min, 0.5 * (forces%ustar(i,j) + forces%ustar(i,j+1)))
1659 idecay_len_tke(i) = ((absf / u_star) * cs%TKE_decay) * gv%H_to_Z
1664 if (do_any)
then ;
do k=1,nz
1667 if (use_eos .and. (k==nkml+1))
then
1670 press(i) = (gv%H_to_RZ * gv%g_Earth) * htot(i)
1671 if (
associated(tv%p_surf)) press(i) = press(i) + 0.5*(tv%p_surf(i,j)+tv%p_surf(i,j+1))
1673 i_2hlay = 1.0 / (h(i,j,k2) + h(i,j+1,k2) + h_neglect)
1674 t_eos(i) = (h(i,j,k2)*tv%T(i,j,k2) + h(i,j+1,k2)*tv%T(i,j+1,k2)) * i_2hlay
1675 s_eos(i) = (h(i,j,k2)*tv%S(i,j,k2) + h(i,j+1,k2)*tv%S(i,j+1,k2)) * i_2hlay
1678 tv%eqn_of_state, (/is-g%IsdB+1,ie-g%IsdB+1/) )
1681 do i=is,ie ;
if (do_i(i))
then
1683 hlay = 0.5*(h(i,j,k) + h(i,j+1,k))
1684 if (hlay > h_tiny)
then
1685 i_2hlay = 1.0 / (h(i,j,k) + h(i,j+1,k))
1686 u_at_v = 0.5 * (h(i,j,k) * (u(i-1,j,k) + u(i,j,k)) + &
1687 h(i,j+1,k) * (u(i-1,j+1,k) + u(i,j+1,k))) * i_2hlay
1688 uh2 = ((uhtot(i) - htot(i)*u_at_v)**2 + (vhtot(i) - htot(i)*v(i,j,k))**2)
1691 t_lay = (h(i,j,k)*tv%T(i,j,k) + h(i,j+1,k)*tv%T(i,j+1,k)) * i_2hlay
1692 s_lay = (h(i,j,k)*tv%S(i,j,k) + h(i,j+1,k)*tv%S(i,j+1,k)) * i_2hlay
1693 ghprime = g_h_rho0 * (dr_dt(i) * (t_lay*htot(i) - thtot(i)) + &
1694 dr_ds(i) * (s_lay*htot(i) - shtot(i)))
1696 ghprime = g_h_rho0 * (gv%Rlay(k)*htot(i) - rhtot(i))
1699 if (ghprime > 0.0)
then
1700 ribulk = cs%bulk_Ri_ML * exp(-htot(i) * idecay_len_tke(i))
1701 if (ribulk * uh2 <= htot(i)**2 * ghprime)
then
1702 visc%nkml_visc_v(i,j) = real(k_massive(i))
1704 elseif (ribulk * uh2 <= (htot(i) + hlay)**2 * ghprime)
then
1705 visc%nkml_visc_v(i,j) = real(k-1) + &
1706 ( sqrt(ribulk * uh2 / ghprime) - htot(i) ) / hlay
1713 if (do_i(i)) do_any = .true.
1716 if (.not.do_any)
exit
1719 do i=is,ie ;
if (do_i(i))
then
1720 htot(i) = htot(i) + 0.5 * (h(i,j,k) + h(i,j+1,k))
1721 vhtot(i) = vhtot(i) + 0.5 * (h(i,j,k) + h(i,j+1,k)) * v(i,j,k)
1722 uhtot(i) = uhtot(i) + 0.25 * (h(i,j,k) * (u(i-1,j,k) + u(i,j,k)) + &
1723 h(i,j+1,k) * (u(i-1,j+1,k) + u(i,j+1,k)))
1725 thtot(i) = thtot(i) + 0.5 * (h(i,j,k)*tv%T(i,j,k) + h(i,j+1,k)*tv%T(i,j+1,k))
1726 shtot(i) = shtot(i) + 0.5 * (h(i,j,k)*tv%S(i,j,k) + h(i,j+1,k)*tv%S(i,j+1,k))
1728 rhtot(i) = rhtot(i) + 0.5 * (h(i,j,k) + h(i,j+1,k)) * gv%Rlay(k)
1733 if (do_any)
then ;
do i=is,ie ;
if (do_i(i))
then
1734 visc%nkml_visc_v(i,j) = k_massive(i)
1735 endif ;
enddo ;
endif
1738 do_any_shelf = .false.
1739 if (
associated(forces%frac_shelf_v))
then
1741 if (forces%frac_shelf_v(i,j)*g%mask2dCv(i,j) == 0.0)
then
1743 visc%tbl_thick_shelf_v(i,j) = 0.0 ; visc%kv_tbl_shelf_v(i,j) = 0.0
1745 do_i(i) = .true. ; do_any_shelf = .true.
1750 if (do_any_shelf)
then
1751 do k=1,nz ;
do i=is,ie ;
if (do_i(i))
then
1752 if (v(i,j,k) * (h(i,j+1,k) - h(i,j,k)) >= 0)
then
1753 h_at_vel(i,k) = 2.0*h(i,j,k)*h(i,j+1,k) / &
1754 (h(i,j,k) + h(i,j+1,k) + h_neglect)
1756 h_at_vel(i,k) = 0.5 * (h(i,j,k) + h(i,j+1,k))
1759 h_at_vel(i,k) = 0.0 ; ustar(i) = 0.0
1760 endif ;
enddo ;
enddo
1762 do i=is,ie ;
if (do_i(i))
then
1763 htot_vel = 0.0 ; hwtot = 0.0 ; hutot = 0.0
1764 thtot(i) = 0.0 ; shtot(i) = 0.0
1765 if (use_eos .or. .not.cs%linear_drag)
then ;
do k=1,nz
1766 if (htot_vel>=cs%Htbl_shelf)
exit
1767 hweight = min(cs%Htbl_shelf - htot_vel, h_at_vel(i,k))
1768 if (hweight <= 1.5*gv%Angstrom_H + h_neglect) cycle
1770 htot_vel = htot_vel + h_at_vel(i,k)
1771 hwtot = hwtot + hweight
1773 if (.not.cs%linear_drag)
then
1774 u_at_v = set_u_at_v(u, h, g, i, j, k, mask_u, obc)
1775 hutot = hutot + hweight * sqrt(v(i,j,k)**2 + &
1776 u_at_v**2 + u_bg_sq)
1779 thtot(i) = thtot(i) + hweight * 0.5 * (tv%T(i,j,k) + tv%T(i,j+1,k))
1780 shtot(i) = shtot(i) + hweight * 0.5 * (tv%S(i,j,k) + tv%S(i,j+1,k))
1784 if (.not.cs%linear_drag)
then ;
if (hwtot > 0.0)
then
1785 ustar(i) = cdrag_sqrt_z * hutot/hwtot
1787 ustar(i) = cdrag_sqrt_z * cs%drag_bg_vel
1790 if (use_eos)
then ;
if (hwtot > 0.0)
then
1791 t_eos(i) = thtot(i)/hwtot ; s_eos(i) = shtot(i)/hwtot
1793 t_eos(i) = 0.0 ; s_eos(i) = 0.0
1799 tv%eqn_of_state, (/is-g%IsdB+1,ie-g%IsdB+1/) )
1802 do i=is,ie ;
if (do_i(i))
then
1805 ustarsq = rho0x400_g * ustar(i)**2
1808 thtot(i) = 0.0 ; shtot(i) = 0.0
1810 if (h_at_vel(i,k) <= 0.0) cycle
1811 t_lay = 0.5 * (tv%T(i,j,k) + tv%T(i,j+1,k))
1812 s_lay = 0.5 * (tv%S(i,j,k) + tv%S(i,j+1,k))
1813 oldfn = dr_dt(i)*(t_lay*htot(i) - thtot(i)) + dr_ds(i)*(s_lay*htot(i) - shtot(i))
1814 if (oldfn >= ustarsq)
exit
1816 dfn = (dr_dt(i)*(0.5*(tv%T(i,j,k+1)+tv%T(i,j+1,k+1)) - t_lay) + &
1817 dr_ds(i)*(0.5*(tv%S(i,j,k+1)+tv%S(i,j+1,k+1)) - s_lay)) * &
1818 (h_at_vel(i,k)+htot(i))
1819 if ((oldfn + dfn) <= ustarsq)
then
1822 dh = h_at_vel(i,k) * sqrt((ustarsq-oldfn) / (dfn))
1825 htot(i) = htot(i) + dh
1826 thtot(i) = thtot(i) + t_lay*dh ; shtot(i) = shtot(i) + s_lay*dh
1828 if ((oldfn < ustarsq) .and. (h_at_vel(i,nz) > 0.0))
then
1829 t_lay = 0.5*(tv%T(i,j,nz) + tv%T(i,j+1,nz))
1830 s_lay = 0.5*(tv%S(i,j,nz) + tv%S(i,j+1,nz))
1831 if (dr_dt(i)*(t_lay*htot(i) - thtot(i)) + &
1832 dr_ds(i)*(s_lay*htot(i) - shtot(i)) < ustarsq) &
1833 htot(i) = htot(i) + h_at_vel(i,nz)
1838 rlay = gv%Rlay(k) ; rlb = gv%Rlay(k+1)
1840 oldfn = rlay*htot(i) - rhtot(i)
1841 if (oldfn >= ustarsq)
exit
1843 dfn = (rlb - rlay)*(h_at_vel(i,k)+htot(i))
1844 if ((oldfn + dfn) <= ustarsq)
then
1847 dh = h_at_vel(i,k) * sqrt((ustarsq-oldfn) / (dfn))
1850 htot(i) = htot(i) + dh
1851 rhtot = rhtot + rlay*dh
1853 if (gv%Rlay(nz)*htot(i) - rhtot(i) < ustarsq) &
1854 htot(i) = htot(i) + h_at_vel(i,nz)
1861 ustar1 = ustar(i)*gv%Z_to_H
1862 h2f2 = (htot(i)*(g%CoriolisBu(i-1,j)+g%CoriolisBu(i,j)) + h_neglect*cs%omega)**2
1863 tbl_thick_z = gv%H_to_Z * max(cs%Htbl_shelf_min, &
1864 ( htot(i)*ustar1 ) / ( 0.5*ustar1 + sqrt((0.5*ustar1)**2 + h2f2 ) ) )
1865 visc%tbl_thick_shelf_v(i,j) = tbl_thick_z
1866 visc%Kv_tbl_shelf_v(i,j) = max(cs%Kv_TBL_min, cdrag_sqrt*ustar(i)*tbl_thick_z)
1874 if (
associated(visc%nkml_visc_u) .and.
associated(visc%nkml_visc_v)) &
1875 call uvchksum(
"nkml_visc_[uv]", visc%nkml_visc_u, visc%nkml_visc_v, &
1876 g%HI, haloshift=0, scalar_pair=.true.)
1878 if (cs%id_nkml_visc_u > 0) &
1879 call post_data(cs%id_nkml_visc_u, visc%nkml_visc_u, cs%diag)
1880 if (cs%id_nkml_visc_v > 0) &
1881 call post_data(cs%id_nkml_visc_v, visc%nkml_visc_v, cs%diag)
1883 end subroutine set_viscous_ml
1886 subroutine set_visc_register_restarts(HI, GV, param_file, visc, restart_CS)
1896 logical :: use_kappa_shear, ks_at_vertex
1897 logical :: adiabatic, usekpp, useepbl
1898 logical :: use_cvmix_shear, mle_use_pbl_mld, use_cvmix_conv
1899 integer :: isd, ied, jsd, jed, nz
1901 character(len=40) :: mdl =
"MOM_set_visc"
1902 isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
1904 call get_param(param_file, mdl,
"ADIABATIC", adiabatic, default=.false., &
1907 use_kappa_shear = .false. ; ks_at_vertex = .false. ; use_cvmix_shear = .false.
1908 usekpp = .false. ; useepbl = .false. ; use_cvmix_conv = .false.
1910 if (.not.adiabatic)
then
1911 use_kappa_shear = kappa_shear_is_used(param_file)
1912 ks_at_vertex = kappa_shear_at_vertex(param_file)
1913 use_cvmix_shear = cvmix_shear_is_used(param_file)
1914 use_cvmix_conv = cvmix_conv_is_used(param_file)
1915 call get_param(param_file, mdl,
"USE_KPP", usekpp, &
1916 "If true, turns on the [CVMix] KPP scheme of Large et al., 1994, "//&
1917 "to calculate diffusivities and non-local transport in the OBL.", &
1918 default=.false., do_not_log=.true.)
1919 call get_param(param_file, mdl,
"ENERGETICS_SFC_PBL", useepbl, &
1920 "If true, use an implied energetics planetary boundary "//&
1921 "layer scheme to determine the diffusivity and viscosity "//&
1922 "in the surface boundary layer.", default=.false., do_not_log=.true.)
1925 if (use_kappa_shear .or. usekpp .or. useepbl .or. use_cvmix_shear .or. use_cvmix_conv)
then
1928 "Shear-driven turbulent diffusivity at interfaces",
"m2 s-1", z_grid=
'i')
1930 if (usekpp .or. useepbl .or. use_cvmix_shear .or. use_cvmix_conv .or. &
1931 (use_kappa_shear .and. .not.ks_at_vertex ))
then
1934 "Shear-driven turbulent viscosity at interfaces",
"m2 s-1", z_grid=
'i')
1936 if (use_kappa_shear .and. ks_at_vertex)
then
1937 call safe_alloc_ptr(visc%TKE_turb, hi%IsdB, hi%IedB, hi%JsdB, hi%JedB, nz+1)
1938 call safe_alloc_ptr(visc%Kv_shear_Bu, hi%IsdB, hi%IedB, hi%JsdB, hi%JedB, nz+1)
1940 "Shear-driven turbulent viscosity at vertex interfaces",
"m2 s-1", &
1941 hor_grid=
"Bu", z_grid=
'i')
1942 elseif (use_kappa_shear)
then
1952 call get_param(param_file, mdl,
"MLE_USE_PBL_MLD", mle_use_pbl_mld, &
1953 default=.false., do_not_log=.true.)
1955 call get_param(param_file, mdl,
"HFREEZE", hfreeze, &
1956 default=-1.0, do_not_log=.true.)
1958 if (mle_use_pbl_mld)
then
1961 "Instantaneous active mixing layer depth",
"m")
1964 if (hfreeze >= 0.0 .and. .not.mle_use_pbl_mld)
then
1969 end subroutine set_visc_register_restarts
1972 subroutine set_visc_init(Time, G, GV, US, param_file, diag, visc, CS, restart_CS, OBC)
1973 type(time_type),
target,
intent(in) :: time
1979 type(
diag_ctrl),
target,
intent(inout) :: diag
1989 real :: csmag_chan_dflt, smag_const1, tke_decay_dflt, bulk_ri_ml_dflt
1990 real :: kv_background
1991 real :: omega_frac_dflt
1996 real :: z2_t_rescale
1998 integer :: i, j, k, is, ie, js, je, n
1999 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, nz
2000 logical :: default_2018_answers
2001 logical :: use_kappa_shear, adiabatic, use_omega, mle_use_pbl_mld
2003 logical :: use_regridding
2004 character(len=200) :: filename, tideamp_file
2007 # include "version_variable.h"
2008 character(len=40) :: mdl =
"MOM_set_visc"
2010 if (
associated(cs))
then
2011 call mom_error(warning,
"set_visc_init called with an associated "// &
2012 "control structure.")
2019 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
2020 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = gv%ke
2021 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
2027 cs%RiNo_mix = .false.
2028 call get_param(param_file, mdl,
"INPUTDIR", cs%inputdir, default=
".")
2029 cs%inputdir = slasher(cs%inputdir)
2030 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
2031 "This sets the default value for the various _2018_ANSWERS parameters.", &
2033 call get_param(param_file, mdl,
"SET_VISC_2018_ANSWERS", cs%answers_2018, &
2034 "If true, use the order of arithmetic and expressions that recover the "//&
2035 "answers from the end of 2018. Otherwise, use updated and more robust "//&
2036 "forms of the same expressions.", default=default_2018_answers)
2037 call get_param(param_file, mdl,
"BOTTOMDRAGLAW", cs%bottomdraglaw, &
2038 "If true, the bottom stress is calculated with a drag "//&
2039 "law of the form c_drag*|u|*u. The velocity magnitude "//&
2040 "may be an assumed value or it may be based on the "//&
2041 "actual velocity in the bottommost HBBL, depending on "//&
2042 "LINEAR_DRAG.", default=.true.)
2043 call get_param(param_file, mdl,
"CHANNEL_DRAG", cs%Channel_drag, &
2044 "If true, the bottom drag is exerted directly on each "//&
2045 "layer proportional to the fraction of the bottom it "//&
2046 "overlies.", default=.false.)
2047 call get_param(param_file, mdl,
"LINEAR_DRAG", cs%linear_drag, &
2048 "If LINEAR_DRAG and BOTTOMDRAGLAW are defined the drag "//&
2049 "law is cdrag*DRAG_BG_VEL*u.", default=.false.)
2050 call get_param(param_file, mdl,
"ADIABATIC", adiabatic, default=.false., &
2053 call log_param(param_file, mdl,
"ADIABATIC",adiabatic, &
2054 "There are no diapycnal mass fluxes if ADIABATIC is "//&
2055 "true. This assumes that KD = KDML = 0.0 and that "//&
2056 "there is no buoyancy forcing, but makes the model "//&
2057 "faster by eliminating subroutine calls.", default=.false.)
2060 if (.not.adiabatic)
then
2061 cs%RiNo_mix = kappa_shear_is_used(param_file)
2064 call get_param(param_file, mdl,
"PRANDTL_TURB", visc%Prandtl_turb, &
2065 "The turbulent Prandtl number applied to shear "//&
2066 "instability.", units=
"nondim", default=1.0)
2067 call get_param(param_file, mdl,
"DEBUG", cs%debug, default=.false.)
2069 call get_param(param_file, mdl,
"DYNAMIC_VISCOUS_ML", cs%dynamic_viscous_ML, &
2070 "If true, use a bulk Richardson number criterion to "//&
2071 "determine the mixed layer thickness for viscosity.", &
2073 if (cs%dynamic_viscous_ML)
then
2074 call get_param(param_file, mdl,
"BULK_RI_ML", bulk_ri_ml_dflt, default=0.0)
2075 call get_param(param_file, mdl,
"BULK_RI_ML_VISC", cs%bulk_Ri_ML, &
2076 "The efficiency with which mean kinetic energy released "//&
2077 "by mechanically forced entrainment of the mixed layer "//&
2078 "is converted to turbulent kinetic energy. By default, "//&
2079 "BULK_RI_ML_VISC = BULK_RI_ML or 0.", units=
"nondim", &
2080 default=bulk_ri_ml_dflt)
2081 call get_param(param_file, mdl,
"TKE_DECAY", tke_decay_dflt, default=0.0)
2082 call get_param(param_file, mdl,
"TKE_DECAY_VISC", cs%TKE_decay, &
2083 "TKE_DECAY_VISC relates the vertical rate of decay of "//&
2084 "the TKE available for mechanical entrainment to the "//&
2085 "natural Ekman depth for use in calculating the dynamic "//&
2086 "mixed layer viscosity. By default, "//&
2087 "TKE_DECAY_VISC = TKE_DECAY or 0.", units=
"nondim", &
2088 default=tke_decay_dflt)
2089 call get_param(param_file, mdl,
"ML_USE_OMEGA", use_omega, &
2090 "If true, use the absolute rotation rate instead of the "//&
2091 "vertical component of rotation when setting the decay "//&
2092 "scale for turbulence.", default=.false., do_not_log=.true.)
2093 omega_frac_dflt = 0.0
2095 call mom_error(warning,
"ML_USE_OMEGA is depricated; use ML_OMEGA_FRAC=1.0 instead.")
2096 omega_frac_dflt = 1.0
2098 call get_param(param_file, mdl,
"ML_OMEGA_FRAC", cs%omega_frac, &
2099 "When setting the decay scale for turbulence, use this "//&
2100 "fraction of the absolute rotation rate blended with the "//&
2101 "local value of f, as sqrt((1-of)*f^2 + of*4*omega^2).", &
2102 units=
"nondim", default=omega_frac_dflt)
2103 call get_param(param_file, mdl,
"OMEGA", cs%omega, &
2104 "The rotation rate of the earth.", units=
"s-1", &
2105 default=7.2921e-5, scale=us%T_to_s)
2107 cs%ustar_min = 2e-4*cs%omega*(gv%Angstrom_Z + gv%H_to_Z*gv%H_subroundoff)
2109 call get_param(param_file, mdl,
"OMEGA", cs%omega, &
2110 "The rotation rate of the earth.", units=
"s-1", &
2111 default=7.2921e-5, scale=us%T_to_s)
2114 call get_param(param_file, mdl,
"HBBL", cs%Hbbl, &
2115 "The thickness of a bottom boundary layer with a "//&
2116 "viscosity of KVBBL if BOTTOMDRAGLAW is not defined, or "//&
2117 "the thickness over which near-bottom velocities are "//&
2118 "averaged for the drag law if BOTTOMDRAGLAW is defined "//&
2119 "but LINEAR_DRAG is not.", units=
"m", fail_if_missing=.true.)
2120 if (cs%bottomdraglaw)
then
2121 call get_param(param_file, mdl,
"CDRAG", cs%cdrag, &
2122 "CDRAG is the drag coefficient relating the magnitude of "//&
2123 "the velocity field to the bottom stress. CDRAG is only "//&
2124 "used if BOTTOMDRAGLAW is defined.", units=
"nondim", &
2126 call get_param(param_file, mdl,
"BBL_USE_TIDAL_BG", cs%BBL_use_tidal_bg, &
2127 "Flag to use the tidal RMS amplitude in place of constant "//&
2128 "background velocity for computing u* in the BBL. "//&
2129 "This flag is only used when BOTTOMDRAGLAW is true and "//&
2130 "LINEAR_DRAG is false.", default=.false.)
2131 if (cs%BBL_use_tidal_bg)
then
2132 call get_param(param_file, mdl,
"TIDEAMP_FILE", tideamp_file, &
2133 "The path to the file containing the spatially varying "//&
2134 "tidal amplitudes with INT_TIDE_DISSIPATION.", default=
"tideamp.nc")
2136 call get_param(param_file, mdl,
"DRAG_BG_VEL", cs%drag_bg_vel, &
2137 "DRAG_BG_VEL is either the assumed bottom velocity (with "//&
2138 "LINEAR_DRAG) or an unresolved velocity that is "//&
2139 "combined with the resolved velocity to estimate the "//&
2140 "velocity magnitude. DRAG_BG_VEL is only used when "//&
2141 "BOTTOMDRAGLAW is defined.", units=
"m s-1", default=0.0, scale=us%m_s_to_L_T)
2143 call get_param(param_file, mdl,
"USE_REGRIDDING", use_regridding, &
2144 do_not_log = .true., default = .false. )
2145 call get_param(param_file, mdl,
"BBL_USE_EOS", cs%BBL_use_EOS, &
2146 "If true, use the equation of state in determining the "//&
2147 "properties of the bottom boundary layer. Otherwise use "//&
2148 "the layer target potential densities. The default of "//&
2149 "this is determined by USE_REGRIDDING.", default=use_regridding)
2150 if (use_regridding .and. (.not. cs%BBL_use_EOS)) &
2151 call mom_error(fatal,
"When using MOM6 in ALE mode it is required to "//&
2152 "set BBL_USE_EOS to True")
2154 call get_param(param_file, mdl,
"BBL_THICK_MIN", cs%BBL_thick_min, &
2155 "The minimum bottom boundary layer thickness that can be "//&
2156 "used with BOTTOMDRAGLAW. This might be "//&
2157 "Kv/(cdrag*drag_bg_vel) to give Kv as the minimum "//&
2158 "near-bottom viscosity.", units=
"m", default=0.0)
2159 call get_param(param_file, mdl,
"HTBL_SHELF_MIN", cs%Htbl_shelf_min, &
2160 "The minimum top boundary layer thickness that can be "//&
2161 "used with BOTTOMDRAGLAW. This might be "//&
2162 "Kv/(cdrag*drag_bg_vel) to give Kv as the minimum "//&
2163 "near-top viscosity.", units=
"m", default=cs%BBL_thick_min, scale=gv%m_to_H)
2164 call get_param(param_file, mdl,
"HTBL_SHELF", cs%Htbl_shelf, &
2165 "The thickness over which near-surface velocities are "//&
2166 "averaged for the drag law under an ice shelf. By "//&
2167 "default this is the same as HBBL", units=
"m", default=cs%Hbbl, scale=gv%m_to_H)
2169 cs%Hbbl = cs%Hbbl * gv%m_to_H
2170 cs%BBL_thick_min = cs%BBL_thick_min * gv%m_to_H
2172 call get_param(param_file, mdl,
"KV", kv_background, &
2173 "The background kinematic viscosity in the interior. "//&
2174 "The molecular value, ~1e-6 m2 s-1, may be used.", &
2175 units=
"m2 s-1", fail_if_missing=.true.)
2177 call get_param(param_file, mdl,
"USE_KPP", use_kpp, &
2178 "If true, turns on the [CVMix] KPP scheme of Large et al., 1994, "//&
2179 "to calculate diffusivities and non-local transport in the OBL.", &
2180 do_not_log=.true., default=.false.)
2182 call get_param(param_file, mdl,
"KV_BBL_MIN", cs%KV_BBL_min, &
2183 "The minimum viscosities in the bottom boundary layer.", &
2184 units=
"m2 s-1", default=kv_background, scale=us%m2_s_to_Z2_T)
2185 call get_param(param_file, mdl,
"KV_TBL_MIN", cs%KV_TBL_min, &
2186 "The minimum viscosities in the top boundary layer.", &
2187 units=
"m2 s-1", default=kv_background, scale=us%m2_s_to_Z2_T)
2188 call get_param(param_file, mdl,
"CORRECT_BBL_BOUNDS", cs%correct_BBL_bounds, &
2189 "If true, uses the correct bounds on the BBL thickness and "//&
2190 "viscosity so that the bottom layer feels the intended drag.", &
2193 if (cs%Channel_drag)
then
2194 call get_param(param_file, mdl,
"SMAG_LAP_CONST", smag_const1, default=-1.0)
2196 csmag_chan_dflt = 0.15
2197 if (smag_const1 >= 0.0) csmag_chan_dflt = smag_const1
2199 call get_param(param_file, mdl,
"SMAG_CONST_CHANNEL", cs%c_Smag, &
2200 "The nondimensional Laplacian Smagorinsky constant used "//&
2201 "in calculating the channel drag if it is enabled. The "//&
2202 "default is to use the same value as SMAG_LAP_CONST if "//&
2203 "it is defined, or 0.15 if it is not. The value used is "//&
2204 "also 0.15 if the specified value is negative.", &
2205 units=
"nondim", default=csmag_chan_dflt)
2206 if (cs%c_Smag < 0.0) cs%c_Smag = 0.15
2209 call get_param(param_file, mdl,
"MLE_USE_PBL_MLD", mle_use_pbl_mld, &
2210 default=.false., do_not_log=.true.)
2212 if (cs%RiNo_mix .and. kappa_shear_at_vertex(param_file))
then
2214 call pass_var(visc%Kv_shear_Bu, g%Domain, position=corner, complete=.true.)
2217 if (cs%bottomdraglaw)
then
2218 allocate(visc%bbl_thick_u(isdb:iedb,jsd:jed)) ; visc%bbl_thick_u(:,:) = 0.0
2219 allocate(visc%kv_bbl_u(isdb:iedb,jsd:jed)) ; visc%kv_bbl_u(:,:) = 0.0
2220 allocate(visc%bbl_thick_v(isd:ied,jsdb:jedb)) ; visc%bbl_thick_v(:,:) = 0.0
2221 allocate(visc%kv_bbl_v(isd:ied,jsdb:jedb)) ; visc%kv_bbl_v(:,:) = 0.0
2222 allocate(visc%ustar_bbl(isd:ied,jsd:jed)) ; visc%ustar_bbl(:,:) = 0.0
2223 allocate(visc%TKE_bbl(isd:ied,jsd:jed)) ; visc%TKE_bbl(:,:) = 0.0
2225 cs%id_bbl_thick_u = register_diag_field(
'ocean_model',
'bbl_thick_u', &
2226 diag%axesCu1, time,
'BBL thickness at u points',
'm', conversion=us%Z_to_m)
2227 cs%id_kv_bbl_u = register_diag_field(
'ocean_model',
'kv_bbl_u', diag%axesCu1, &
2228 time,
'BBL viscosity at u points',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
2229 cs%id_bbl_u = register_diag_field(
'ocean_model',
'bbl_u', diag%axesCu1, &
2230 time,
'BBL mean u current',
'm s-1', conversion=us%L_T_to_m_s)
2231 if (cs%id_bbl_u>0)
then
2232 allocate(cs%bbl_u(isdb:iedb,jsd:jed)) ; cs%bbl_u(:,:) = 0.0
2234 cs%id_bbl_thick_v = register_diag_field(
'ocean_model',
'bbl_thick_v', &
2235 diag%axesCv1, time,
'BBL thickness at v points',
'm', conversion=us%Z_to_m)
2236 cs%id_kv_bbl_v = register_diag_field(
'ocean_model',
'kv_bbl_v', diag%axesCv1, &
2237 time,
'BBL viscosity at v points',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
2238 cs%id_bbl_v = register_diag_field(
'ocean_model',
'bbl_v', diag%axesCv1, &
2239 time,
'BBL mean v current',
'm s-1', conversion=us%L_T_to_m_s)
2240 if (cs%id_bbl_v>0)
then
2241 allocate(cs%bbl_v(isd:ied,jsdb:jedb)) ; cs%bbl_v(:,:) = 0.0
2243 if (cs%BBL_use_tidal_bg)
then
2244 allocate(cs%tideamp(isd:ied,jsd:jed)) ; cs%tideamp(:,:) = 0.0
2245 filename = trim(cs%inputdir) // trim(tideamp_file)
2246 call log_param(param_file, mdl,
"INPUTDIR/TIDEAMP_FILE", filename)
2247 call mom_read_data(filename,
'tideamp', cs%tideamp, g%domain, timelevel=1, scale=us%m_to_Z*us%T_to_s)
2251 if (cs%Channel_drag)
then
2252 allocate(visc%Ray_u(isdb:iedb,jsd:jed,nz)) ; visc%Ray_u(:,:,:) = 0.0
2253 allocate(visc%Ray_v(isd:ied,jsdb:jedb,nz)) ; visc%Ray_v(:,:,:) = 0.0
2254 cs%id_Ray_u = register_diag_field(
'ocean_model',
'Rayleigh_u', diag%axesCuL, &
2255 time,
'Rayleigh drag velocity at u points',
'm s-1', conversion=us%Z_to_m*us%s_to_T)
2256 cs%id_Ray_v = register_diag_field(
'ocean_model',
'Rayleigh_v', diag%axesCvL, &
2257 time,
'Rayleigh drag velocity at v points',
'm s-1', conversion=us%Z_to_m*us%s_to_T)
2261 if (cs%dynamic_viscous_ML)
then
2262 allocate(visc%nkml_visc_u(isdb:iedb,jsd:jed)) ; visc%nkml_visc_u(:,:) = 0.0
2263 allocate(visc%nkml_visc_v(isd:ied,jsdb:jedb)) ; visc%nkml_visc_v(:,:) = 0.0
2264 cs%id_nkml_visc_u = register_diag_field(
'ocean_model',
'nkml_visc_u', &
2265 diag%axesCu1, time,
'Number of layers in viscous mixed layer at u points',
'm')
2266 cs%id_nkml_visc_v = register_diag_field(
'ocean_model',
'nkml_visc_v', &
2267 diag%axesCv1, time,
'Number of layers in viscous mixed layer at v points',
'm')
2270 call register_restart_field_as_obsolete(
'Kd_turb',
'Kd_shear', restart_cs)
2271 call register_restart_field_as_obsolete(
'Kv_turb',
'Kv_shear', restart_cs)
2276 if ((us%m_to_Z_restart /= 0.0) .and. (us%m_to_Z_restart /= us%m_to_Z)) &
2277 z_rescale = us%m_to_Z / us%m_to_Z_restart
2279 if ((us%s_to_T_restart /= 0.0) .and. (us%s_to_T_restart /= us%s_to_T)) &
2280 i_t_rescale = us%s_to_T_restart / us%s_to_T
2281 z2_t_rescale = z_rescale**2*i_t_rescale
2283 if (z2_t_rescale /= 1.0)
then
2284 if (
associated(visc%Kd_shear))
then ;
if (
query_initialized(visc%Kd_shear,
"Kd_shear", restart_cs))
then
2285 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2286 visc%Kd_shear(i,j,k) = z2_t_rescale * visc%Kd_shear(i,j,k)
2287 enddo ;
enddo ;
enddo
2290 if (
associated(visc%Kv_shear))
then ;
if (
query_initialized(visc%Kv_shear,
"Kv_shear", restart_cs))
then
2291 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2292 visc%Kv_shear(i,j,k) = z2_t_rescale * visc%Kv_shear(i,j,k)
2293 enddo ;
enddo ;
enddo
2296 if (
associated(visc%Kv_shear_Bu))
then ;
if (
query_initialized(visc%Kv_shear_Bu,
"Kv_shear_Bu", restart_cs))
then
2297 do k=1,nz+1 ;
do j=js-1,je ;
do i=is-1,ie
2298 visc%Kv_shear_Bu(i,j,k) = z2_t_rescale * visc%Kv_shear_Bu(i,j,k)
2299 enddo ;
enddo ;
enddo
2303 if (mle_use_pbl_mld .and. (z_rescale /= 1.0))
then
2304 if (
associated(visc%MLD))
then ;
if (
query_initialized(visc%MLD,
"MLD", restart_cs))
then
2305 do j=js,je ;
do i=is,ie
2306 visc%MLD(i,j) = z_rescale * visc%MLD(i,j)
2311 end subroutine set_visc_init
2314 subroutine set_visc_end(visc, CS)
2319 if (cs%bottomdraglaw)
then
2320 deallocate(visc%bbl_thick_u) ;
deallocate(visc%bbl_thick_v)
2321 deallocate(visc%kv_bbl_u) ;
deallocate(visc%kv_bbl_v)
2322 if (
allocated(cs%bbl_u))
deallocate(cs%bbl_u)
2323 if (
allocated(cs%bbl_v))
deallocate(cs%bbl_v)
2325 if (cs%Channel_drag)
then
2326 deallocate(visc%Ray_u) ;
deallocate(visc%Ray_v)
2328 if (cs%dynamic_viscous_ML)
then
2329 deallocate(visc%nkml_visc_u) ;
deallocate(visc%nkml_visc_v)
2331 if (
associated(visc%Kd_shear))
deallocate(visc%Kd_shear)
2332 if (
associated(visc%Kv_slow))
deallocate(visc%Kv_slow)
2333 if (
associated(visc%TKE_turb))
deallocate(visc%TKE_turb)
2334 if (
associated(visc%Kv_shear))
deallocate(visc%Kv_shear)
2335 if (
associated(visc%Kv_shear_Bu))
deallocate(visc%Kv_shear_Bu)
2336 if (
associated(visc%ustar_bbl))
deallocate(visc%ustar_bbl)
2337 if (
associated(visc%TKE_bbl))
deallocate(visc%TKE_bbl)
2338 if (
associated(visc%taux_shelf))
deallocate(visc%taux_shelf)
2339 if (
associated(visc%tauy_shelf))
deallocate(visc%tauy_shelf)
2340 if (
associated(visc%tbl_thick_shelf_u))
deallocate(visc%tbl_thick_shelf_u)
2341 if (
associated(visc%tbl_thick_shelf_v))
deallocate(visc%tbl_thick_shelf_v)
2342 if (
associated(visc%kv_tbl_shelf_u))
deallocate(visc%kv_tbl_shelf_u)
2343 if (
associated(visc%kv_tbl_shelf_v))
deallocate(visc%kv_tbl_shelf_v)
2346 end subroutine set_visc_end