16 use mom_pointaccel,
only : write_u_accel, write_v_accel, pointaccel_init
25 implicit none ;
private 27 #include <MOM_memory.h> 29 public vertvisc, vertvisc_remnant, vertvisc_coef
30 public vertvisc_limit_vel, vertvisc_init, vertvisc_end
31 public updatecfltruncationvalue
52 logical :: cfl_based_trunc
64 logical :: cflrampingisactivated = .false.
65 type(time_type) :: rampstarttime
67 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_,NK_INTERFACE_) :: &
69 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_,NKMEM_) :: &
71 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_,NK_INTERFACE_) :: &
73 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_,NKMEM_) :: &
75 real,
pointer,
dimension(:,:) :: a1_shelf_u => null()
77 real,
pointer,
dimension(:,:) :: a1_shelf_v => null()
81 logical :: bottomdraglaw
86 logical :: channel_drag
89 logical :: harmonic_visc
97 logical :: direct_stress
99 logical :: dynamic_viscous_ml
103 logical :: answers_2018
109 integer,
pointer :: ntrunc
111 character(len=200) :: u_trunc_file
113 character(len=200) :: v_trunc_file
115 logical :: stokesmixing
123 integer :: id_du_dt_visc = -1, id_dv_dt_visc = -1, id_au_vv = -1, id_av_vv = -1
124 integer :: id_h_u = -1, id_h_v = -1, id_hml_u = -1 , id_hml_v = -1
125 integer :: id_taux_bot = -1, id_tauy_bot = -1
126 integer :: id_kv_slow = -1, id_kv_u = -1, id_kv_v = -1
128 integer :: id_hf_du_dt_visc_2d = -1, id_hf_dv_dt_visc_2d = -1
157 subroutine vertvisc(u, v, h, forces, visc, dt, OBC, ADp, CDp, G, GV, US, CS, &
158 taux_bot, tauy_bot, Waves)
162 real,
dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
164 real,
dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
166 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
170 real,
intent(in) :: dt
176 real,
dimension(SZIB_(G),SZJ_(G)), &
177 optional,
intent(out) :: taux_bot
179 real,
dimension(SZI_(G),SZJB_(G)), &
180 optional,
intent(out) :: tauy_bot
183 optional,
pointer :: waves
192 real :: c1(szib_(g),szk_(g))
194 real :: ray(szib_(g),szk_(g))
209 real :: zds, hfr, h_a
210 real :: surface_stress(szib_(g))
213 real,
allocatable,
dimension(:,:) :: hf_du_dt_visc_2d
214 real,
allocatable,
dimension(:,:) :: hf_dv_dt_visc_2d
216 logical :: do_i(szib_(g))
217 logical :: dostokesmixing
219 integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz, n
220 is = g%isc ; ie = g%iec; js = g%jsc; je = g%jec
221 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB ; nz = g%ke
223 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_vert_friction(visc): "// &
224 "Module must be initialized before it is used.")
226 if (cs%direct_stress)
then 227 hmix = cs%Hmix_stress
230 dt_rho0 = dt / gv%H_to_RZ
231 dt_z_to_h = dt*gv%Z_to_H
232 h_neglect = gv%H_subroundoff
236 dostokesmixing=.false.
237 if (cs%StokesMixing)
then 238 if (
present(waves)) dostokesmixing =
associated(waves)
239 if (.not. dostokesmixing) &
240 call mom_error(fatal,
"Stokes Mixing called without allocated"//&
241 "Waves Control Structure")
244 do k=1,nz ;
do i=isq,ieq ; ray(i,k) = 0.0 ;
enddo ;
enddo 253 do i=isq,ieq ; do_i(i) = (g%mask2dCu(i,j) > 0) ;
enddo 256 if (dostokesmixing)
then ;
do k=1,nz ;
do i=isq,ieq
257 if (do_i(i)) u(i,j,k) = u(i,j,k) + us%m_s_to_L_T*waves%Us_x(i,j,k)
258 enddo ;
enddo ;
endif 260 if (
associated(adp%du_dt_visc))
then ;
do k=1,nz ;
do i=isq,ieq
261 adp%du_dt_visc(i,j,k) = u(i,j,k)
262 enddo ;
enddo ;
endif 267 if (cs%direct_stress)
then 268 do i=isq,ieq ;
if (do_i(i))
then 269 surface_stress(i) = 0.0
271 stress = dt_rho0 * forces%taux(i,j)
273 h_a = 0.5 * (h(i,j,k) + h(i+1,j,k)) + h_neglect
274 hfr = 1.0 ;
if ((zds+h_a) > hmix) hfr = (hmix - zds) / h_a
275 u(i,j,k) = u(i,j,k) + i_hmix * hfr * stress
276 zds = zds + h_a ;
if (zds >= hmix)
exit 280 surface_stress(i) = dt_rho0 * (g%mask2dCu(i,j)*forces%taux(i,j))
283 if (cs%Channel_drag)
then ;
do k=1,nz ;
do i=isq,ieq
284 ray(i,k) = visc%Ray_u(i,j,k)
285 enddo ;
enddo ;
endif 310 do i=isq,ieq ;
if (do_i(i))
then 311 b_denom_1 = cs%h_u(i,j,1) + dt_z_to_h * (ray(i,1) + cs%a_u(i,j,1))
312 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h*cs%a_u(i,j,2))
313 d1(i) = b_denom_1 * b1(i)
314 u(i,j,1) = b1(i) * (cs%h_u(i,j,1) * u(i,j,1) + surface_stress(i))
316 do k=2,nz ;
do i=isq,ieq ;
if (do_i(i))
then 317 c1(i,k) = dt_z_to_h * cs%a_u(i,j,k) * b1(i)
318 b_denom_1 = cs%h_u(i,j,k) + dt_z_to_h * (ray(i,k) + cs%a_u(i,j,k)*d1(i))
319 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h * cs%a_u(i,j,k+1))
320 d1(i) = b_denom_1 * b1(i)
321 u(i,j,k) = (cs%h_u(i,j,k) * u(i,j,k) + &
322 dt_z_to_h * cs%a_u(i,j,k) * u(i,j,k-1)) * b1(i)
323 endif ;
enddo ;
enddo 327 do k=nz-1,1,-1 ;
do i=isq,ieq ;
if (do_i(i))
then 328 u(i,j,k) = u(i,j,k) + c1(i,k+1) * u(i,j,k+1)
329 endif ;
enddo ;
enddo 331 if (
associated(adp%du_dt_visc))
then ;
do k=1,nz ;
do i=isq,ieq
332 adp%du_dt_visc(i,j,k) = (u(i,j,k) - adp%du_dt_visc(i,j,k))*idt
333 enddo ;
enddo ;
endif 335 if (
associated(visc%taux_shelf))
then ;
do i=isq,ieq
336 visc%taux_shelf(i,j) = -gv%Rho0*cs%a1_shelf_u(i,j)*u(i,j,1)
339 if (
PRESENT(taux_bot))
then 341 taux_bot(i,j) = gv%Rho0 * (u(i,j,nz)*cs%a_u(i,j,nz+1))
343 if (cs%Channel_drag)
then ;
do k=1,nz ;
do i=isq,ieq
344 taux_bot(i,j) = taux_bot(i,j) + gv%Rho0 * (ray(i,k)*u(i,j,k))
345 enddo ;
enddo ;
endif 349 if (dostokesmixing)
then ;
do k=1,nz ;
do i=isq,ieq
350 if (do_i(i)) u(i,j,k) = u(i,j,k) - us%m_s_to_L_T*waves%Us_x(i,j,k)
351 enddo ;
enddo ;
endif 361 do i=is,ie ; do_i(i) = (g%mask2dCv(i,j) > 0) ;
enddo 364 if (dostokesmixing)
then ;
do k=1,nz ;
do i=is,ie
365 if (do_i(i)) v(i,j,k) = v(i,j,k) + us%m_s_to_L_T*waves%Us_y(i,j,k)
366 enddo ;
enddo ;
endif 368 if (
associated(adp%dv_dt_visc))
then ;
do k=1,nz ;
do i=is,ie
369 adp%dv_dt_visc(i,j,k) = v(i,j,k)
370 enddo ;
enddo ;
endif 375 if (cs%direct_stress)
then 376 do i=is,ie ;
if (do_i(i))
then 377 surface_stress(i) = 0.0
379 stress = dt_rho0 * forces%tauy(i,j)
381 h_a = 0.5 * (h(i,j,k) + h(i,j+1,k)) + h_neglect
382 hfr = 1.0 ;
if ((zds+h_a) > hmix) hfr = (hmix - zds) / h_a
383 v(i,j,k) = v(i,j,k) + i_hmix * hfr * stress
384 zds = zds + h_a ;
if (zds >= hmix)
exit 388 surface_stress(i) = dt_rho0 * (g%mask2dCv(i,j)*forces%tauy(i,j))
391 if (cs%Channel_drag)
then ;
do k=1,nz ;
do i=is,ie
392 ray(i,k) = visc%Ray_v(i,j,k)
393 enddo ;
enddo ;
endif 395 do i=is,ie ;
if (do_i(i))
then 396 b_denom_1 = cs%h_v(i,j,1) + dt_z_to_h * (ray(i,1) + cs%a_v(i,j,1))
397 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h*cs%a_v(i,j,2))
398 d1(i) = b_denom_1 * b1(i)
399 v(i,j,1) = b1(i) * (cs%h_v(i,j,1) * v(i,j,1) + surface_stress(i))
401 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 402 c1(i,k) = dt_z_to_h * cs%a_v(i,j,k) * b1(i)
403 b_denom_1 = cs%h_v(i,j,k) + dt_z_to_h * (ray(i,k) + cs%a_v(i,j,k)*d1(i))
404 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h * cs%a_v(i,j,k+1))
405 d1(i) = b_denom_1 * b1(i)
406 v(i,j,k) = (cs%h_v(i,j,k) * v(i,j,k) + dt_z_to_h * cs%a_v(i,j,k) * v(i,j,k-1)) * b1(i)
407 endif ;
enddo ;
enddo 408 do k=nz-1,1,-1 ;
do i=is,ie ;
if (do_i(i))
then 409 v(i,j,k) = v(i,j,k) + c1(i,k+1) * v(i,j,k+1)
410 endif ;
enddo ;
enddo 412 if (
associated(adp%dv_dt_visc))
then ;
do k=1,nz ;
do i=is,ie
413 adp%dv_dt_visc(i,j,k) = (v(i,j,k) - adp%dv_dt_visc(i,j,k))*idt
414 enddo ;
enddo ;
endif 416 if (
associated(visc%tauy_shelf))
then ;
do i=is,ie
417 visc%tauy_shelf(i,j) = -gv%Rho0*cs%a1_shelf_v(i,j)*v(i,j,1)
420 if (
present(tauy_bot))
then 422 tauy_bot(i,j) = gv%Rho0 * (v(i,j,nz)*cs%a_v(i,j,nz+1))
424 if (cs%Channel_drag)
then ;
do k=1,nz ;
do i=is,ie
425 tauy_bot(i,j) = tauy_bot(i,j) + gv%Rho0 * (ray(i,k)*v(i,j,k))
426 enddo ;
enddo ;
endif 430 if (dostokesmixing)
then ;
do k=1,nz ;
do i=is,ie
431 if (do_i(i)) v(i,j,k) = v(i,j,k) - us%m_s_to_L_T*waves%Us_y(i,j,k)
432 enddo ;
enddo ;
endif 436 call vertvisc_limit_vel(u, v, h, adp, cdp, forces, visc, dt, g, gv, us, cs)
439 if (
associated(obc))
then 440 do n=1,obc%number_of_segments
441 if (obc%segment(n)%specified)
then 442 if (obc%segment(n)%is_N_or_S)
then 443 j = obc%segment(n)%HI%JsdB
444 do k=1,nz ;
do i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
445 v(i,j,k) = obc%segment(n)%normal_vel(i,j,k)
447 elseif (obc%segment(n)%is_E_or_W)
then 448 i = obc%segment(n)%HI%IsdB
449 do k=1,nz ;
do j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
450 u(i,j,k) = obc%segment(n)%normal_vel(i,j,k)
458 if (cs%id_du_dt_visc > 0) &
459 call post_data(cs%id_du_dt_visc, adp%du_dt_visc, cs%diag)
460 if (cs%id_dv_dt_visc > 0) &
461 call post_data(cs%id_dv_dt_visc, adp%dv_dt_visc, cs%diag)
462 if (
present(taux_bot) .and. (cs%id_taux_bot > 0)) &
463 call post_data(cs%id_taux_bot, taux_bot, cs%diag)
464 if (
present(tauy_bot) .and. (cs%id_tauy_bot > 0)) &
465 call post_data(cs%id_tauy_bot, tauy_bot, cs%diag)
483 if (cs%id_hf_du_dt_visc_2d > 0)
then 484 allocate(hf_du_dt_visc_2d(g%IsdB:g%IedB,g%jsd:g%jed))
485 hf_du_dt_visc_2d(:,:) = 0.0
486 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
487 hf_du_dt_visc_2d(i,j) = hf_du_dt_visc_2d(i,j) + adp%du_dt_visc(i,j,k) * adp%diag_hfrac_u(i,j,k)
488 enddo ;
enddo ;
enddo 489 call post_data(cs%id_hf_du_dt_visc_2d, hf_du_dt_visc_2d, cs%diag)
490 deallocate(hf_du_dt_visc_2d)
492 if (cs%id_hf_dv_dt_visc_2d > 0)
then 493 allocate(hf_dv_dt_visc_2d(g%isd:g%ied,g%JsdB:g%JedB))
494 hf_dv_dt_visc_2d(:,:) = 0.0
495 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
496 hf_dv_dt_visc_2d(i,j) = hf_dv_dt_visc_2d(i,j) + adp%dv_dt_visc(i,j,k) * adp%diag_hfrac_v(i,j,k)
497 enddo ;
enddo ;
enddo 498 call post_data(cs%id_hf_dv_dt_visc_2d, hf_dv_dt_visc_2d, cs%diag)
499 deallocate(hf_dv_dt_visc_2d)
502 end subroutine vertvisc
508 subroutine vertvisc_remnant(visc, visc_rem_u, visc_rem_v, dt, G, GV, US, CS)
512 real,
dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
513 intent(inout) :: visc_rem_u
516 real,
dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
517 intent(inout) :: visc_rem_v
520 real,
intent(in) :: dt
527 real :: c1(szib_(g),szk_(g))
529 real :: ray(szib_(g),szk_(g))
533 logical :: do_i(szib_(g))
535 integer :: i, j, k, is, ie, isq, ieq, jsq, jeq, nz
536 is = g%isc ; ie = g%iec
537 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB ; nz = g%ke
539 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_vert_friction(visc): "// &
540 "Module must be initialized before it is used.")
542 dt_z_to_h = dt*gv%Z_to_H
544 do k=1,nz ;
do i=isq,ieq ; ray(i,k) = 0.0 ;
enddo ;
enddo 551 do i=isq,ieq ; do_i(i) = (g%mask2dCu(i,j) > 0) ;
enddo 553 if (cs%Channel_drag)
then ;
do k=1,nz ;
do i=isq,ieq
554 ray(i,k) = visc%Ray_u(i,j,k)
555 enddo ;
enddo ;
endif 557 do i=isq,ieq ;
if (do_i(i))
then 558 b_denom_1 = cs%h_u(i,j,1) + dt_z_to_h * (ray(i,1) + cs%a_u(i,j,1))
559 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h*cs%a_u(i,j,2))
560 d1(i) = b_denom_1 * b1(i)
561 visc_rem_u(i,j,1) = b1(i) * cs%h_u(i,j,1)
563 do k=2,nz ;
do i=isq,ieq ;
if (do_i(i))
then 564 c1(i,k) = dt_z_to_h * cs%a_u(i,j,k)*b1(i)
565 b_denom_1 = cs%h_u(i,j,k) + dt_z_to_h * (ray(i,k) + cs%a_u(i,j,k)*d1(i))
566 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h * cs%a_u(i,j,k+1))
567 d1(i) = b_denom_1 * b1(i)
568 visc_rem_u(i,j,k) = (cs%h_u(i,j,k) + dt_z_to_h * cs%a_u(i,j,k) * visc_rem_u(i,j,k-1)) * b1(i)
569 endif ;
enddo ;
enddo 570 do k=nz-1,1,-1 ;
do i=isq,ieq ;
if (do_i(i))
then 571 visc_rem_u(i,j,k) = visc_rem_u(i,j,k) + c1(i,k+1)*visc_rem_u(i,j,k+1)
573 endif ;
enddo ;
enddo 582 do i=is,ie ; do_i(i) = (g%mask2dCv(i,j) > 0) ;
enddo 584 if (cs%Channel_drag)
then ;
do k=1,nz ;
do i=is,ie
585 ray(i,k) = visc%Ray_v(i,j,k)
586 enddo ;
enddo ;
endif 588 do i=is,ie ;
if (do_i(i))
then 589 b_denom_1 = cs%h_v(i,j,1) + dt_z_to_h * (ray(i,1) + cs%a_v(i,j,1))
590 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h*cs%a_v(i,j,2))
591 d1(i) = b_denom_1 * b1(i)
592 visc_rem_v(i,j,1) = b1(i) * cs%h_v(i,j,1)
594 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 595 c1(i,k) = dt_z_to_h * cs%a_v(i,j,k)*b1(i)
596 b_denom_1 = cs%h_v(i,j,k) + dt_z_to_h * (ray(i,k) + cs%a_v(i,j,k)*d1(i))
597 b1(i) = 1.0 / (b_denom_1 + dt_z_to_h * cs%a_v(i,j,k+1))
598 d1(i) = b_denom_1 * b1(i)
599 visc_rem_v(i,j,k) = (cs%h_v(i,j,k) + dt_z_to_h * cs%a_v(i,j,k) * visc_rem_v(i,j,k-1)) * b1(i)
600 endif ;
enddo ;
enddo 601 do k=nz-1,1,-1 ;
do i=is,ie ;
if (do_i(i))
then 602 visc_rem_v(i,j,k) = visc_rem_v(i,j,k) + c1(i,k+1)*visc_rem_v(i,j,k+1)
603 endif ;
enddo ;
enddo 607 call uvchksum(
"visc_rem_[uv]", visc_rem_u, visc_rem_v, g%HI, haloshift=0, &
611 end subroutine vertvisc_remnant
617 subroutine vertvisc_coef(u, v, h, forces, visc, dt, G, GV, US, CS, OBC)
621 real,
dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
623 real,
dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
625 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
629 real,
intent(in) :: dt
639 real,
dimension(SZIB_(G),SZK_(G)) :: &
646 real,
dimension(SZIB_(G),SZK_(G)+1) :: &
653 real,
dimension(SZIB_(G)) :: &
666 real,
allocatable,
dimension(:,:) :: hml_u
667 real,
allocatable,
dimension(:,:) :: hml_v
668 real,
allocatable,
dimension(:,:,:) :: kv_u
669 real,
allocatable,
dimension(:,:,:) :: kv_v
670 real :: zcol(szi_(g))
686 logical,
dimension(SZIB_(G)) :: do_i, do_i_shelf
687 logical :: do_any_shelf
688 integer,
dimension(SZIB_(G)) :: &
691 integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz
692 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
693 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB ; nz = g%ke
695 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_vert_friction(coef): "// &
696 "Module must be initialized before it is used.")
698 h_neglect = gv%H_subroundoff
699 a_cpl_max = 1.0e37 * us%m_to_Z * us%T_to_s
700 i_hbbl(:) = 1.0 / (cs%Hbbl + h_neglect)
701 i_valbl = 0.0 ;
if (cs%harm_BL_val > 0.0) i_valbl = 1.0 / cs%harm_BL_val
703 if (cs%id_Kv_u > 0)
then 704 allocate(kv_u(g%IsdB:g%IedB,g%jsd:g%jed,g%ke)) ; kv_u(:,:,:) = 0.0
707 if (cs%id_Kv_v > 0)
then 708 allocate(kv_v(g%isd:g%ied,g%JsdB:g%JedB,g%ke)) ; kv_v(:,:,:) = 0.0
711 if (cs%debug .or. (cs%id_hML_u > 0))
then 712 allocate(hml_u(g%IsdB:g%IedB,g%jsd:g%jed)) ; hml_u(:,:) = 0.0
714 if (cs%debug .or. (cs%id_hML_v > 0))
then 715 allocate(hml_v(g%isd:g%ied,g%JsdB:g%JedB)) ; hml_v(:,:) = 0.0
718 if ((
associated(visc%taux_shelf) .or.
associated(forces%frac_shelf_u)) .and. &
719 .not.
associated(cs%a1_shelf_u))
then 720 allocate(cs%a1_shelf_u(g%IsdB:g%IedB,g%jsd:g%jed)) ; cs%a1_shelf_u(:,:)=0.0
722 if ((
associated(visc%tauy_shelf) .or.
associated(forces%frac_shelf_v)) .and. &
723 .not.
associated(cs%a1_shelf_v))
then 724 allocate(cs%a1_shelf_v(g%isd:g%ied,g%JsdB:g%JedB)) ; cs%a1_shelf_v(:,:)=0.0
731 do i=isq,ieq ; do_i(i) = (g%mask2dCu(i,j) > 0) ;
enddo 733 if (cs%bottomdraglaw)
then ;
do i=isq,ieq
734 kv_bbl(i) = visc%Kv_bbl_u(i,j)
735 bbl_thick(i) = visc%bbl_thick_u(i,j) * gv%Z_to_H + h_neglect
736 if (do_i(i)) i_hbbl(i) = 1.0 / bbl_thick(i)
739 do k=1,nz ;
do i=isq,ieq ;
if (do_i(i))
then 740 h_harm(i,k) = 2.0*h(i,j,k)*h(i+1,j,k) / (h(i,j,k)+h(i+1,j,k)+h_neglect)
741 h_arith(i,k) = 0.5*(h(i+1,j,k)+h(i,j,k))
742 h_delta(i,k) = h(i+1,j,k) - h(i,j,k)
743 endif ;
enddo ;
enddo 745 dmin(i) = min(g%bathyT(i,j), g%bathyT(i+1,j)) * gv%Z_to_H
750 if (
associated(obc))
then ;
if (obc%number_of_segments > 0)
then 751 do i=isq,ieq ;
if (do_i(i) .and. (obc%segnum_u(i,j) /= obc_none))
then 752 if (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_e)
then 753 do k=1,nz ; h_harm(i,k) = h(i,j,k) ; h_arith(i,k) = h(i,j,k) ; h_delta(i,k) = 0. ;
enddo 754 dmin(i) = g%bathyT(i,j) * gv%Z_to_H
756 elseif (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_w)
then 757 do k=1,nz ; h_harm(i,k) = h(i+1,j,k) ; h_arith(i,k) = h(i+1,j,k) ; h_delta(i,k) = 0. ;
enddo 758 dmin(i) = g%bathyT(i+1,j) * gv%Z_to_H
769 if (cs%harmonic_visc)
then 770 do i=isq,ieq ; z_i(i,nz+1) = 0.0 ;
enddo 771 do k=nz,1,-1 ;
do i=isq,ieq ;
if (do_i(i))
then 772 hvel(i,k) = h_harm(i,k)
773 if (u(i,j,k) * h_delta(i,k) < 0)
then 774 z2 = z_i(i,k+1) ; botfn = 1.0 / (1.0 + 0.09*z2*z2*z2*z2*z2*z2)
775 hvel(i,k) = (1.0-botfn)*h_harm(i,k) + botfn*h_arith(i,k)
777 z_i(i,k) = z_i(i,k+1) + h_harm(i,k)*i_hbbl(i)
778 endif ;
enddo ;
enddo 780 do i=isq,ieq ; zh(i) = 0.0 ; z_i(i,nz+1) = 0.0 ;
enddo 781 do i=isq,ieq+1 ; zcol(i) = -g%bathyT(i,j) * gv%Z_to_H ;
enddo 783 do i=isq,ieq+1 ; zcol(i) = zcol(i) + h(i,j,k) ;
enddo 784 do i=isq,ieq ;
if (do_i(i))
then 785 zh(i) = zh(i) + h_harm(i,k)
787 z_clear = max(zcol(i),zcol(i+1)) + dmin(i)
788 if (zi_dir(i) < 0) z_clear = zcol(i) + dmin(i)
789 if (zi_dir(i) > 0) z_clear = zcol(i+1) + dmin(i)
791 z_i(i,k) = max(zh(i), z_clear) * i_hbbl(i)
793 hvel(i,k) = h_arith(i,k)
794 if (u(i,j,k) * h_delta(i,k) > 0)
then 795 if (zh(i) * i_hbbl(i) < cs%harm_BL_val)
then 796 hvel(i,k) = h_harm(i,k)
798 z2_wt = 1.0 ;
if (zh(i) * i_hbbl(i) < 2.0*cs%harm_BL_val) &
799 z2_wt = max(0.0, min(1.0, zh(i) * i_hbbl(i) * i_valbl - 1.0))
800 z2 = z2_wt * (max(zh(i), z_clear) * i_hbbl(i))
801 botfn = 1.0 / (1.0 + 0.09*z2*z2*z2*z2*z2*z2)
802 hvel(i,k) = (1.0-botfn)*h_arith(i,k) + botfn*h_harm(i,k)
810 call find_coupling_coef(a_cpl, hvel, do_i, h_harm, bbl_thick, kv_bbl, z_i, h_ml, &
811 dt, j, g, gv, us, cs, visc, forces, work_on_u=.true., obc=obc)
812 if (
allocated(hml_u))
then 813 do i=isq,ieq ;
if (do_i(i))
then ; hml_u(i,j) = h_ml(i) ;
endif ;
enddo 816 do_any_shelf = .false.
817 if (
associated(forces%frac_shelf_u))
then 819 cs%a1_shelf_u(i,j) = 0.0
820 do_i_shelf(i) = (do_i(i) .and. forces%frac_shelf_u(i,j) > 0.0)
821 if (do_i_shelf(i)) do_any_shelf = .true.
823 if (do_any_shelf)
then 824 if (cs%harmonic_visc)
then 825 call find_coupling_coef(a_shelf, hvel, do_i_shelf, h_harm, bbl_thick, &
826 kv_bbl, z_i, h_ml, dt, j, g, gv, us, cs, &
827 visc, forces, work_on_u=.true., obc=obc, shelf=.true.)
830 do i=isq,ieq ;
if (do_i_shelf(i))
then 831 zh(i) = 0.0 ; ztop_min(i) = min(zcol(i), zcol(i+1))
832 i_htbl(i) = 1.0 / (visc%tbl_thick_shelf_u(i,j)*gv%Z_to_H + h_neglect)
835 do i=isq,ieq+1 ; zcol(i) = zcol(i) - h(i,j,k) ;
enddo 836 do i=isq,ieq ;
if (do_i_shelf(i))
then 837 zh(i) = zh(i) + h_harm(i,k)
839 hvel_shelf(i,k) = hvel(i,k)
840 if (u(i,j,k) * h_delta(i,k) > 0)
then 841 if (zh(i) * i_htbl(i) < cs%harm_BL_val)
then 842 hvel_shelf(i,k) = min(hvel(i,k), h_harm(i,k))
844 z2_wt = 1.0 ;
if (zh(i) * i_htbl(i) < 2.0*cs%harm_BL_val) &
845 z2_wt = max(0.0, min(1.0, zh(i) * i_htbl(i) * i_valbl - 1.0))
846 z2 = z2_wt * (max(zh(i), ztop_min(i) - min(zcol(i),zcol(i+1))) * i_htbl(i))
847 topfn = 1.0 / (1.0 + 0.09*z2**6)
848 hvel_shelf(i,k) = min(hvel(i,k), (1.0-topfn)*h_arith(i,k) + topfn*h_harm(i,k))
853 call find_coupling_coef(a_shelf, hvel_shelf, do_i_shelf, h_harm, &
854 bbl_thick, kv_bbl, z_i, h_ml, dt, j, g, gv, us, cs, &
855 visc, forces, work_on_u=.true., obc=obc, shelf=.true.)
857 do i=isq,ieq ;
if (do_i_shelf(i)) cs%a1_shelf_u(i,j) = a_shelf(i,1) ;
enddo 861 if (do_any_shelf)
then 862 do k=1,nz+1 ;
do i=isq,ieq ;
if (do_i_shelf(i))
then 863 cs%a_u(i,j,k) = min(a_cpl_max, forces%frac_shelf_u(i,j) * a_shelf(i,k) + &
864 (1.0-forces%frac_shelf_u(i,j)) * a_cpl(i,k))
868 elseif (do_i(i))
then 869 cs%a_u(i,j,k) = min(a_cpl_max, a_cpl(i,k))
870 endif ;
enddo ;
enddo 871 do k=1,nz ;
do i=isq,ieq ;
if (do_i_shelf(i))
then 873 cs%h_u(i,j,k) = forces%frac_shelf_u(i,j) * hvel_shelf(i,k) + &
874 (1.0-forces%frac_shelf_u(i,j)) * hvel(i,k) + h_neglect
875 elseif (do_i(i))
then 876 cs%h_u(i,j,k) = hvel(i,k) + h_neglect
877 endif ;
enddo ;
enddo 879 do k=1,nz+1 ;
do i=isq,ieq ;
if (do_i(i)) cs%a_u(i,j,k) = min(a_cpl_max, a_cpl(i,k)) ;
enddo ;
enddo 880 do k=1,nz ;
do i=isq,ieq ;
if (do_i(i)) cs%h_u(i,j,k) = hvel(i,k) + h_neglect ;
enddo ;
enddo 884 if (cs%id_Kv_u > 0)
then 885 do k=1,nz ;
do i=isq,ieq
886 if (do_i(i)) kv_u(i,j,k) = 0.5 * gv%H_to_Z*(cs%a_u(i,j,k)+cs%a_u(i,j,k+1)) * cs%h_u(i,j,k)
898 do i=is,ie ; do_i(i) = (g%mask2dCv(i,j) > 0) ;
enddo 900 if (cs%bottomdraglaw)
then ;
do i=is,ie
901 kv_bbl(i) = visc%Kv_bbl_v(i,j)
902 bbl_thick(i) = visc%bbl_thick_v(i,j) * gv%Z_to_H + h_neglect
903 if (do_i(i)) i_hbbl(i) = 1.0 / bbl_thick(i)
906 do k=1,nz ;
do i=is,ie ;
if (do_i(i))
then 907 h_harm(i,k) = 2.0*h(i,j,k)*h(i,j+1,k) / (h(i,j,k)+h(i,j+1,k)+h_neglect)
908 h_arith(i,k) = 0.5*(h(i,j+1,k)+h(i,j,k))
909 h_delta(i,k) = h(i,j+1,k) - h(i,j,k)
910 endif ;
enddo ;
enddo 912 dmin(i) = min(g%bathyT(i,j), g%bathyT(i,j+1)) * gv%Z_to_H
917 if (
associated(obc))
then ;
if (obc%number_of_segments > 0)
then 918 do i=is,ie ;
if (do_i(i) .and. (obc%segnum_v(i,j) /= obc_none))
then 919 if (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_n)
then 920 do k=1,nz ; h_harm(i,k) = h(i,j,k) ; h_arith(i,k) = h(i,j,k) ; h_delta(i,k) = 0. ;
enddo 921 dmin(i) = g%bathyT(i,j) * gv%Z_to_H
923 elseif (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_s)
then 924 do k=1,nz ; h_harm(i,k) = h(i,j+1,k) ; h_arith(i,k) = h(i,j+1,k) ; h_delta(i,k) = 0. ;
enddo 925 dmin(i) = g%bathyT(i,j+1) * gv%Z_to_H
936 if (cs%harmonic_visc)
then 937 do i=is,ie ; z_i(i,nz+1) = 0.0 ;
enddo 939 do k=nz,1,-1 ;
do i=is,ie ;
if (do_i(i))
then 940 hvel(i,k) = h_harm(i,k)
941 if (v(i,j,k) * h_delta(i,k) < 0)
then 942 z2 = z_i(i,k+1) ; botfn = 1.0 / (1.0 + 0.09*z2*z2*z2*z2*z2*z2)
943 hvel(i,k) = (1.0-botfn)*h_harm(i,k) + botfn*h_arith(i,k)
945 z_i(i,k) = z_i(i,k+1) + h_harm(i,k)*i_hbbl(i)
946 endif ;
enddo ;
enddo 949 zh(i) = 0.0 ; z_i(i,nz+1) = 0.0
950 zcol1(i) = -g%bathyT(i,j) * gv%Z_to_H
951 zcol2(i) = -g%bathyT(i,j+1) * gv%Z_to_H
953 do k=nz,1,-1 ;
do i=is,ie ;
if (do_i(i))
then 954 zh(i) = zh(i) + h_harm(i,k)
955 zcol1(i) = zcol1(i) + h(i,j,k) ; zcol2(i) = zcol2(i) + h(i,j+1,k)
957 z_clear = max(zcol1(i),zcol2(i)) + dmin(i)
958 if (zi_dir(i) < 0) z_clear = zcol1(i) + dmin(i)
959 if (zi_dir(i) > 0) z_clear = zcol2(i) + dmin(i)
961 z_i(i,k) = max(zh(i), z_clear) * i_hbbl(i)
963 hvel(i,k) = h_arith(i,k)
964 if (v(i,j,k) * h_delta(i,k) > 0)
then 965 if (zh(i) * i_hbbl(i) < cs%harm_BL_val)
then 966 hvel(i,k) = h_harm(i,k)
968 z2_wt = 1.0 ;
if (zh(i) * i_hbbl(i) < 2.0*cs%harm_BL_val) &
969 z2_wt = max(0.0, min(1.0, zh(i) * i_hbbl(i) * i_valbl - 1.0))
970 z2 = z2_wt * (max(zh(i), max(zcol1(i),zcol2(i)) + dmin(i)) * i_hbbl(i))
971 botfn = 1.0 / (1.0 + 0.09*z2*z2*z2*z2*z2*z2)
972 hvel(i,k) = (1.0-botfn)*h_arith(i,k) + botfn*h_harm(i,k)
976 endif ;
enddo ;
enddo 979 call find_coupling_coef(a_cpl, hvel, do_i, h_harm, bbl_thick, kv_bbl, z_i, h_ml, &
980 dt, j, g, gv, us, cs, visc, forces, work_on_u=.false., obc=obc)
981 if (
allocated(hml_v))
then 982 do i=is,ie ;
if (do_i(i))
then ; hml_v(i,j) = h_ml(i) ;
endif ;
enddo 984 do_any_shelf = .false.
985 if (
associated(forces%frac_shelf_v))
then 987 cs%a1_shelf_v(i,j) = 0.0
988 do_i_shelf(i) = (do_i(i) .and. forces%frac_shelf_v(i,j) > 0.0)
989 if (do_i_shelf(i)) do_any_shelf = .true.
991 if (do_any_shelf)
then 992 if (cs%harmonic_visc)
then 993 call find_coupling_coef(a_shelf, hvel, do_i_shelf, h_harm, bbl_thick, &
994 kv_bbl, z_i, h_ml, dt, j, g, gv, us, cs, visc, &
995 forces, work_on_u=.false., obc=obc, shelf=.true.)
998 do i=is,ie ;
if (do_i_shelf(i))
then 999 zh(i) = 0.0 ; ztop_min(i) = min(zcol1(i), zcol2(i))
1000 i_htbl(i) = 1.0 / (visc%tbl_thick_shelf_v(i,j)*gv%Z_to_H + h_neglect)
1003 do i=is,ie ;
if (do_i_shelf(i))
then 1004 zcol1(i) = zcol1(i) - h(i,j,k) ; zcol2(i) = zcol2(i) - h(i,j+1,k)
1005 zh(i) = zh(i) + h_harm(i,k)
1007 hvel_shelf(i,k) = hvel(i,k)
1008 if (v(i,j,k) * h_delta(i,k) > 0)
then 1009 if (zh(i) * i_htbl(i) < cs%harm_BL_val)
then 1010 hvel_shelf(i,k) = min(hvel(i,k), h_harm(i,k))
1012 z2_wt = 1.0 ;
if (zh(i) * i_htbl(i) < 2.0*cs%harm_BL_val) &
1013 z2_wt = max(0.0, min(1.0, zh(i) * i_htbl(i) * i_valbl - 1.0))
1014 z2 = z2_wt * (max(zh(i), ztop_min(i) - min(zcol1(i),zcol2(i))) * i_htbl(i))
1015 topfn = 1.0 / (1.0 + 0.09*z2**6)
1016 hvel_shelf(i,k) = min(hvel(i,k), (1.0-topfn)*h_arith(i,k) + topfn*h_harm(i,k))
1021 call find_coupling_coef(a_shelf, hvel_shelf, do_i_shelf, h_harm, &
1022 bbl_thick, kv_bbl, z_i, h_ml, dt, j, g, gv, us, cs, &
1023 visc, forces, work_on_u=.false., obc=obc, shelf=.true.)
1025 do i=is,ie ;
if (do_i_shelf(i)) cs%a1_shelf_v(i,j) = a_shelf(i,1) ;
enddo 1029 if (do_any_shelf)
then 1030 do k=1,nz+1 ;
do i=is,ie ;
if (do_i_shelf(i))
then 1031 cs%a_v(i,j,k) = min(a_cpl_max, forces%frac_shelf_v(i,j) * a_shelf(i,k) + &
1032 (1.0-forces%frac_shelf_v(i,j)) * a_cpl(i,k))
1036 elseif (do_i(i))
then 1037 cs%a_v(i,j,k) = min(a_cpl_max, a_cpl(i,k))
1038 endif ;
enddo ;
enddo 1039 do k=1,nz ;
do i=is,ie ;
if (do_i_shelf(i))
then 1041 cs%h_v(i,j,k) = forces%frac_shelf_v(i,j) * hvel_shelf(i,k) + &
1042 (1.0-forces%frac_shelf_v(i,j)) * hvel(i,k) + h_neglect
1043 elseif (do_i(i))
then 1044 cs%h_v(i,j,k) = hvel(i,k) + h_neglect
1045 endif ;
enddo ;
enddo 1047 do k=1,nz+1 ;
do i=is,ie ;
if (do_i(i)) cs%a_v(i,j,k) = min(a_cpl_max, a_cpl(i,k)) ;
enddo ;
enddo 1048 do k=1,nz ;
do i=is,ie ;
if (do_i(i)) cs%h_v(i,j,k) = hvel(i,k) + h_neglect ;
enddo ;
enddo 1052 if (cs%id_Kv_v > 0)
then 1053 do k=1,nz ;
do i=is,ie
1054 if (do_i(i)) kv_v(i,j,k) = 0.5 * gv%H_to_Z*(cs%a_v(i,j,k)+cs%a_v(i,j,k+1)) * cs%h_v(i,j,k)
1061 call uvchksum(
"vertvisc_coef h_[uv]", cs%h_u, cs%h_v, g%HI, haloshift=0, &
1062 scale=gv%H_to_m, scalar_pair=.true.)
1063 call uvchksum(
"vertvisc_coef a_[uv]", cs%a_u, cs%a_v, g%HI, haloshift=0, &
1064 scale=us%Z_to_m*us%s_to_T, scalar_pair=.true.)
1065 if (
allocated(hml_u) .and.
allocated(hml_v)) &
1066 call uvchksum(
"vertvisc_coef hML_[uv]", hml_u, hml_v, g%HI, &
1067 haloshift=0, scale=gv%H_to_m, scalar_pair=.true.)
1071 if (
associated(visc%Kv_slow) .and. (cs%id_Kv_slow > 0)) &
1072 call post_data(cs%id_Kv_slow, visc%Kv_slow, cs%diag)
1073 if (cs%id_Kv_u > 0)
call post_data(cs%id_Kv_u, kv_u, cs%diag)
1074 if (cs%id_Kv_v > 0)
call post_data(cs%id_Kv_v, kv_v, cs%diag)
1075 if (cs%id_au_vv > 0)
call post_data(cs%id_au_vv, cs%a_u, cs%diag)
1076 if (cs%id_av_vv > 0)
call post_data(cs%id_av_vv, cs%a_v, cs%diag)
1077 if (cs%id_h_u > 0)
call post_data(cs%id_h_u, cs%h_u, cs%diag)
1078 if (cs%id_h_v > 0)
call post_data(cs%id_h_v, cs%h_v, cs%diag)
1079 if (cs%id_hML_u > 0)
call post_data(cs%id_hML_u, hml_u, cs%diag)
1080 if (cs%id_hML_v > 0)
call post_data(cs%id_hML_v, hml_v, cs%diag)
1082 if (
allocated(hml_u))
deallocate(hml_u)
1083 if (
allocated(hml_v))
deallocate(hml_v)
1085 end subroutine vertvisc_coef
1090 subroutine find_coupling_coef(a_cpl, hvel, do_i, h_harm, bbl_thick, kv_bbl, z_i, h_ml, &
1091 dt, j, G, GV, US, CS, visc, forces, work_on_u, OBC, shelf)
1092 type(ocean_grid_type),
intent(in) :: G
1093 type(verticalGrid_type),
intent(in) :: GV
1094 type(unit_scale_type),
intent(in) :: US
1095 real,
dimension(SZIB_(G),SZK_(GV)+1), &
1096 intent(out) :: a_cpl
1097 real,
dimension(SZIB_(G),SZK_(GV)), &
1099 logical,
dimension(SZIB_(G)), &
1101 real,
dimension(SZIB_(G),SZK_(GV)), &
1102 intent(in) :: h_harm
1104 real,
dimension(SZIB_(G)),
intent(in) :: bbl_thick
1105 real,
dimension(SZIB_(G)),
intent(in) :: kv_bbl
1106 real,
dimension(SZIB_(G),SZK_(GV)+1), &
1109 real,
dimension(SZIB_(G)),
intent(out) :: h_ml
1110 integer,
intent(in) :: j
1111 real,
intent(in) :: dt
1112 type(vertvisc_CS),
pointer :: CS
1113 type(vertvisc_type),
intent(in) :: visc
1114 type(mech_forcing),
intent(in) :: forces
1115 logical,
intent(in) :: work_on_u
1117 type(ocean_OBC_type),
pointer :: OBC
1118 logical,
optional,
intent(in) :: shelf
1123 real,
dimension(SZIB_(G)) :: &
1132 real,
dimension(SZIB_(G),SZK_(GV)+1) :: &
1149 logical :: do_shelf, do_OBCs
1150 integer :: i, k, is, ie, max_nk
1156 if (work_on_u)
then ; is = g%IscB ; ie = g%IecB
1157 else ; is = g%isc ; ie = g%iec ;
endif 1159 h_neglect = gv%H_subroundoff
1161 if (cs%answers_2018)
then 1165 i_amax = (1.0e-10*us%Z_to_m) * dt
1170 do_shelf = .false. ;
if (
present(shelf)) do_shelf = shelf
1172 if (
associated(obc))
then ; do_obcs = (obc%number_of_segments > 0) ;
endif 1177 do i=is,ie ; kv_tot(i,1) = 0.0 ;
enddo 1178 if ((gv%nkml>0) .or. do_shelf)
then ;
do k=2,nz ;
do i=is,ie
1179 if (do_i(i)) kv_tot(i,k) = cs%Kv
1180 enddo ;
enddo ;
else 1181 i_hmix = 1.0 / (cs%Hmix + h_neglect)
1182 do i=is,ie ; z_t(i) = h_neglect*i_hmix ;
enddo 1183 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 1184 z_t(i) = z_t(i) + h_harm(i,k-1)*i_hmix
1185 kv_tot(i,k) = cs%Kv + cs%Kvml / ((z_t(i)*z_t(i)) * &
1186 (1.0 + 0.09*z_t(i)*z_t(i)*z_t(i)*z_t(i)*z_t(i)*z_t(i)))
1187 endif ;
enddo ;
enddo 1190 do i=is,ie ;
if (do_i(i))
then 1191 if (cs%bottomdraglaw)
then 1193 if (r < bbl_thick(i))
then 1194 a_cpl(i,nz+1) = kv_bbl(i) / (i_amax*kv_bbl(i) + (r+h_neglect)*gv%H_to_Z)
1196 a_cpl(i,nz+1) = kv_bbl(i) / (i_amax*kv_bbl(i) + (bbl_thick(i)+h_neglect)*gv%H_to_Z)
1199 a_cpl(i,nz+1) = cs%Kvbbl / ((0.5*hvel(i,nz)+h_neglect)*gv%H_to_Z + i_amax*cs%Kvbbl)
1203 if (
associated(visc%Kv_shear))
then 1206 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 1207 kv_add(i,k) = 0.5*(visc%Kv_shear(i,j,k) + visc%Kv_shear(i+1,j,k))
1208 endif ;
enddo ;
enddo 1210 do i=is,ie ;
if (do_i(i) .and. (obc%segnum_u(i,j) /= obc_none))
then 1211 if (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_e)
then 1212 do k=2,nz ; kv_add(i,k) = visc%Kv_shear(i,j,k) ;
enddo 1213 elseif (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_w)
then 1214 do k=2,nz ; kv_add(i,k) = visc%Kv_shear(i+1,j,k) ;
enddo 1218 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 1219 kv_tot(i,k) = kv_tot(i,k) + kv_add(i,k)
1220 endif ;
enddo ;
enddo 1222 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 1223 kv_add(i,k) = 0.5*(visc%Kv_shear(i,j,k) + visc%Kv_shear(i,j+1,k))
1224 endif ;
enddo ;
enddo 1226 do i=is,ie ;
if (do_i(i) .and. (obc%segnum_v(i,j) /= obc_none))
then 1227 if (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_n)
then 1228 do k=2,nz ; kv_add(i,k) = visc%Kv_shear(i,j,k) ;
enddo 1229 elseif (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_s)
then 1230 do k=2,nz ; kv_add(i,k) = visc%Kv_shear(i,j+1,k) ;
enddo 1234 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 1235 kv_tot(i,k) = kv_tot(i,k) + kv_add(i,k)
1236 endif ;
enddo ;
enddo 1240 if (
associated(visc%Kv_shear_Bu))
then 1242 do k=2,nz ;
do i=is,ie ;
If (do_i(i))
then 1243 kv_tot(i,k) = kv_tot(i,k) + (0.5)*(visc%Kv_shear_Bu(i,j-1,k) + visc%Kv_shear_Bu(i,j,k))
1244 endif ;
enddo ;
enddo 1246 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 1247 kv_tot(i,k) = kv_tot(i,k) + (0.5)*(visc%Kv_shear_Bu(i-1,j,k) + visc%Kv_shear_Bu(i,j,k))
1248 endif ;
enddo ;
enddo 1252 do k=nz,2,-1 ;
do i=is,ie ;
if (do_i(i))
then 1256 botfn = 1.0 / (1.0 + 0.09*z2*z2*z2*z2*z2*z2)
1258 if (cs%bottomdraglaw)
then 1259 kv_tot(i,k) = kv_tot(i,k) + (kv_bbl(i) - cs%Kv)*botfn
1260 r = 0.5*(hvel(i,k) + hvel(i,k-1))
1261 if (r > bbl_thick(i))
then 1262 h_shear = ((1.0 - botfn) * r + botfn*bbl_thick(i)) + h_neglect
1264 h_shear = r + h_neglect
1267 kv_tot(i,k) = kv_tot(i,k) + (cs%Kvbbl-cs%Kv)*botfn
1268 h_shear = 0.5*(hvel(i,k) + hvel(i,k-1) + h_neglect)
1272 a_cpl(i,k) = kv_tot(i,k) / (h_shear*gv%H_to_Z + i_amax*kv_tot(i,k))
1273 endif ;
enddo ;
enddo 1277 do i=is,ie ;
if (do_i(i))
then 1279 kv_tbl(i) = visc%Kv_tbl_shelf_u(i,j)
1280 tbl_thick(i) = visc%tbl_thick_shelf_u(i,j) * gv%Z_to_H + h_neglect
1282 kv_tbl(i) = visc%Kv_tbl_shelf_v(i,j)
1283 tbl_thick(i) = visc%tbl_thick_shelf_v(i,j) * gv%Z_to_H + h_neglect
1288 if (0.5*hvel(i,1) > tbl_thick(i))
then 1289 a_cpl(i,1) = kv_tbl(i) / (tbl_thick(i)*gv%H_to_Z + i_amax*kv_tbl(i))
1291 a_cpl(i,1) = kv_tbl(i) / ((0.5*hvel(i,1)+h_neglect)*gv%H_to_Z + i_amax*kv_tbl(i))
1295 do k=2,nz ;
do i=is,ie ;
if (do_i(i))
then 1296 z_t(i) = z_t(i) + hvel(i,k-1) / tbl_thick(i)
1297 topfn = 1.0 / (1.0 + 0.09 * z_t(i)**6)
1299 r = 0.5*(hvel(i,k)+hvel(i,k-1))
1300 if (r > tbl_thick(i))
then 1301 h_shear = ((1.0 - topfn) * r + topfn*tbl_thick(i)) + h_neglect
1303 h_shear = r + h_neglect
1306 kv_top = topfn * kv_tbl(i)
1307 a_cpl(i,k) = a_cpl(i,k) + kv_top / (h_shear*gv%H_to_Z + i_amax*kv_top)
1308 endif ;
enddo ;
enddo 1309 elseif (cs%dynamic_viscous_ML .or. (gv%nkml>0))
then 1311 do i=is,ie ;
if (do_i(i))
then 1312 if (gv%nkml>0) nk_visc(i) =
real(GV%nkml+1)
1314 u_star(i) = 0.5*(forces%ustar(i,j) + forces%ustar(i+1,j))
1315 absf(i) = 0.5*(abs(g%CoriolisBu(i,j-1)) + abs(g%CoriolisBu(i,j)))
1316 if (cs%dynamic_viscous_ML) nk_visc(i) = visc%nkml_visc_u(i,j) + 1
1318 u_star(i) = 0.5*(forces%ustar(i,j) + forces%ustar(i,j+1))
1319 absf(i) = 0.5*(abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j)))
1320 if (cs%dynamic_viscous_ML) nk_visc(i) = visc%nkml_visc_v(i,j) + 1
1322 h_ml(i) = h_neglect ; z_t(i) = 0.0
1323 max_nk = max(max_nk,ceiling(nk_visc(i) - 1.0))
1326 if (do_obcs)
then ;
if (work_on_u)
then 1327 do i=is,ie ;
if (do_i(i) .and. (obc%segnum_u(i,j) /= obc_none))
then 1328 if (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_e) &
1329 u_star(i) = forces%ustar(i,j)
1330 if (obc%segment(obc%segnum_u(i,j))%direction == obc_direction_w) &
1331 u_star(i) = forces%ustar(i+1,j)
1334 do i=is,ie ;
if (do_i(i) .and. (obc%segnum_v(i,j) /= obc_none))
then 1335 if (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_n) &
1336 u_star(i) = forces%ustar(i,j)
1337 if (obc%segment(obc%segnum_v(i,j))%direction == obc_direction_s) &
1338 u_star(i) = forces%ustar(i,j+1)
1342 do k=1,max_nk ;
do i=is,ie ;
if (do_i(i))
then 1343 if (k+1 <= nk_visc(i))
then 1344 h_ml(i) = h_ml(i) + hvel(i,k)
1345 elseif (k < nk_visc(i))
then 1346 h_ml(i) = h_ml(i) + (nk_visc(i) - k) * hvel(i,k)
1348 endif ;
enddo ;
enddo 1350 do k=2,max_nk ;
do i=is,ie ;
if (do_i(i))
then ;
if (k < nk_visc(i))
then 1352 z_t(i) = z_t(i) + hvel(i,k-1)
1353 temp1 = (z_t(i)*h_ml(i) - z_t(i)*z_t(i))*gv%H_to_Z
1356 visc_ml = u_star(i) * 0.41 * (temp1*u_star(i)) / (absf(i)*temp1 + (h_ml(i)+h_neglect)*u_star(i))
1357 a_ml = visc_ml / (0.25*(hvel(i,k)+hvel(i,k-1) + h_neglect) * gv%H_to_Z + 0.5*i_amax*visc_ml)
1359 if (a_ml > a_cpl(i,k)) a_cpl(i,k) = a_ml
1360 endif ;
endif ;
enddo ;
enddo 1363 end subroutine find_coupling_coef
1368 subroutine vertvisc_limit_vel(u, v, h, ADp, CDp, forces, visc, dt, G, GV, US, CS)
1372 real,
dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
1374 real,
dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
1376 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1382 real,
intent(in) :: dt
1392 real :: vel_report(szib_(g),szjb_(g))
1393 real :: u_old(szib_(g),szj_(g),szk_(g))
1394 real :: v_old(szi_(g),szjb_(g),szk_(g))
1395 logical :: trunc_any, dowrite(szib_(g),szjb_(g))
1396 integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz
1397 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1398 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1401 truncvel = 0.9*maxvel
1402 h_report = 6.0 * gv%Angstrom_H
1403 dt_rho0 = (us%L_T_to_m_s*us%Z_to_m) * dt / gv%Rho0
1405 if (len_trim(cs%u_trunc_file) > 0)
then 1409 do i=isq,ieq ; dowrite(i,j) = .false. ;
enddo 1410 if (cs%CFL_based_trunc)
then 1411 do i=isq,ieq ; vel_report(i,j) = 3.0e8*us%m_s_to_L_T ;
enddo 1412 do k=1,nz ;
do i=isq,ieq
1413 if (abs(u(i,j,k)) < cs%vel_underflow) u(i,j,k) = 0.0
1414 if (u(i,j,k) < 0.0)
then 1415 cfl = (-u(i,j,k) * dt) * (g%dy_Cu(i,j) * g%IareaT(i+1,j))
1417 cfl = (u(i,j,k) * dt) * (g%dy_Cu(i,j) * g%IareaT(i,j))
1419 if (cfl > cs%CFL_trunc) trunc_any = .true.
1420 if (cfl > cs%CFL_report)
then 1421 dowrite(i,j) = .true.
1422 vel_report(i,j) = min(vel_report(i,j), abs(u(i,j,k)))
1426 do i=isq,ieq; vel_report(i,j) = maxvel;
enddo 1427 do k=1,nz ;
do i=isq,ieq
1428 if (abs(u(i,j,k)) < cs%vel_underflow)
then ; u(i,j,k) = 0.0
1429 elseif (abs(u(i,j,k)) > maxvel)
then 1430 dowrite(i,j) = .true. ; trunc_any = .true.
1435 do i=isq,ieq ;
if (dowrite(i,j))
then 1436 u_old(i,j,:) = u(i,j,:)
1439 if (trunc_any)
then ;
if (cs%CFL_based_trunc)
then 1440 do k=1,nz ;
do i=isq,ieq
1441 if ((u(i,j,k) * (dt * g%dy_Cu(i,j))) * g%IareaT(i+1,j) < -cs%CFL_trunc)
then 1442 u(i,j,k) = (-0.9*cs%CFL_trunc) * (g%areaT(i+1,j) / (dt * g%dy_Cu(i,j)))
1443 if (h(i,j,k) + h(i+1,j,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1444 elseif ((u(i,j,k) * (dt * g%dy_Cu(i,j))) * g%IareaT(i,j) > cs%CFL_trunc)
then 1445 u(i,j,k) = (0.9*cs%CFL_trunc) * (g%areaT(i,j) / (dt * g%dy_Cu(i,j)))
1446 if (h(i,j,k) + h(i+1,j,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1450 do k=1,nz ;
do i=isq,ieq ;
if (abs(u(i,j,k)) > maxvel)
then 1451 u(i,j,k) = sign(truncvel,u(i,j,k))
1452 if (h(i,j,k) + h(i+1,j,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1453 endif ;
enddo ;
enddo 1457 if (cs%CFL_based_trunc)
then 1459 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1460 if (abs(u(i,j,k)) < cs%vel_underflow)
then ; u(i,j,k) = 0.0
1461 elseif ((u(i,j,k) * (dt * g%dy_Cu(i,j))) * g%IareaT(i+1,j) < -cs%CFL_trunc)
then 1462 u(i,j,k) = (-0.9*cs%CFL_trunc) * (g%areaT(i+1,j) / (dt * g%dy_Cu(i,j)))
1463 if (h(i,j,k) + h(i+1,j,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1464 elseif ((u(i,j,k) * (dt * g%dy_Cu(i,j))) * g%IareaT(i,j) > cs%CFL_trunc)
then 1465 u(i,j,k) = (0.9*cs%CFL_trunc) * (g%areaT(i,j) / (dt * g%dy_Cu(i,j)))
1466 if (h(i,j,k) + h(i+1,j,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1468 enddo ;
enddo ;
enddo 1471 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1472 if (abs(u(i,j,k)) < cs%vel_underflow)
then ; u(i,j,k) = 0.0
1473 elseif (abs(u(i,j,k)) > maxvel)
then 1474 u(i,j,k) = sign(truncvel, u(i,j,k))
1475 if (h(i,j,k) + h(i+1,j,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1477 enddo ;
enddo ;
enddo 1481 if (len_trim(cs%u_trunc_file) > 0)
then 1482 do j=js,je;
do i=isq,ieq ;
if (dowrite(i,j))
then 1485 call write_u_accel(i, j, u_old, h, adp, cdp, dt, g, gv, us, cs%PointAccel_CSp, &
1486 vel_report(i,j), forces%taux(i,j)*dt_rho0, a=cs%a_u, hv=cs%h_u)
1487 endif ;
enddo ;
enddo 1490 if (len_trim(cs%v_trunc_file) > 0)
then 1494 do i=is,ie ; dowrite(i,j) = .false. ;
enddo 1495 if (cs%CFL_based_trunc)
then 1496 do i=is,ie ; vel_report(i,j) = 3.0e8*us%m_s_to_L_T ;
enddo 1497 do k=1,nz ;
do i=is,ie
1498 if (abs(v(i,j,k)) < cs%vel_underflow) v(i,j,k) = 0.0
1499 if (v(i,j,k) < 0.0)
then 1500 cfl = (-v(i,j,k) * dt) * (g%dx_Cv(i,j) * g%IareaT(i,j+1))
1502 cfl = (v(i,j,k) * dt) * (g%dx_Cv(i,j) * g%IareaT(i,j))
1504 if (cfl > cs%CFL_trunc) trunc_any = .true.
1505 if (cfl > cs%CFL_report)
then 1506 dowrite(i,j) = .true.
1507 vel_report(i,j) = min(vel_report(i,j), abs(v(i,j,k)))
1511 do i=is,ie ; vel_report(i,j) = maxvel ;
enddo 1512 do k=1,nz ;
do i=is,ie
1513 if (abs(v(i,j,k)) < cs%vel_underflow)
then ; v(i,j,k) = 0.0
1514 elseif (abs(v(i,j,k)) > maxvel)
then 1515 dowrite(i,j) = .true. ; trunc_any = .true.
1520 do i=is,ie ;
if (dowrite(i,j))
then 1521 v_old(i,j,:) = v(i,j,:)
1524 if (trunc_any)
then ;
if (cs%CFL_based_trunc)
then 1525 do k=1,nz;
do i=is,ie
1526 if ((v(i,j,k) * (dt * g%dx_Cv(i,j))) * g%IareaT(i,j+1) < -cs%CFL_trunc)
then 1527 v(i,j,k) = (-0.9*cs%CFL_trunc) * (g%areaT(i,j+1) / (dt * g%dx_Cv(i,j)))
1528 if (h(i,j,k) + h(i,j+1,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1529 elseif ((v(i,j,k) * (dt * g%dx_Cv(i,j))) * g%IareaT(i,j) > cs%CFL_trunc)
then 1530 v(i,j,k) = (0.9*cs%CFL_trunc) * (g%areaT(i,j) / (dt * g%dx_Cv(i,j)))
1531 if (h(i,j,k) + h(i,j+1,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1535 do k=1,nz ;
do i=is,ie ;
if (abs(v(i,j,k)) > maxvel)
then 1536 v(i,j,k) = sign(truncvel,v(i,j,k))
1537 if (h(i,j,k) + h(i,j+1,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1538 endif ;
enddo ;
enddo 1542 if (cs%CFL_based_trunc)
then 1544 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1545 if (abs(v(i,j,k)) < cs%vel_underflow)
then ; v(i,j,k) = 0.0
1546 elseif ((v(i,j,k) * (dt * g%dx_Cv(i,j))) * g%IareaT(i,j+1) < -cs%CFL_trunc)
then 1547 v(i,j,k) = (-0.9*cs%CFL_trunc) * (g%areaT(i,j+1) / (dt * g%dx_Cv(i,j)))
1548 if (h(i,j,k) + h(i,j+1,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1549 elseif ((v(i,j,k) * (dt * g%dx_Cv(i,j))) * g%IareaT(i,j) > cs%CFL_trunc)
then 1550 v(i,j,k) = (0.9*cs%CFL_trunc) * (g%areaT(i,j) / (dt * g%dx_Cv(i,j)))
1551 if (h(i,j,k) + h(i,j+1,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1553 enddo ;
enddo ;
enddo 1556 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1557 if (abs(v(i,j,k)) < cs%vel_underflow)
then ; v(i,j,k) = 0.0
1558 elseif (abs(v(i,j,k)) > maxvel)
then 1559 v(i,j,k) = sign(truncvel, v(i,j,k))
1560 if (h(i,j,k) + h(i,j+1,k) > h_report) cs%ntrunc = cs%ntrunc + 1
1562 enddo ;
enddo ;
enddo 1566 if (len_trim(cs%v_trunc_file) > 0)
then 1567 do j=jsq,jeq;
do i=is,ie ;
if (dowrite(i,j))
then 1570 call write_v_accel(i, j, v_old, h, adp, cdp, dt, g, gv, us, cs%PointAccel_CSp, &
1571 vel_report(i,j), forces%tauy(i,j)*dt_rho0, a=cs%a_v, hv=cs%h_v)
1572 endif ;
enddo ;
enddo 1575 end subroutine vertvisc_limit_vel
1578 subroutine vertvisc_init(MIS, Time, G, GV, US, param_file, diag, ADp, dirs, &
1581 target,
intent(in) :: mis
1584 type(time_type),
target,
intent(in) :: time
1589 type(
diag_ctrl),
target,
intent(inout) :: diag
1592 integer,
target,
intent(inout) :: ntrunc
1597 real :: hmix_str_dflt
1600 logical :: default_2018_answers
1601 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, nz
1603 #include "version_variable.h" 1604 character(len=40) :: mdl =
"MOM_vert_friction" 1605 character(len=40) :: thickness_units
1607 if (
associated(cs))
then 1608 call mom_error(warning,
"vertvisc_init called with an associated "// &
1609 "control structure.")
1614 if (gv%Boussinesq) then; thickness_units =
"m" 1615 else; thickness_units =
"kg m-2";
endif 1617 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
1618 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1620 cs%diag => diag ; cs%ntrunc => ntrunc ; ntrunc = 0
1623 call log_version(param_file, mdl, version,
"", log_to_all=.true., debugging=.true.)
1624 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
1625 "This sets the default value for the various _2018_ANSWERS parameters.", &
1627 call get_param(param_file, mdl,
"VERT_FRICTION_2018_ANSWERS", cs%answers_2018, &
1628 "If true, use the order of arithmetic and expressions that recover the answers "//&
1629 "from the end of 2018. Otherwise, use expressions that do not use an arbitrary "//&
1630 "hard-coded maximum viscous coupling coefficient between layers.", &
1631 default=default_2018_answers)
1632 call get_param(param_file, mdl,
"BOTTOMDRAGLAW", cs%bottomdraglaw, &
1633 "If true, the bottom stress is calculated with a drag "//&
1634 "law of the form c_drag*|u|*u. The velocity magnitude "//&
1635 "may be an assumed value or it may be based on the "//&
1636 "actual velocity in the bottommost HBBL, depending on "//&
1637 "LINEAR_DRAG.", default=.true.)
1638 call get_param(param_file, mdl,
"CHANNEL_DRAG", cs%Channel_drag, &
1639 "If true, the bottom drag is exerted directly on each "//&
1640 "layer proportional to the fraction of the bottom it "//&
1641 "overlies.", default=.false.)
1642 call get_param(param_file, mdl,
"DIRECT_STRESS", cs%direct_stress, &
1643 "If true, the wind stress is distributed over the "//&
1644 "topmost HMIX_STRESS of fluid (like in HYCOM), and KVML "//&
1645 "may be set to a very small value.", default=.false.)
1646 call get_param(param_file, mdl,
"DYNAMIC_VISCOUS_ML", cs%dynamic_viscous_ML, &
1647 "If true, use a bulk Richardson number criterion to "//&
1648 "determine the mixed layer thickness for viscosity.", &
1650 call get_param(param_file, mdl,
"U_TRUNC_FILE", cs%u_trunc_file, &
1651 "The absolute path to a file into which the accelerations "//&
1652 "leading to zonal velocity truncations are written. "//&
1653 "Undefine this for efficiency if this diagnostic is not "//&
1654 "needed.", default=
" ", debuggingparam=.true.)
1655 call get_param(param_file, mdl,
"V_TRUNC_FILE", cs%v_trunc_file, &
1656 "The absolute path to a file into which the accelerations "//&
1657 "leading to meridional velocity truncations are written. "//&
1658 "Undefine this for efficiency if this diagnostic is not "//&
1659 "needed.", default=
" ", debuggingparam=.true.)
1660 call get_param(param_file, mdl,
"HARMONIC_VISC", cs%harmonic_visc, &
1661 "If true, use the harmonic mean thicknesses for "//&
1662 "calculating the vertical viscosity.", default=.false.)
1663 call get_param(param_file, mdl,
"HARMONIC_BL_SCALE", cs%harm_BL_val, &
1664 "A scale to determine when water is in the boundary "//&
1665 "layers based solely on harmonic mean thicknesses for "//&
1666 "the purpose of determining the extent to which the "//&
1667 "thicknesses used in the viscosities are upwinded.", &
1668 default=0.0, units=
"nondim")
1669 call get_param(param_file, mdl,
"DEBUG", cs%debug, default=.false.)
1672 call get_param(param_file, mdl,
"HMIX_FIXED", cs%Hmix, &
1673 "The prescribed depth over which the near-surface "//&
1674 "viscosity and diffusivity are elevated when the bulk "//&
1675 "mixed layer is not used.", units=
"m", scale=gv%m_to_H, &
1676 unscaled=hmix_m, fail_if_missing=.true.)
1677 if (cs%direct_stress)
then 1678 if (gv%nkml < 1)
then 1679 call get_param(param_file, mdl,
"HMIX_STRESS", cs%Hmix_stress, &
1680 "The depth over which the wind stress is applied if "//&
1681 "DIRECT_STRESS is true.", units=
"m", default=hmix_m, scale=gv%m_to_H)
1683 call get_param(param_file, mdl,
"HMIX_STRESS", cs%Hmix_stress, &
1684 "The depth over which the wind stress is applied if "//&
1685 "DIRECT_STRESS is true.", units=
"m", fail_if_missing=.true., scale=gv%m_to_H)
1687 if (cs%Hmix_stress <= 0.0)
call mom_error(fatal,
"vertvisc_init: " // &
1688 "HMIX_STRESS must be set to a positive value if DIRECT_STRESS is true.")
1690 call get_param(param_file, mdl,
"KV", cs%Kv, &
1691 "The background kinematic viscosity in the interior. "//&
1692 "The molecular value, ~1e-6 m2 s-1, may be used.", &
1693 units=
"m2 s-1", fail_if_missing=.true., scale=us%m2_s_to_Z2_T, unscaled=kv_dflt)
1695 if (gv%nkml < 1)
call get_param(param_file, mdl,
"KVML", cs%Kvml, &
1696 "The kinematic viscosity in the mixed layer. A typical "//&
1697 "value is ~1e-2 m2 s-1. KVML is not used if "//&
1698 "BULKMIXEDLAYER is true. The default is set by KV.", &
1699 units=
"m2 s-1", default=kv_dflt, scale=us%m2_s_to_Z2_T)
1700 if (.not.cs%bottomdraglaw)
call get_param(param_file, mdl,
"KVBBL", cs%Kvbbl, &
1701 "The kinematic viscosity in the benthic boundary layer. "//&
1702 "A typical value is ~1e-2 m2 s-1. KVBBL is not used if "//&
1703 "BOTTOMDRAGLAW is true. The default is set by KV.", &
1704 units=
"m2 s-1", default=kv_dflt, scale=us%m2_s_to_Z2_T)
1705 call get_param(param_file, mdl,
"HBBL", cs%Hbbl, &
1706 "The thickness of a bottom boundary layer with a "//&
1707 "viscosity of KVBBL if BOTTOMDRAGLAW is not defined, or "//&
1708 "the thickness over which near-bottom velocities are "//&
1709 "averaged for the drag law if BOTTOMDRAGLAW is defined "//&
1710 "but LINEAR_DRAG is not.", units=
"m", fail_if_missing=.true., scale=gv%m_to_H)
1711 call get_param(param_file, mdl,
"MAXVEL", cs%maxvel, &
1712 "The maximum velocity allowed before the velocity "//&
1713 "components are truncated.", units=
"m s-1", default=3.0e8, scale=us%m_s_to_L_T)
1714 call get_param(param_file, mdl,
"CFL_BASED_TRUNCATIONS", cs%CFL_based_trunc, &
1715 "If true, base truncations on the CFL number, and not an "//&
1716 "absolute speed.", default=.true.)
1717 call get_param(param_file, mdl,
"CFL_TRUNCATE", cs%CFL_trunc, &
1718 "The value of the CFL number that will cause velocity "//&
1719 "components to be truncated; instability can occur past 0.5.", &
1720 units=
"nondim", default=0.5)
1721 call get_param(param_file, mdl,
"CFL_REPORT", cs%CFL_report, &
1722 "The value of the CFL number that causes accelerations "//&
1723 "to be reported; the default is CFL_TRUNCATE.", &
1724 units=
"nondim", default=cs%CFL_trunc)
1725 call get_param(param_file, mdl,
"CFL_TRUNCATE_RAMP_TIME", cs%truncRampTime, &
1726 "The time over which the CFL truncation value is ramped "//&
1727 "up at the beginning of the run.", &
1728 units=
"s", default=0.)
1729 cs%CFL_truncE = cs%CFL_trunc
1730 call get_param(param_file, mdl,
"CFL_TRUNCATE_START", cs%CFL_truncS, &
1731 "The start value of the truncation CFL number used when "//&
1732 "ramping up CFL_TRUNC.", &
1733 units=
"nondim", default=0.)
1734 call get_param(param_file, mdl,
"STOKES_MIXING_COMBINED", cs%StokesMixing, &
1735 "Flag to use Stokes drift Mixing via the Lagrangian "//&
1736 " current (Eulerian plus Stokes drift). "//&
1737 " Still needs work and testing, so not recommended for use.",&
1748 if (cs%StokesMixing)
then 1749 call mom_error(fatal,
"Stokes mixing requires user intervention in the code.\n"//&
1750 " Model now exiting. See MOM_vert_friction.F90 for \n"//&
1751 " details (search 'BGR 04/04/2018' to locate comment).")
1753 call get_param(param_file, mdl,
"VEL_UNDERFLOW", cs%vel_underflow, &
1754 "A negligibly small velocity magnitude below which velocity "//&
1755 "components are set to 0. A reasonable value might be "//&
1756 "1e-30 m/s, which is less than an Angstrom divided by "//&
1757 "the age of the universe.", units=
"m s-1", default=0.0, scale=us%m_s_to_L_T)
1759 alloc_(cs%a_u(isdb:iedb,jsd:jed,nz+1)) ; cs%a_u(:,:,:) = 0.0
1760 alloc_(cs%h_u(isdb:iedb,jsd:jed,nz)) ; cs%h_u(:,:,:) = 0.0
1761 alloc_(cs%a_v(isd:ied,jsdb:jedb,nz+1)) ; cs%a_v(:,:,:) = 0.0
1762 alloc_(cs%h_v(isd:ied,jsdb:jedb,nz)) ; cs%h_v(:,:,:) = 0.0
1764 cs%id_Kv_slow = register_diag_field(
'ocean_model',
'Kv_slow', diag%axesTi, time, &
1765 'Slow varying vertical viscosity',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
1767 cs%id_Kv_u = register_diag_field(
'ocean_model',
'Kv_u', diag%axesCuL, time, &
1768 'Total vertical viscosity at u-points',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
1770 cs%id_Kv_v = register_diag_field(
'ocean_model',
'Kv_v', diag%axesCvL, time, &
1771 'Total vertical viscosity at v-points',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
1773 cs%id_au_vv = register_diag_field(
'ocean_model',
'au_visc', diag%axesCui, time, &
1774 'Zonal Viscous Vertical Coupling Coefficient',
'm s-1', conversion=us%Z_to_m*us%s_to_T)
1776 cs%id_av_vv = register_diag_field(
'ocean_model',
'av_visc', diag%axesCvi, time, &
1777 'Meridional Viscous Vertical Coupling Coefficient',
'm s-1', conversion=us%Z_to_m*us%s_to_T)
1779 cs%id_h_u = register_diag_field(
'ocean_model',
'Hu_visc', diag%axesCuL, time, &
1780 'Thickness at Zonal Velocity Points for Viscosity', thickness_units, &
1781 conversion=gv%H_to_m)
1783 cs%id_h_v = register_diag_field(
'ocean_model',
'Hv_visc', diag%axesCvL, time, &
1784 'Thickness at Meridional Velocity Points for Viscosity', thickness_units, &
1785 conversion=gv%H_to_m)
1787 cs%id_hML_u = register_diag_field(
'ocean_model',
'HMLu_visc', diag%axesCu1, time, &
1788 'Mixed Layer Thickness at Zonal Velocity Points for Viscosity', thickness_units, &
1789 conversion=gv%H_to_m)
1791 cs%id_hML_v = register_diag_field(
'ocean_model',
'HMLv_visc', diag%axesCv1, time, &
1792 'Mixed Layer Thickness at Meridional Velocity Points for Viscosity', thickness_units, &
1793 conversion=gv%H_to_m)
1795 cs%id_du_dt_visc = register_diag_field(
'ocean_model',
'du_dt_visc', diag%axesCuL, &
1796 time,
'Zonal Acceleration from Vertical Viscosity',
'm s-2', conversion=us%L_T2_to_m_s2)
1797 if (cs%id_du_dt_visc > 0)
call safe_alloc_ptr(adp%du_dt_visc,isdb,iedb,jsd,jed,nz)
1798 cs%id_dv_dt_visc = register_diag_field(
'ocean_model',
'dv_dt_visc', diag%axesCvL, &
1799 time,
'Meridional Acceleration from Vertical Viscosity',
'm s-2', conversion=us%L_T2_to_m_s2)
1800 if (cs%id_dv_dt_visc > 0)
call safe_alloc_ptr(adp%dv_dt_visc,isd,ied,jsdb,jedb,nz)
1802 cs%id_taux_bot = register_diag_field(
'ocean_model',
'taux_bot', diag%axesCu1, &
1803 time,
'Zonal Bottom Stress from Ocean to Earth',
'Pa', &
1804 conversion=us%RZ_to_kg_m2*us%L_T2_to_m_s2)
1805 cs%id_tauy_bot = register_diag_field(
'ocean_model',
'tauy_bot', diag%axesCv1, &
1806 time,
'Meridional Bottom Stress from Ocean to Earth',
'Pa', &
1807 conversion=us%RZ_to_kg_m2*us%L_T2_to_m_s2)
1827 cs%id_hf_du_dt_visc_2d = register_diag_field(
'ocean_model',
'hf_du_dt_visc_2d', diag%axesCu1, time, &
1828 'Depth-sum Fractional Thickness-weighted Zonal Acceleration from Vertical Viscosity',
'm s-2', &
1829 conversion=us%L_T2_to_m_s2)
1830 if (cs%id_hf_du_dt_visc_2d > 0)
then 1831 call safe_alloc_ptr(adp%du_dt_visc,isdb,iedb,jsd,jed,nz)
1832 call safe_alloc_ptr(adp%diag_hfrac_u,isdb,iedb,jsd,jed,nz)
1835 cs%id_hf_dv_dt_visc_2d = register_diag_field(
'ocean_model',
'hf_dv_dt_visc_2d', diag%axesCv1, time, &
1836 'Depth-sum Fractional Thickness-weighted Meridional Acceleration from Vertical Viscosity',
'm s-2', &
1837 conversion=us%L_T2_to_m_s2)
1838 if (cs%id_hf_dv_dt_visc_2d > 0)
then 1839 call safe_alloc_ptr(adp%dv_dt_visc,isd,ied,jsdb,jedb,nz)
1840 call safe_alloc_ptr(adp%diag_hfrac_v,isd,ied,jsd,jedb,nz)
1843 if ((len_trim(cs%u_trunc_file) > 0) .or. (len_trim(cs%v_trunc_file) > 0)) &
1844 call pointaccel_init(mis, time, g, param_file, diag, dirs, cs%PointAccel_CSp)
1846 end subroutine vertvisc_init
1851 subroutine updatecfltruncationvalue(Time, CS, activate)
1852 type(time_type),
target,
intent(in) :: time
1854 logical,
optional,
intent(in) :: activate
1858 real :: deltatime, wghta
1859 character(len=12) :: msg
1861 if (cs%truncRampTime==0.)
return 1865 if (
present(activate))
then 1867 cs%rampStartTime = time
1868 cs%CFLrampingIsActivated = .true.
1871 if (.not.cs%CFLrampingIsActivated)
return 1872 deltatime = max( 0., time_type_to_real( time - cs%rampStartTime ) )
1873 if (deltatime >= cs%truncRampTime)
then 1874 cs%CFL_trunc = cs%CFL_truncE
1875 cs%truncRampTime = 0.
1877 wghta = min( 1., deltatime / cs%truncRampTime )
1880 wghta = 1. - ( (1. - wghta)**2 )
1881 cs%CFL_trunc = cs%CFL_truncS + wghta * ( cs%CFL_truncE - cs%CFL_truncS )
1883 write(msg(1:12),
'(es12.3)') cs%CFL_trunc
1884 call mom_error(note,
"MOM_vert_friction: updateCFLtruncationValue set CFL"// &
1885 " limit to "//trim(msg))
1886 end subroutine updatecfltruncationvalue
1889 subroutine vertvisc_end(CS)
1893 dealloc_(cs%a_u) ; dealloc_(cs%h_u)
1894 dealloc_(cs%a_v) ; dealloc_(cs%h_v)
1895 if (
associated(cs%a1_shelf_u))
deallocate(cs%a1_shelf_u)
1896 if (
associated(cs%a1_shelf_v))
deallocate(cs%a1_shelf_v)
1898 end subroutine vertvisc_end
Wraps the FMS time manager functions.
This module implements boundary forcing for MOM6.
Ocean grid type. See mom_grid for details.
A structure that can be parsed to read and document run-time parameters.
Implements vertical viscosity (vertvisc)
Provides the ocean grid type.
Vertical viscosities, drag coefficients, and related fields.
The control structure for the MOM_PointAccel module.
The MOM6 facility to parse input files for runtime parameters.
Pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances.
Structure that contains pointers to the mechanical forcing at the surface used to drive the liquid oc...
Interface for surface waves.
Make a diagnostic available for averaging or output.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances.
Routines for error handling and I/O management.
Pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM.F90.
Container for all surface wave related parameters.
The control structure with parameters and memory for the MOM_vert_friction module.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Debug accelerations at a given point.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Controls where open boundary conditions are applied.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Do a halo update on an array.
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.