25 implicit none ;
private
27 #include <MOM_memory.h>
29 public horizontal_viscosity, hor_visc_init, hor_visc_end
45 logical :: better_bound_kh
49 logical :: better_bound_ah
57 logical :: smagorinsky_kh
59 logical :: smagorinsky_ah
63 logical :: modified_leith
65 logical :: use_beta_in_leith
68 logical :: use_qg_leith_visc
70 logical :: bound_coriolis
73 logical :: use_kh_bg_2d
76 logical :: use_land_mask
81 logical :: anisotropic
82 logical :: add_les_viscosity
85 logical :: dynamic_aniso
87 logical :: res_scale_meke
90 logical :: answers_2018
95 real :: gme_efficiency
102 real allocable_,
dimension(NIMEM_,NJMEM_) :: kh_bg_xx
106 real allocable_,
dimension(NIMEM_,NJMEM_) :: kh_bg_2d
110 real allocable_,
dimension(NIMEM_,NJMEM_) :: ah_bg_xx
114 real allocable_,
dimension(NIMEM_,NJMEM_) :: reduction_xx
117 real allocable_,
dimension(NIMEM_,NJMEM_) :: &
118 kh_max_xx, & !< The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].
119 ah_max_xx, & !< The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].
120 n1n2_h, & !< Factor n1*n2 in the anisotropic direction tensor at h-points
121 n1n1_m_n2n2_h, & !< Factor n1**2-n2**2 in the anisotropic direction tensor at h-points
122 grid_sp_h2, & !< Harmonic mean of the squares of the grid [L2 ~> m2]
124 real allocable_,
dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: kh_bg_xy
128 real allocable_,
dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: ah_bg_xy
132 real allocable_,
dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: reduction_xy
135 real allocable_,
dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
136 kh_max_xy, & !< The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].
137 ah_max_xy, & !< The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].
138 n1n2_q, & !< Factor n1*n2 in the anisotropic direction tensor at q-points
141 real allocable_,
dimension(NIMEM_,NJMEM_) :: &
142 dx2h, & !< Pre-calculated dx^2 at h points [L2 ~> m2]
143 dy2h, & !< Pre-calculated dy^2 at h points [L2 ~> m2]
144 dx_dyt, & !< Pre-calculated dx/dy at h points [nondim]
146 real allocable_,
dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
147 dx2q, & !< Pre-calculated dx^2 at q points [L2 ~> m2]
148 dy2q, & !< Pre-calculated dy^2 at q points [L2 ~> m2]
149 dx_dybu, & !< Pre-calculated dx/dy at q points [nondim]
151 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_) :: &
152 idx2dycu, & !< 1/(dx^2 dy) at u points [L-3 ~> m-3]
154 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_) :: &
155 idx2dycv, & !< 1/(dx^2 dy) at v points [L-3 ~> m-3]
160 real allocable_,
dimension(NIMEM_,NJMEM_) :: &
161 laplac2_const_xx, & !< Laplacian metric-dependent constants [L2 ~> m2]
162 biharm6_const_xx, & !< Biharmonic metric-dependent constants [L6 ~> m6]
163 laplac3_const_xx, & !< Laplacian metric-dependent constants [L3 ~> m3]
164 biharm_const_xx, & !< Biharmonic metric-dependent constants [L4 ~> m4]
165 biharm_const2_xx, & !< Biharmonic metric-dependent constants [T L4 ~> s m4]
168 real allocable_,
dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
169 laplac2_const_xy, & !< Laplacian metric-dependent constants [L2 ~> m2]
170 biharm6_const_xy, & !< Biharmonic metric-dependent constants [L6 ~> m6]
171 laplac3_const_xy, & !< Laplacian metric-dependent constants [L3 ~> m3]
172 biharm_const_xy, & !< Biharmonic metric-dependent constants [L4 ~> m4]
173 biharm_const2_xy, & !< Biharmonic metric-dependent constants [T L4 ~> s m4]
185 integer :: id_grid_re_ah = -1, id_grid_re_kh = -1
186 integer :: id_diffu = -1, id_diffv = -1
188 integer :: id_hf_diffu_2d = -1, id_hf_diffv_2d = -1
189 integer :: id_ah_h = -1, id_ah_q = -1
190 integer :: id_kh_h = -1, id_kh_q = -1
191 integer :: id_gme_coeff_h = -1, id_gme_coeff_q = -1
192 integer :: id_vort_xy_q = -1, id_div_xx_h = -1
193 integer :: id_sh_xy_q = -1, id_sh_xx_h = -1
194 integer :: id_frictwork = -1, id_frictworkintz = -1
195 integer :: id_frictwork_gme = -1
215 subroutine horizontal_viscosity(u, v, h, diffu, diffv, MEKE, VarMix, G, GV, US, &
216 CS, OBC, BT, TD, ADp)
219 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
221 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
223 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
225 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
228 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
246 real,
dimension(SZIB_(G),SZJ_(G)) :: &
252 real,
dimension(SZI_(G),SZJB_(G)) :: &
258 real,
dimension(SZI_(G),SZJ_(G)) :: &
270 grad_vort_mag_h_2d, &
279 real,
allocatable,
dimension(:,:) :: hf_diffu_2d
280 real,
allocatable,
dimension(:,:) :: hf_diffv_2d
282 real,
dimension(SZIB_(G),SZJB_(G)) :: &
284 ddel2vdx, ddel2udy, &
295 grad_vort_mag_q_2d, &
304 real,
dimension(SZIB_(G),SZJB_(G),SZK_(G)) :: &
309 gme_coeff_q, & !< GME coeff. at q-points [L2 T-1 ~> m2 s-1]
312 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)+1) :: &
314 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)+1) :: &
316 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
324 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
325 grid_re_kh, & !< Grid Reynolds number for Laplacian horizontal viscosity at h points [nondim]
326 grid_re_ah, & !< Grid Reynolds number for Biharmonic horizontal viscosity at h points [nondim]
336 real :: vert_vort_mag
346 real :: visc_bound_rem
355 real :: gme_coeff_limiter
362 real :: backscat_subround
367 logical :: rescale_kh, legacy_bound
368 logical :: find_frictwork
369 logical :: apply_obc = .false.
370 logical :: use_meke_ku
371 logical :: use_meke_au
372 integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz
373 integer :: i, j, k, n
374 real :: inv_pi3, inv_pi2, inv_pi6
375 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
376 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
378 h_neglect = gv%H_subroundoff
379 h_neglect3 = h_neglect**3
380 inv_pi3 = 1.0/((4.0*atan(1.0))**3)
381 inv_pi2 = 1.0/((4.0*atan(1.0))**2)
382 inv_pi6 = inv_pi3 * inv_pi3
387 if (
present(obc))
then ;
if (
associated(obc))
then ;
if (obc%OBC_pe)
then
388 apply_obc = obc%Flather_u_BCs_exist_globally .or. obc%Flather_v_BCs_exist_globally
390 endif ;
endif ;
endif
392 if (.not.
associated(cs))
call mom_error(fatal, &
393 "MOM_hor_visc: Module must be initialized before it is used.")
394 if (.not.(cs%Laplacian .or. cs%biharmonic))
return
396 find_frictwork = (cs%id_FrictWork > 0)
397 if (cs%id_FrictWorkIntz > 0) find_frictwork = .true.
398 if (
associated(meke))
then
399 if (
associated(meke%mom_src)) find_frictwork = .true.
400 backscat_subround = 0.0
401 if (find_frictwork .and.
associated(meke%mom_src) .and. (meke%backscatter_Ro_c > 0.0) .and. &
402 (meke%backscatter_Ro_Pow /= 0.0)) &
403 backscat_subround = (1.0e-16/meke%backscatter_Ro_c)**(1.0/meke%backscatter_Ro_Pow)
407 if (
associated(varmix))
then
408 rescale_kh = varmix%Resoln_scaled_Kh
409 if (rescale_kh .and. &
410 (.not.
associated(varmix%Res_fn_h) .or. .not.
associated(varmix%Res_fn_q))) &
411 call mom_error(fatal,
"MOM_hor_visc: VarMix%Res_fn_h and " //&
412 "VarMix%Res_fn_q both need to be associated with Resoln_scaled_Kh.")
414 legacy_bound = (cs%Smagorinsky_Kh .or. cs%Leith_Kh) .and. &
415 (cs%bound_Kh .and. .not.cs%better_bound_Kh)
418 use_meke_ku =
associated(meke%Ku)
419 use_meke_au =
associated(meke%Au)
422 do j=jsq-1,jeq+2 ;
do i=isq-1,ieq+2
423 boundary_mask_h(i,j) = (g%mask2dCu(i,j) * g%mask2dCv(i,j) * g%mask2dCu(i-1,j) * g%mask2dCv(i,j-1))
426 do j=js-2,jeq+1 ;
do i=is-2,ieq+1
427 boundary_mask_q(i,j) = (g%mask2dCv(i,j) * g%mask2dCv(i+1,j) * g%mask2dCu(i,j) * g%mask2dCu(i,j-1))
431 gme_coeff_h(:,:,:) = 0.0
432 gme_coeff_q(:,:,:) = 0.0
433 str_xx_gme(:,:) = 0.0
434 str_xy_gme(:,:) = 0.0
437 call barotropic_get_tav(bt, ubtav, vbtav, g, us)
440 do j=js-1,je+2 ;
do i=is-1,ie+2
441 dudx_bt(i,j) = cs%DY_dxT(i,j)*(g%IdyCu(i,j) * ubtav(i,j) - &
442 g%IdyCu(i-1,j) * ubtav(i-1,j))
443 dvdy_bt(i,j) = cs%DX_dyT(i,j)*(g%IdxCv(i,j) * vbtav(i,j) - &
444 g%IdxCv(i,j-1) * vbtav(i,j-1))
447 do j=jsq,jeq+1 ;
do i=isq,ieq+1
448 sh_xx_bt(i,j) = dudx_bt(i,j) - dvdy_bt(i,j)
452 do j=js-2,jeq+1 ;
do i=is-2,ieq+1
453 dvdx_bt(i,j) = cs%DY_dxBu(i,j)*(vbtav(i+1,j)*g%IdyCv(i+1,j) &
454 - vbtav(i,j)*g%IdyCv(i,j))
455 dudy_bt(i,j) = cs%DX_dyBu(i,j)*(ubtav(i,j+1)*g%IdxCu(i,j+1) &
456 - ubtav(i,j)*g%IdxCu(i,j))
459 call pass_vector(dudx_bt, dvdy_bt, g%Domain, stagger=bgrid_ne)
460 call pass_vector(dvdx_bt, dudy_bt, g%Domain, stagger=agrid)
463 do j=js-1,jeq ;
do i=is-1,ieq
464 sh_xy_bt(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx_bt(i,j) + dudy_bt(i,j) )
467 do j=js-1,jeq ;
do i=is-1,ieq
468 sh_xy_bt(i,j) = g%mask2dBu(i,j) * ( dvdx_bt(i,j) + dudy_bt(i,j) )
472 do j=jsq-1,jeq+2 ;
do i=isq-1,ieq+2
473 grad_vel_mag_bt_h(i,j) = boundary_mask_h(i,j) * (dudx_bt(i,j)**2 + dvdy_bt(i,j)**2 + &
474 (0.25*((dvdx_bt(i,j)+dvdx_bt(i-1,j-1))+(dvdx_bt(i,j-1)+dvdx_bt(i-1,j))))**2 + &
475 (0.25*((dudy_bt(i,j)+dudy_bt(i-1,j-1))+(dudy_bt(i,j-1)+dudy_bt(i-1,j))))**2)
478 do j=js-2,jeq+1 ;
do i=is-2,ieq+1
479 grad_vel_mag_bt_q(i,j) = boundary_mask_q(i,j) * (dvdx_bt(i,j)**2 + dudy_bt(i,j)**2 + &
480 (0.25*((dudx_bt(i,j)+dudx_bt(i+1,j+1))+(dudx_bt(i,j+1)+dudx_bt(i+1,j))))**2 + &
481 (0.25*((dvdy_bt(i,j)+dvdy_bt(i+1,j+1))+(dvdy_bt(i,j+1)+dvdy_bt(i+1,j))))**2)
521 do j=jsq-1,jeq+2 ;
do i=isq-1,ieq+2
522 dudx(i,j) = cs%DY_dxT(i,j)*(g%IdyCu(i,j) * u(i,j,k) - &
523 g%IdyCu(i-1,j) * u(i-1,j,k))
524 dvdy(i,j) = cs%DX_dyT(i,j)*(g%IdxCv(i,j) * v(i,j,k) - &
525 g%IdxCv(i,j-1) * v(i,j-1,k))
526 sh_xx(i,j) = dudx(i,j) - dvdy(i,j)
530 do j=jsq-2,jeq+2 ;
do i=isq-2,ieq+2
531 dvdx(i,j) = cs%DY_dxBu(i,j)*(v(i+1,j,k)*g%IdyCv(i+1,j) - v(i,j,k)*g%IdyCv(i,j))
532 dudy(i,j) = cs%DX_dyBu(i,j)*(u(i,j+1,k)*g%IdxCu(i,j+1) - u(i,j,k)*g%IdxCu(i,j))
540 if (cs%use_land_mask)
then
541 do j=js-2,je+2 ;
do i=isq-1,ieq+1
542 h_u(i,j) = 0.5 * (g%mask2dT(i,j)*h(i,j,k) + g%mask2dT(i+1,j)*h(i+1,j,k))
544 do j=jsq-1,jeq+1 ;
do i=is-2,ie+2
545 h_v(i,j) = 0.5 * (g%mask2dT(i,j)*h(i,j,k) + g%mask2dT(i,j+1)*h(i,j+1,k))
548 do j=js-2,je+2 ;
do i=isq-1,ieq+1
549 h_u(i,j) = 0.5 * (h(i,j,k) + h(i+1,j,k))
551 do j=jsq-1,jeq+1 ;
do i=is-2,ie+2
552 h_v(i,j) = 0.5 * (h(i,j,k) + h(i,j+1,k))
558 if (apply_obc)
then ;
do n=1,obc%number_of_segments
559 j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
560 if (obc%zero_strain .or. obc%freeslip_strain .or. obc%computed_strain)
then
561 if (obc%segment(n)%is_N_or_S .and. (j >= js-2) .and. (j <= jeq+1))
then
562 do i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
563 if (obc%zero_strain)
then
564 dvdx(i,j) = 0. ; dudy(i,j) = 0.
565 elseif (obc%freeslip_strain)
then
567 elseif (obc%computed_strain)
then
568 if (obc%segment(n)%direction == obc_direction_n)
then
569 dudy(i,j) = 2.0*cs%DX_dyBu(i,j)* &
570 (obc%segment(n)%tangential_vel(i,j,k) - u(i,j,k))*g%IdxCu(i,j)
572 dudy(i,j) = 2.0*cs%DX_dyBu(i,j)* &
573 (u(i,j+1,k) - obc%segment(n)%tangential_vel(i,j,k))*g%IdxCu(i,j+1)
575 elseif (obc%specified_strain)
then
576 if (obc%segment(n)%direction == obc_direction_n)
then
577 dudy(i,j) = cs%DX_dyBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdxCu(i,j)*g%dxBu(i,j)
579 dudy(i,j) = cs%DX_dyBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdxCu(i,j+1)*g%dxBu(i,j)
583 elseif (obc%segment(n)%is_E_or_W .and. (i >= is-2) .and. (i <= ieq+1))
then
584 do j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
585 if (obc%zero_strain)
then
586 dvdx(i,j) = 0. ; dudy(i,j) = 0.
587 elseif (obc%freeslip_strain)
then
589 elseif (obc%computed_strain)
then
590 if (obc%segment(n)%direction == obc_direction_e)
then
591 dvdx(i,j) = 2.0*cs%DY_dxBu(i,j)* &
592 (obc%segment(n)%tangential_vel(i,j,k) - v(i,j,k))*g%IdyCv(i,j)
594 dvdx(i,j) = 2.0*cs%DY_dxBu(i,j)* &
595 (v(i+1,j,k) - obc%segment(n)%tangential_vel(i,j,k))*g%IdyCv(i+1,j)
597 elseif (obc%specified_strain)
then
598 if (obc%segment(n)%direction == obc_direction_e)
then
599 dvdx(i,j) = cs%DY_dxBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdyCv(i,j)*g%dxBu(i,j)
601 dvdx(i,j) = cs%DY_dxBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdyCv(i+1,j)*g%dxBu(i,j)
609 if (obc%segment(n)%direction == obc_direction_n)
then
614 if ((j >= jsq-1) .and. (j <= jeq+1))
then
615 do i = max(is-2,obc%segment(n)%HI%isd), min(ie+2,obc%segment(n)%HI%ied)
619 elseif (obc%segment(n)%direction == obc_direction_s)
then
620 if ((j >= jsq-1) .and. (j <= jeq+1))
then
621 do i = max(is-2,obc%segment(n)%HI%isd), min(ie+2,obc%segment(n)%HI%ied)
622 h_v(i,j) = h(i,j+1,k)
625 elseif (obc%segment(n)%direction == obc_direction_e)
then
626 if ((i >= isq-1) .and. (i <= ieq+1))
then
627 do j = max(js-2,obc%segment(n)%HI%jsd), min(je+2,obc%segment(n)%HI%jed)
631 elseif (obc%segment(n)%direction == obc_direction_w)
then
632 if ((i >= isq-1) .and. (i <= ieq+1))
then
633 do j = max(js-2,obc%segment(n)%HI%jsd), min(je+2,obc%segment(n)%HI%jed)
634 h_u(i,j) = h(i+1,j,k)
640 if (apply_obc)
then ;
do n=1,obc%number_of_segments
641 j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
642 if (obc%segment(n)%direction == obc_direction_n)
then
643 if ((j >= js-2) .and. (j <= je))
then
644 do i = max(isq-1,obc%segment(n)%HI%IsdB), min(ieq+1,obc%segment(n)%HI%IedB)
645 h_u(i,j+1) = h_u(i,j)
648 elseif (obc%segment(n)%direction == obc_direction_s)
then
649 if ((j >= js-1) .and. (j <= je+1))
then
650 do i = max(isq-1,obc%segment(n)%HI%isd), min(ieq+1,obc%segment(n)%HI%ied)
651 h_u(i,j) = h_u(i,j+1)
654 elseif (obc%segment(n)%direction == obc_direction_e)
then
655 if ((i >= is-2) .and. (i <= ie))
then
656 do j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+1,obc%segment(n)%HI%jed)
657 h_v(i+1,j) = h_v(i,j)
660 elseif (obc%segment(n)%direction == obc_direction_w)
then
661 if ((i >= is-1) .and. (i <= ie+1))
then
662 do j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+1,obc%segment(n)%HI%jed)
663 h_v(i,j) = h_v(i+1,j)
672 do j=js-2,jeq+1 ;
do i=is-2,ieq+1
673 sh_xy(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx(i,j) + dudy(i,j) )
676 do j=js-2,jeq+1 ;
do i=is-2,ieq+1
677 sh_xy(i,j) = g%mask2dBu(i,j) * ( dvdx(i,j) + dudy(i,j) )
682 if (cs%biharmonic)
then
683 do j=js-1,jeq+1 ;
do i=isq-1,ieq+1
684 del2u(i,j) = cs%Idxdy2u(i,j)*(cs%dy2h(i+1,j)*sh_xx(i+1,j) - cs%dy2h(i,j)*sh_xx(i,j)) + &
685 cs%Idx2dyCu(i,j)*(cs%dx2q(i,j)*sh_xy(i,j) - cs%dx2q(i,j-1)*sh_xy(i,j-1))
687 do j=jsq-1,jeq+1 ;
do i=is-1,ieq+1
688 del2v(i,j) = cs%Idxdy2v(i,j)*(cs%dy2q(i,j)*sh_xy(i,j) - cs%dy2q(i-1,j)*sh_xy(i-1,j)) - &
689 cs%Idx2dyCv(i,j)*(cs%dx2h(i,j+1)*sh_xx(i,j+1) - cs%dx2h(i,j)*sh_xx(i,j))
691 if (apply_obc) then;
if (obc%zero_biharmonic)
then
692 do n=1,obc%number_of_segments
693 i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
694 if (obc%segment(n)%is_N_or_S .and. (j >= jsq-1) .and. (j <= jeq+1))
then
695 do i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
698 elseif (obc%segment(n)%is_E_or_W .and. (i >= isq-1) .and. (i <= ieq+1))
then
699 do j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
709 do j=jsq-2,jeq+2 ;
do i=isq-2,ieq+2
710 vort_xy(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx(i,j) - dudy(i,j) )
713 do j=jsq-2,jeq+2 ;
do i=isq-2,ieq+2
714 vort_xy(i,j) = g%mask2dBu(i,j) * ( dvdx(i,j) - dudy(i,j) )
719 do j=jsq-1,jeq+2 ;
do i=isq-1,ieq+2
720 div_xx(i,j) = dudx(i,j) + dvdy(i,j)
723 if ((cs%Leith_Kh) .or. (cs%Leith_Ah))
then
726 do j=jsq-1,jeq+1 ;
do i=isq-1,ieq+2
727 dy_dxbu = g%dyBu(i,j) * g%IdxBu(i,j)
728 vort_xy_dx(i,j) = dy_dxbu * (vort_xy(i,j) * g%IdyCu(i,j) - vort_xy(i-1,j) * g%IdyCu(i-1,j))
731 do j=jsq-1,jeq+2 ;
do i=isq-1,ieq+1
732 dx_dybu = g%dxBu(i,j) * g%IdyBu(i,j)
733 vort_xy_dy(i,j) = dx_dybu * (vort_xy(i,j) * g%IdxCv(i,j) - vort_xy(i,j-1) * g%IdxCv(i,j-1))
737 do j=jsq-1,jeq+1 ;
do i=isq-1,ieq+1
738 dy_dxbu = g%dyBu(i,j) * g%IdxBu(i,j)
739 dx_dybu = g%dxBu(i,j) * g%IdyBu(i,j)
741 del2vort_q(i,j) = dy_dxbu * (vort_xy_dx(i+1,j) * g%IdyCv(i+1,j) - vort_xy_dx(i,j) * g%IdyCv(i,j)) + &
742 dx_dybu * (vort_xy_dy(i,j+1) * g%IdyCu(i,j+1) - vort_xy_dy(i,j) * g%IdyCu(i,j))
744 do j=jsq,jeq+1 ;
do i=isq,ieq+1
745 del2vort_h(i,j) = 0.25*(del2vort_q(i,j) + del2vort_q(i-1,j) + del2vort_q(i,j-1) + del2vort_q(i-1,j-1))
748 if (cs%modified_Leith)
then
751 do j=jsq-1,jeq+2 ;
do i=isq-1,ieq+1
752 div_xx_dx(i,j) = g%IdxCu(i,j)*(div_xx(i+1,j) - div_xx(i,j))
754 do j=jsq-1,jeq+1 ;
do i=isq-1,ieq+2
755 div_xx_dy(i,j) = g%IdyCv(i,j)*(div_xx(i,j+1) - div_xx(i,j))
759 do j=jsq,jeq+1 ;
do i=isq,ieq+1
760 grad_div_mag_h(i,j) =sqrt((0.5*(div_xx_dx(i,j) + div_xx_dx(i-1,j)))**2 + &
761 (0.5 * (div_xx_dy(i,j) + div_xx_dy(i,j-1)))**2)
764 do j=jsq-1,jeq+1 ;
do i=isq-1,ieq+1
765 grad_div_mag_q(i,j) =sqrt((0.5*(div_xx_dx(i,j) + div_xx_dx(i,j+1)))**2 + &
766 (0.5 * (div_xx_dy(i,j) + div_xx_dy(i+1,j)))**2)
771 do j=jsq-1,jeq+2 ;
do i=is-2,ieq+1
774 do j=jsq-1,jeq+1 ;
do i=isq-1,ieq+2
777 do j=jsq,jeq+1 ;
do i=isq,ieq+1
778 grad_div_mag_h(i,j) = 0.0
780 do j=jsq-1,jeq+1 ;
do i=isq-1,ieq+1
781 grad_div_mag_q(i,j) = 0.0
787 if (cs%use_beta_in_Leith)
then
788 do j=js-2,jeq+1 ;
do i=is-1,ieq+1
789 vort_xy_dx(i,j) = vort_xy_dx(i,j) + 0.5 * ( g%dF_dx(i,j) + g%dF_dx(i,j+1))
791 do j=js-1,jeq+1 ;
do i=is-2,ieq+1
792 vort_xy_dy(i,j) = vort_xy_dy(i,j) + 0.5 * ( g%dF_dy(i,j) + g%dF_dy(i+1,j))
796 if (cs%use_QG_Leith_visc)
then
798 do j=jsq,jeq+1 ;
do i=isq,ieq+1
799 grad_vort_mag_h_2d(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i,j-1)))**2 + &
800 (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i-1,j)))**2 )
802 do j=js-1,jeq ;
do i=is-1,ieq
803 grad_vort_mag_q_2d(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i+1,j)))**2 + &
804 (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i,j+1)))**2 )
808 call calc_qg_leith_viscosity(varmix, g, gv, us, h, k, div_xx_dx, div_xx_dy, &
809 vort_xy_dx, vort_xy_dy)
813 do j=jsq,jeq+1 ;
do i=isq,ieq+1
814 grad_vort_mag_h(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i,j-1)))**2 + &
815 (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i-1,j)))**2 )
817 do j=js-1,jeq ;
do i=is-1,ieq
818 grad_vort_mag_q(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i+1,j)))**2 + &
819 (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i,j+1)))**2 )
826 do j=jsq,jeq+1 ;
do i=isq,ieq+1
827 if ((cs%Smagorinsky_Kh) .or. (cs%Smagorinsky_Ah))
then
828 shear_mag = sqrt(sh_xx(i,j)*sh_xx(i,j) + &
829 0.25*((sh_xy(i-1,j-1)*sh_xy(i-1,j-1) + sh_xy(i,j)*sh_xy(i,j)) + &
830 (sh_xy(i-1,j)*sh_xy(i-1,j) + sh_xy(i,j-1)*sh_xy(i,j-1))))
832 if ((cs%Leith_Kh) .or. (cs%Leith_Ah))
then
833 if (cs%use_QG_Leith_visc)
then
834 vert_vort_mag = min(grad_vort_mag_h(i,j) + grad_div_mag_h(i,j),3.*grad_vort_mag_h_2d(i,j))
836 vert_vort_mag = (grad_vort_mag_h(i,j) + grad_div_mag_h(i,j))
839 if (cs%better_bound_Ah .or. cs%better_bound_Kh)
then
840 hrat_min = min(1.0, min(h_u(i,j), h_u(i-1,j), h_v(i,j), h_v(i,j-1)) / &
841 (h(i,j,k) + h_neglect) )
845 if (cs%Laplacian)
then
848 kh = cs%Kh_bg_xx(i,j)
849 if (cs%add_LES_viscosity)
then
850 if (cs%Smagorinsky_Kh) kh = kh + cs%Laplac2_const_xx(i,j) * shear_mag
851 if (cs%Leith_Kh) kh = kh + cs%Laplac3_const_xx(i,j) * vert_vort_mag*inv_pi3
853 if (cs%Smagorinsky_Kh) kh = max( kh, cs%Laplac2_const_xx(i,j) * shear_mag )
854 if (cs%Leith_Kh) kh = max( kh, cs%Laplac3_const_xx(i,j) * vert_vort_mag*inv_pi3)
857 if (rescale_kh) kh = varmix%Res_fn_h(i,j) * kh
858 if (cs%res_scale_MEKE) meke_res_fn = varmix%Res_fn_h(i,j)
860 if (legacy_bound) kh = min(kh, cs%Kh_Max_xx(i,j))
861 kh = max( kh, cs%Kh_bg_min )
863 kh = kh + meke%Ku(i,j) * meke_res_fn
864 if (cs%anisotropic) kh = kh + cs%Kh_aniso * ( 1. - cs%n1n2_h(i,j)**2 )
868 if (cs%better_bound_Kh)
then
869 if (kh >= hrat_min*cs%Kh_Max_xx(i,j))
then
871 kh = hrat_min*cs%Kh_Max_xx(i,j)
873 visc_bound_rem = 1.0 - kh / (hrat_min*cs%Kh_Max_xx(i,j))
877 if ((cs%id_Kh_h>0) .or. find_frictwork .or. cs%debug) kh_h(i,j,k) = kh
879 if (cs%id_grid_Re_Kh>0)
then
880 ke = 0.125*((u(i,j,k)+u(i-1,j,k))**2 + (v(i,j,k)+v(i,j-1,k))**2)
881 grid_re_kh(i,j,k) = (sqrt(ke) * sqrt(cs%grid_sp_h2(i,j))) &
882 / max(kh, cs%min_grid_Kh)
885 if (cs%id_div_xx_h>0) div_xx_h(i,j,k) = div_xx(i,j)
886 if (cs%id_sh_xx_h>0) sh_xx_h(i,j,k) = sh_xx(i,j)
888 str_xx(i,j) = -kh * sh_xx(i,j)
893 if (cs%anisotropic)
then
895 local_strain = 0.25 * ( (sh_xy(i,j) + sh_xy(i-1,j-1)) + (sh_xy(i-1,j) + sh_xy(i,j-1)) )
897 str_xx(i,j) = str_xx(i,j) - cs%Kh_aniso * cs%n1n2_h(i,j) * cs%n1n1_m_n2n2_h(i,j) * local_strain
900 if (cs%biharmonic)
then
903 ahsm = 0.0; ahlth = 0.0
904 if ((cs%Smagorinsky_Ah) .or. (cs%Leith_Ah))
then
905 if (cs%Smagorinsky_Ah)
then
906 if (cs%bound_Coriolis)
then
907 ahsm = shear_mag * (cs%Biharm_const_xx(i,j) + &
908 cs%Biharm_const2_xx(i,j)*shear_mag)
910 ahsm = cs%Biharm_const_xx(i,j) * shear_mag
913 if (cs%Leith_Ah) ahlth = cs%Biharm6_const_xx(i,j) * abs(del2vort_h(i,j)) * inv_pi6
914 ah = max(max(cs%Ah_bg_xx(i,j), ahsm), ahlth)
915 if (cs%bound_Ah .and. .not.cs%better_bound_Ah) &
916 ah = min(ah, cs%Ah_Max_xx(i,j))
918 ah = cs%Ah_bg_xx(i,j)
921 if (use_meke_au) ah = ah + meke%Au(i,j)
923 if (cs%Re_Ah > 0.0)
then
924 ke = 0.125*((u(i,j,k)+u(i-1,j,k))**2 + (v(i,j,k)+v(i,j-1,k))**2)
925 ah = sqrt(ke) * cs%Re_Ah_const_xx(i,j)
928 if (cs%better_bound_Ah)
then
929 ah = min(ah, visc_bound_rem*hrat_min*cs%Ah_Max_xx(i,j))
932 if ((cs%id_Ah_h>0) .or. find_frictwork .or. cs%debug) ah_h(i,j,k) = ah
934 if (cs%id_grid_Re_Ah>0)
then
935 ke = 0.125*((u(i,j,k)+u(i-1,j,k))**2 + (v(i,j,k)+v(i,j-1,k))**2)
936 grid_re_ah(i,j,k) = (sqrt(ke) * cs%grid_sp_h3(i,j)) &
937 / max(ah, cs%min_grid_Ah)
940 str_xx(i,j) = str_xx(i,j) + ah * &
941 (cs%DY_dxT(i,j) * (g%IdyCu(i,j)*del2u(i,j) - g%IdyCu(i-1,j)*del2u(i-1,j)) - &
942 cs%DX_dyT(i,j) * (g%IdxCv(i,j)*del2v(i,j) - g%IdxCv(i,j-1)*del2v(i,j-1)))
945 bhstr_xx(i,j) = ah * &
946 (cs%DY_dxT(i,j) * (g%IdyCu(i,j)*del2u(i,j) - g%IdyCu(i-1,j)*del2u(i-1,j)) - &
947 cs%DX_dyT(i,j) * (g%IdxCv(i,j)*del2v(i,j) - g%IdxCv(i,j-1)*del2v(i,j-1)))
948 bhstr_xx(i,j) = bhstr_xx(i,j) * (h(i,j,k) * cs%reduction_xx(i,j))
954 if (cs%biharmonic)
then
956 do j=js-1,jeq ;
do i=is-1,ieq
957 ddel2vdx(i,j) = cs%DY_dxBu(i,j)*(del2v(i+1,j)*g%IdyCv(i+1,j) - del2v(i,j)*g%IdyCv(i,j))
958 ddel2udy(i,j) = cs%DX_dyBu(i,j)*(del2u(i,j+1)*g%IdxCu(i,j+1) - del2u(i,j)*g%IdxCu(i,j))
961 if (apply_obc)
then ;
if (obc%zero_strain .or. obc%freeslip_strain)
then
962 do n=1,obc%number_of_segments
963 j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
964 if (obc%segment(n)%is_N_or_S .and. (j >= js-1) .and. (j <= jeq))
then
965 do i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
966 if (obc%zero_strain)
then
967 ddel2vdx(i,j) = 0. ; ddel2udy(i,j) = 0.
968 elseif (obc%freeslip_strain)
then
972 elseif (obc%segment(n)%is_E_or_W .and. (i >= is-1) .and. (i <= ieq))
then
973 do j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
974 if (obc%zero_strain)
then
975 ddel2vdx(i,j) = 0. ; ddel2udy(i,j) = 0.
976 elseif (obc%freeslip_strain)
then
987 do j=js-1,jeq ;
do i=is-1,ieq
988 if ((cs%Smagorinsky_Kh) .or. (cs%Smagorinsky_Ah))
then
989 shear_mag = sqrt(sh_xy(i,j)*sh_xy(i,j) + &
990 0.25*((sh_xx(i,j)*sh_xx(i,j) + sh_xx(i+1,j+1)*sh_xx(i+1,j+1)) + &
991 (sh_xx(i,j+1)*sh_xx(i,j+1) + sh_xx(i+1,j)*sh_xx(i+1,j))))
993 if ((cs%Leith_Kh) .or. (cs%Leith_Ah))
then
994 if (cs%use_QG_Leith_visc)
then
995 vert_vort_mag = min(grad_vort_mag_q(i,j) + grad_div_mag_q(i,j), 3.*grad_vort_mag_q_2d(i,j))
997 vert_vort_mag = (grad_vort_mag_q(i,j) + grad_div_mag_q(i,j))
1000 h2uq = 4.0 * h_u(i,j) * h_u(i,j+1)
1001 h2vq = 4.0 * h_v(i,j) * h_v(i+1,j)
1002 hq(i,j) = 2.0 * h2uq * h2vq / (h_neglect3 + (h2uq + h2vq) * &
1003 ((h_u(i,j) + h_u(i,j+1)) + (h_v(i,j) + h_v(i+1,j))))
1005 if (cs%better_bound_Ah .or. cs%better_bound_Kh)
then
1006 hrat_min = min(1.0, min(h_u(i,j), h_u(i,j+1), h_v(i,j), h_v(i+1,j)) / &
1007 (hq(i,j) + h_neglect) )
1008 visc_bound_rem = 1.0
1011 if (cs%no_slip .and. (g%mask2dBu(i,j) < 0.5))
then
1012 if ((g%mask2dCu(i,j) + g%mask2dCu(i,j+1)) + &
1013 (g%mask2dCv(i,j) + g%mask2dCv(i+1,j)) > 0.0)
then
1016 hu = g%mask2dCu(i,j) * h_u(i,j) + g%mask2dCu(i,j+1) * h_u(i,j+1)
1017 hv = g%mask2dCv(i,j) * h_v(i,j) + g%mask2dCv(i+1,j) * h_v(i+1,j)
1018 if ((g%mask2dCu(i,j) + g%mask2dCu(i,j+1)) * &
1019 (g%mask2dCv(i,j) + g%mask2dCv(i+1,j)) == 0.0)
then
1025 hq(i,j) = 2.0 * (hu * hv) / ((hu + hv) + h_neglect)
1026 hrat_min = min(1.0, min(hu, hv) / (hq(i,j) + h_neglect) )
1031 if (cs%Laplacian)
then
1034 kh = cs%Kh_bg_xy(i,j)
1035 if (cs%add_LES_viscosity)
then
1036 if (cs%Smagorinsky_Kh) kh = kh + cs%Laplac2_const_xx(i,j) * shear_mag
1037 if (cs%Leith_Kh) kh = kh + cs%Laplac3_const_xx(i,j) * vert_vort_mag*inv_pi3
1039 if (cs%Smagorinsky_Kh) kh = max( kh, cs%Laplac2_const_xy(i,j) * shear_mag )
1040 if (cs%Leith_Kh) kh = max( kh, cs%Laplac3_const_xy(i,j) * vert_vort_mag*inv_pi3)
1043 if (rescale_kh) kh = varmix%Res_fn_q(i,j) * kh
1044 if (cs%res_scale_MEKE) meke_res_fn = varmix%Res_fn_q(i,j)
1046 if (legacy_bound) kh = min(kh, cs%Kh_Max_xy(i,j))
1047 kh = max( kh, cs%Kh_bg_min )
1048 if (use_meke_ku)
then
1049 kh = kh + 0.25*( (meke%Ku(i,j) + meke%Ku(i+1,j+1)) + &
1050 (meke%Ku(i+1,j) + meke%Ku(i,j+1)) ) * meke_res_fn
1053 if (cs%anisotropic) kh = kh + cs%Kh_aniso * cs%n1n2_q(i,j)**2
1057 if (cs%better_bound_Kh)
then
1058 if (kh >= hrat_min*cs%Kh_Max_xy(i,j))
then
1059 visc_bound_rem = 0.0
1060 kh = hrat_min*cs%Kh_Max_xy(i,j)
1061 elseif (cs%Kh_Max_xy(i,j)>0.)
then
1062 visc_bound_rem = 1.0 - kh / (hrat_min*cs%Kh_Max_xy(i,j))
1066 if (cs%id_Kh_q>0 .or. cs%debug) kh_q(i,j,k) = kh
1067 if (cs%id_vort_xy_q>0) vort_xy_q(i,j,k) = vort_xy(i,j)
1068 if (cs%id_sh_xy_q>0) sh_xy_q(i,j,k) = sh_xy(i,j)
1070 str_xy(i,j) = -kh * sh_xy(i,j)
1075 if (cs%anisotropic)
then
1077 local_strain = 0.25 * ( (sh_xx(i,j) + sh_xx(i+1,j+1)) + (sh_xx(i+1,j) + sh_xx(i,j+1)) )
1079 str_xy(i,j) = str_xy(i,j) - cs%Kh_aniso * cs%n1n2_q(i,j) * cs%n1n1_m_n2n2_q(i,j) * local_strain
1082 if (cs%biharmonic)
then
1085 ahsm = 0.0 ; ahlth = 0.0
1086 if (cs%Smagorinsky_Ah .or. cs%Leith_Ah)
then
1087 if (cs%Smagorinsky_Ah)
then
1088 if (cs%bound_Coriolis)
then
1089 ahsm = shear_mag * (cs%Biharm_const_xy(i,j) + &
1090 cs%Biharm_const2_xy(i,j)*shear_mag)
1092 ahsm = cs%Biharm_const_xy(i,j) * shear_mag
1095 if (cs%Leith_Ah) ahlth = cs%Biharm6_const_xy(i,j) * abs(del2vort_q(i,j)) * inv_pi6
1096 ah = max(max(cs%Ah_bg_xy(i,j), ahsm), ahlth)
1097 if (cs%bound_Ah .and. .not.cs%better_bound_Ah) &
1098 ah = min(ah, cs%Ah_Max_xy(i,j))
1100 ah = cs%Ah_bg_xy(i,j)
1103 if (use_meke_au)
then
1104 ah = ah + 0.25*( (meke%Au(i,j) + meke%Au(i+1,j+1)) + &
1105 (meke%Au(i+1,j) + meke%Au(i,j+1)) )
1108 if (cs%Re_Ah > 0.0)
then
1109 ke = 0.125*((u(i,j,k)+u(i,j+1,k))**2 + (v(i,j,k)+v(i+1,j,k))**2)
1110 ah = sqrt(ke) * cs%Re_Ah_const_xy(i,j)
1113 if (cs%better_bound_Ah)
then
1114 ah = min(ah, visc_bound_rem*hrat_min*cs%Ah_Max_xy(i,j))
1117 if (cs%id_Ah_q>0 .or. cs%debug) ah_q(i,j,k) = ah
1119 str_xy(i,j) = str_xy(i,j) + ah * ( ddel2vdx(i,j) + ddel2udy(i,j) )
1122 bhstr_xy(i,j) = ah * ( ddel2vdx(i,j) + ddel2udy(i,j) ) * &
1123 (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1129 if (cs%use_GME)
then
1130 call thickness_diffuse_get_kh(td, kh_u_gme, kh_v_gme, g)
1133 do j=jsq,jeq+1 ;
do i=isq,ieq+1
1134 if (grad_vel_mag_bt_h(i,j)>0)
then
1135 gme_coeff = cs%GME_efficiency * (min(g%bathyT(i,j)/cs%GME_h0,1.0)**2) * &
1136 (0.25*(kh_u_gme(i,j,k)+kh_u_gme(i-1,j,k)+kh_v_gme(i,j,k)+kh_v_gme(i,j-1,k)))
1142 gme_coeff = gme_coeff * boundary_mask_h(i,j)
1143 gme_coeff = min(gme_coeff, cs%GME_limiter)
1145 if ((cs%id_GME_coeff_h>0) .or. find_frictwork) gme_coeff_h(i,j,k) = gme_coeff
1146 str_xx_gme(i,j) = gme_coeff * sh_xx_bt(i,j)
1150 do j=js-1,jeq ;
do i=is-1,ieq
1151 if (grad_vel_mag_bt_q(i,j)>0)
then
1152 gme_coeff = cs%GME_efficiency * (min(g%bathyT(i,j)/cs%GME_h0,1.0)**2) * &
1153 (0.25*(kh_u_gme(i,j,k)+kh_u_gme(i,j+1,k)+kh_v_gme(i,j,k)+kh_v_gme(i+1,j,k)))
1159 gme_coeff = gme_coeff * boundary_mask_q(i,j)
1160 gme_coeff = min(gme_coeff, cs%GME_limiter)
1162 if (cs%id_GME_coeff_q>0) gme_coeff_q(i,j,k) = gme_coeff
1163 str_xy_gme(i,j) = gme_coeff * sh_xy_bt(i,j)
1168 call smooth_gme(cs,g,gme_flux_h=str_xx_gme)
1169 call smooth_gme(cs,g,gme_flux_q=str_xy_gme)
1171 do j=jsq,jeq+1 ;
do i=isq,ieq+1
1172 str_xx(i,j) = (str_xx(i,j) + str_xx_gme(i,j)) * (h(i,j,k) * cs%reduction_xx(i,j))
1175 do j=js-1,jeq ;
do i=is-1,ieq
1177 if (cs%no_slip)
then
1178 str_xy(i,j) = (str_xy(i,j) + str_xy_gme(i,j)) * (hq(i,j) * cs%reduction_xy(i,j))
1180 str_xy(i,j) = (str_xy(i,j) + str_xy_gme(i,j)) * (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1184 if (
associated(meke%GME_snk))
then
1185 do j=js,je ;
do i=is,ie
1186 frictwork_gme(i,j,k) = gme_coeff_h(i,j,k) * h(i,j,k) * gv%H_to_kg_m2 * grad_vel_mag_bt_h(i,j)
1191 do j=jsq,jeq+1 ;
do i=isq,ieq+1
1192 str_xx(i,j) = str_xx(i,j) * (h(i,j,k) * cs%reduction_xx(i,j))
1195 do j=js-1,jeq ;
do i=is-1,ieq
1196 if (cs%no_slip)
then
1197 str_xy(i,j) = str_xy(i,j) * (hq(i,j) * cs%reduction_xy(i,j))
1199 str_xy(i,j) = str_xy(i,j) * (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1206 do j=js,je ;
do i=isq,ieq
1207 diffu(i,j,k) = ((g%IdyCu(i,j)*(cs%dy2h(i,j) *str_xx(i,j) - &
1208 cs%dy2h(i+1,j)*str_xx(i+1,j)) + &
1209 g%IdxCu(i,j)*(cs%dx2q(i,j-1)*str_xy(i,j-1) - &
1210 cs%dx2q(i,j) *str_xy(i,j))) * &
1211 g%IareaCu(i,j)) / (h_u(i,j) + h_neglect)
1217 do n=1,obc%number_of_segments
1218 if (obc%segment(n)%is_E_or_W)
then
1219 i = obc%segment(n)%HI%IsdB
1220 do j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
1228 do j=jsq,jeq ;
do i=is,ie
1229 diffv(i,j,k) = ((g%IdyCv(i,j)*(cs%dy2q(i-1,j)*str_xy(i-1,j) - &
1230 cs%dy2q(i,j) *str_xy(i,j)) - &
1231 g%IdxCv(i,j)*(cs%dx2h(i,j) *str_xx(i,j) - &
1232 cs%dx2h(i,j+1)*str_xx(i,j+1))) * &
1233 g%IareaCv(i,j)) / (h_v(i,j) + h_neglect)
1238 do n=1,obc%number_of_segments
1239 if (obc%segment(n)%is_N_or_S)
then
1240 j = obc%segment(n)%HI%JsdB
1241 do i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
1248 if (find_frictwork)
then ;
do j=js,je ;
do i=is,ie
1251 frictwork(i,j,k) = gv%H_to_RZ * ( &
1252 (str_xx(i,j)*(u(i,j,k)-u(i-1,j,k))*g%IdxT(i,j) &
1253 -str_xx(i,j)*(v(i,j,k)-v(i,j-1,k))*g%IdyT(i,j)) &
1254 +0.25*((str_xy(i,j)*( &
1255 (u(i,j+1,k)-u(i,j,k))*g%IdyBu(i,j) &
1256 +(v(i+1,j,k)-v(i,j,k))*g%IdxBu(i,j) ) &
1257 +str_xy(i-1,j-1)*( &
1258 (u(i-1,j,k)-u(i-1,j-1,k))*g%IdyBu(i-1,j-1) &
1259 +(v(i,j-1,k)-v(i-1,j-1,k))*g%IdxBu(i-1,j-1) )) &
1261 (u(i-1,j+1,k)-u(i-1,j,k))*g%IdyBu(i-1,j) &
1262 +(v(i,j,k)-v(i-1,j,k))*g%IdxBu(i-1,j) ) &
1264 (u(i,j,k)-u(i,j-1,k))*g%IdyBu(i,j-1) &
1265 +(v(i+1,j-1,k)-v(i,j-1,k))*g%IdxBu(i,j-1) )) ) )
1266 enddo ;
enddo ;
endif
1271 if (find_frictwork .and.
associated(meke))
then ;
if (
associated(meke%mom_src))
then
1273 do j=js,je ;
do i=is,ie
1274 meke%mom_src(i,j) = 0.
1276 if (
associated(meke%GME_snk))
then
1277 do j=js,je ;
do i=is,ie
1278 meke%GME_snk(i,j) = 0.
1282 if (meke%backscatter_Ro_c /= 0.)
then
1283 do j=js,je ;
do i=is,ie
1284 fath = 0.25*( (abs(g%CoriolisBu(i-1,j-1)) + abs(g%CoriolisBu(i,j))) + &
1285 (abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j-1))) )
1286 shear_mag = sqrt(sh_xx(i,j)*sh_xx(i,j) + &
1287 0.25*((sh_xy(i-1,j-1)*sh_xy(i-1,j-1) + sh_xy(i,j)*sh_xy(i,j)) + &
1288 (sh_xy(i-1,j)*sh_xy(i-1,j) + sh_xy(i,j-1)*sh_xy(i,j-1))))
1289 if (cs%answers_2018)
then
1290 fath = (us%s_to_T*fath)**meke%backscatter_Ro_pow
1293 shear_mag = ( ( (us%s_to_T*shear_mag)**meke%backscatter_Ro_pow ) + 1.e-30 ) &
1294 * meke%backscatter_Ro_c
1297 roscl = shear_mag / ( fath + shear_mag )
1299 if (fath <= backscat_subround*shear_mag)
then
1302 sh_f_pow = meke%backscatter_Ro_c * (shear_mag / fath)**meke%backscatter_Ro_pow
1303 roscl = sh_f_pow / (1.0 + sh_f_pow)
1308 meke%mom_src(i,j) = meke%mom_src(i,j) + gv%H_to_RZ * ( &
1309 ((str_xx(i,j)-roscl*bhstr_xx(i,j))*(u(i,j,k)-u(i-1,j,k))*g%IdxT(i,j) &
1310 -(str_xx(i,j)-roscl*bhstr_xx(i,j))*(v(i,j,k)-v(i,j-1,k))*g%IdyT(i,j)) &
1311 +0.25*(((str_xy(i,j)-roscl*bhstr_xy(i,j))*( &
1312 (u(i,j+1,k)-u(i,j,k))*g%IdyBu(i,j) &
1313 +(v(i+1,j,k)-v(i,j,k))*g%IdxBu(i,j) ) &
1314 +(str_xy(i-1,j-1)-roscl*bhstr_xy(i-1,j-1))*( &
1315 (u(i-1,j,k)-u(i-1,j-1,k))*g%IdyBu(i-1,j-1) &
1316 +(v(i,j-1,k)-v(i-1,j-1,k))*g%IdxBu(i-1,j-1) )) &
1317 +((str_xy(i-1,j)-roscl*bhstr_xy(i-1,j))*( &
1318 (u(i-1,j+1,k)-u(i-1,j,k))*g%IdyBu(i-1,j) &
1319 +(v(i,j,k)-v(i-1,j,k))*g%IdxBu(i-1,j) ) &
1320 +(str_xy(i,j-1)-roscl*bhstr_xy(i,j-1))*( &
1321 (u(i,j,k)-u(i,j-1,k))*g%IdyBu(i,j-1) &
1322 +(v(i+1,j-1,k)-v(i,j-1,k))*g%IdxBu(i,j-1) )) ) )
1326 do j=js,je ;
do i=is,ie
1327 meke%mom_src(i,j) = meke%mom_src(i,j) + frictwork(i,j,k)
1330 if (cs%use_GME .and.
associated(meke))
then ;
if (
associated(meke%GME_snk))
then
1331 do j=js,je ;
do i=is,ie
1332 meke%GME_snk(i,j) = meke%GME_snk(i,j) + frictwork_gme(i,j,k)
1341 if (cs%id_diffu>0)
call post_data(cs%id_diffu, diffu, cs%diag)
1342 if (cs%id_diffv>0)
call post_data(cs%id_diffv, diffv, cs%diag)
1343 if (cs%id_FrictWork>0)
call post_data(cs%id_FrictWork, frictwork, cs%diag)
1344 if (cs%id_FrictWork_GME>0)
call post_data(cs%id_FrictWork_GME, frictwork_gme, cs%diag)
1345 if (cs%id_Ah_h>0)
call post_data(cs%id_Ah_h, ah_h, cs%diag)
1346 if (cs%id_grid_Re_Ah>0)
call post_data(cs%id_grid_Re_Ah, grid_re_ah, cs%diag)
1347 if (cs%id_div_xx_h>0)
call post_data(cs%id_div_xx_h, div_xx_h, cs%diag)
1348 if (cs%id_vort_xy_q>0)
call post_data(cs%id_vort_xy_q, vort_xy_q, cs%diag)
1349 if (cs%id_sh_xx_h>0)
call post_data(cs%id_sh_xx_h, sh_xx_h, cs%diag)
1350 if (cs%id_sh_xy_q>0)
call post_data(cs%id_sh_xy_q, sh_xy_q, cs%diag)
1351 if (cs%id_Ah_q>0)
call post_data(cs%id_Ah_q, ah_q, cs%diag)
1352 if (cs%id_Kh_h>0)
call post_data(cs%id_Kh_h, kh_h, cs%diag)
1353 if (cs%id_grid_Re_Kh>0)
call post_data(cs%id_grid_Re_Kh, grid_re_kh, cs%diag)
1354 if (cs%id_Kh_q>0)
call post_data(cs%id_Kh_q, kh_q, cs%diag)
1355 if (cs%id_GME_coeff_h > 0)
call post_data(cs%id_GME_coeff_h, gme_coeff_h, cs%diag)
1356 if (cs%id_GME_coeff_q > 0)
call post_data(cs%id_GME_coeff_q, gme_coeff_q, cs%diag)
1359 if (cs%Laplacian)
then
1360 call hchksum(kh_h,
"Kh_h", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
1361 call bchksum(kh_q,
"Kh_q", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
1363 if (cs%biharmonic)
call hchksum(ah_h,
"Ah_h", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
1364 if (cs%biharmonic)
call bchksum(ah_q,
"Ah_q", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
1367 if (cs%id_FrictWorkIntz > 0)
then
1369 do i=is,ie ; frictworkintz(i,j) = frictwork(i,j,1) ;
enddo
1370 do k=2,nz ;
do i=is,ie
1371 frictworkintz(i,j) = frictworkintz(i,j) + frictwork(i,j,k)
1374 call post_data(cs%id_FrictWorkIntz, frictworkintz, cs%diag)
1393 if (
present(adp) .and. (cs%id_hf_diffu_2d > 0))
then
1394 allocate(hf_diffu_2d(g%IsdB:g%IedB,g%jsd:g%jed))
1395 hf_diffu_2d(:,:) = 0.0
1396 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1397 hf_diffu_2d(i,j) = hf_diffu_2d(i,j) + diffu(i,j,k) * adp%diag_hfrac_u(i,j,k)
1398 enddo ;
enddo ;
enddo
1399 call post_data(cs%id_hf_diffu_2d, hf_diffu_2d, cs%diag)
1400 deallocate(hf_diffu_2d)
1402 if (
present(adp) .and. (cs%id_hf_diffv_2d > 0))
then
1403 allocate(hf_diffv_2d(g%isd:g%ied,g%JsdB:g%JedB))
1404 hf_diffv_2d(:,:) = 0.0
1405 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1406 hf_diffv_2d(i,j) = hf_diffv_2d(i,j) + diffv(i,j,k) * adp%diag_hfrac_v(i,j,k)
1407 enddo ;
enddo ;
enddo
1408 call post_data(cs%id_hf_diffv_2d, hf_diffv_2d, cs%diag)
1409 deallocate(hf_diffv_2d)
1412 end subroutine horizontal_viscosity
1417 subroutine hor_visc_init(Time, G, US, param_file, diag, CS, MEKE, ADp)
1418 type(time_type),
intent(in) :: time
1423 type(
diag_ctrl),
target,
intent(inout) :: diag
1428 real,
dimension(SZIB_(G),SZJ_(G)) :: u0u, u0v
1429 real,
dimension(SZI_(G),SZJB_(G)) :: v0u, v0v
1436 real :: min_grid_sp_h2
1437 real :: min_grid_sp_h4
1442 real :: boundcorconst
1448 real :: kh_vel_scale
1449 real :: ah_vel_scale
1450 real :: ah_time_scale
1451 real :: smag_lap_const
1452 real :: smag_bi_const
1453 real :: leith_lap_const
1454 real :: leith_bi_const
1459 real :: bound_cor_vel
1463 real :: kh_pwr_of_sine
1464 logical :: bound_cor_def
1470 logical :: default_2018_answers
1471 character(len=64) :: inputdir, filename
1474 real :: aniso_grid_dir(2)
1475 integer :: aniso_mode
1476 integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz
1477 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
1480 #include "version_variable.h"
1481 character(len=40) :: mdl =
"MOM_hor_visc"
1482 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1483 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1484 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1485 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1486 if (
associated(cs))
then
1487 call mom_error(warning,
"hor_visc_init called with an associated "// &
1488 "control structure.")
1497 cs%bound_Kh = .false. ; cs%better_bound_Kh = .false. ; cs%Smagorinsky_Kh = .false. ; cs%Leith_Kh = .false.
1498 cs%bound_Ah = .false. ; cs%better_bound_Ah = .false. ; cs%Smagorinsky_Ah = .false. ; cs%Leith_Ah = .false.
1499 cs%use_QG_Leith_visc = .false.
1500 cs%bound_Coriolis = .false.
1501 cs%Modified_Leith = .false.
1502 cs%anisotropic = .false.
1503 cs%dynamic_aniso = .false.
1507 call get_param(param_file, mdl,
"GET_ALL_PARAMS", get_all, default=.false.)
1508 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
1509 "This sets the default value for the various _2018_ANSWERS parameters.", &
1511 call get_param(param_file, mdl,
"HOR_VISC_2018_ANSWERS", cs%answers_2018, &
1512 "If true, use the order of arithmetic and expressions that recover the "//&
1513 "answers from the end of 2018. Otherwise, use updated and more robust "//&
1514 "forms of the same expressions.", default=default_2018_answers)
1515 call get_param(param_file, mdl,
"DEBUG", cs%debug, default=.false.)
1516 call get_param(param_file, mdl,
"LAPLACIAN", cs%Laplacian, &
1517 "If true, use a Laplacian horizontal viscosity.", &
1519 if (cs%Laplacian .or. get_all)
then
1520 call get_param(param_file, mdl,
"KH", kh, &
1521 "The background Laplacian horizontal viscosity.", &
1522 units =
"m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1523 call get_param(param_file, mdl,
"KH_BG_MIN", cs%Kh_bg_min, &
1524 "The minimum value allowed for Laplacian horizontal viscosity, KH.", &
1525 units =
"m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1526 call get_param(param_file, mdl,
"KH_VEL_SCALE", kh_vel_scale, &
1527 "The velocity scale which is multiplied by the grid "//&
1528 "spacing to calculate the Laplacian viscosity. "//&
1529 "The final viscosity is the largest of this scaled "//&
1530 "viscosity, the Smagorinsky and Leith viscosities, and KH.", &
1531 units=
"m s-1", default=0.0, scale=us%m_s_to_L_T)
1532 call get_param(param_file, mdl,
"KH_SIN_LAT", kh_sin_lat, &
1533 "The amplitude of a latitudinally-dependent background "//&
1534 "viscosity of the form KH_SIN_LAT*(SIN(LAT)**KH_PWR_OF_SINE).", &
1535 units =
"m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1536 if (kh_sin_lat>0. .or. get_all) &
1537 call get_param(param_file, mdl,
"KH_PWR_OF_SINE", kh_pwr_of_sine, &
1538 "The power used to raise SIN(LAT) when using a latitudinally "//&
1539 "dependent background viscosity.", &
1540 units =
"nondim", default=4.0)
1541 call get_param(param_file, mdl,
"SMAGORINSKY_KH", cs%Smagorinsky_Kh, &
1542 "If true, use a Smagorinsky nonlinear eddy viscosity.", &
1544 if (cs%Smagorinsky_Kh .or. get_all) &
1545 call get_param(param_file, mdl,
"SMAG_LAP_CONST", smag_lap_const, &
1546 "The nondimensional Laplacian Smagorinsky constant, "//&
1547 "often 0.15.", units=
"nondim", default=0.0, &
1548 fail_if_missing = cs%Smagorinsky_Kh)
1549 call get_param(param_file, mdl,
"LEITH_KH", cs%Leith_Kh, &
1550 "If true, use a Leith nonlinear eddy viscosity.", &
1552 call get_param(param_file, mdl,
"MODIFIED_LEITH", cs%Modified_Leith, &
1553 "If true, add a term to Leith viscosity which is "//&
1554 "proportional to the gradient of divergence.", &
1556 call get_param(param_file, mdl,
"RES_SCALE_MEKE_VISC", cs%res_scale_MEKE, &
1557 "If true, the viscosity contribution from MEKE is scaled by "//&
1558 "the resolution function.", default=.false.)
1559 if (cs%Leith_Kh .or. get_all)
then
1560 call get_param(param_file, mdl,
"LEITH_LAP_CONST", leith_lap_const, &
1561 "The nondimensional Laplacian Leith constant, "//&
1562 "often set to 1.0", units=
"nondim", default=0.0, &
1563 fail_if_missing = cs%Leith_Kh)
1564 call get_param(param_file, mdl,
"USE_QG_LEITH_VISC", cs%use_QG_Leith_visc, &
1565 "If true, use QG Leith nonlinear eddy viscosity.", &
1567 if (cs%use_QG_Leith_visc .and. .not. cs%Leith_Kh)
call mom_error(fatal, &
1568 "MOM_lateral_mixing_coeffs.F90, VarMix_init:"//&
1569 "LEITH_KH must be True when USE_QG_LEITH_VISC=True.")
1571 if (cs%Leith_Kh .or. cs%Leith_Ah .or. get_all)
then
1572 call get_param(param_file, mdl,
"USE_BETA_IN_LEITH", cs%use_beta_in_Leith, &
1573 "If true, include the beta term in the Leith nonlinear eddy viscosity.", &
1574 default=cs%Leith_Kh)
1575 call get_param(param_file, mdl,
"MODIFIED_LEITH", cs%modified_Leith, &
1576 "If true, add a term to Leith viscosity which is \n"//&
1577 "proportional to the gradient of divergence.", &
1580 call get_param(param_file, mdl,
"BOUND_KH", cs%bound_Kh, &
1581 "If true, the Laplacian coefficient is locally limited "//&
1582 "to be stable.", default=.true.)
1583 call get_param(param_file, mdl,
"BETTER_BOUND_KH", cs%better_bound_Kh, &
1584 "If true, the Laplacian coefficient is locally limited "//&
1585 "to be stable with a better bounding than just BOUND_KH.", &
1586 default=cs%bound_Kh)
1587 call get_param(param_file, mdl,
"ANISOTROPIC_VISCOSITY", cs%anisotropic, &
1588 "If true, allow anistropic viscosity in the Laplacian "//&
1589 "horizontal viscosity.", default=.false.)
1590 call get_param(param_file, mdl,
"ADD_LES_VISCOSITY", cs%add_LES_viscosity, &
1591 "If true, adds the viscosity from Smagorinsky and Leith to the "//&
1592 "background viscosity instead of taking the maximum.", default=.false.)
1594 if (cs%anisotropic .or. get_all)
then
1595 call get_param(param_file, mdl,
"KH_ANISO", cs%Kh_aniso, &
1596 "The background Laplacian anisotropic horizontal viscosity.", &
1597 units =
"m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1598 call get_param(param_file, mdl,
"ANISOTROPIC_MODE", aniso_mode, &
1599 "Selects the mode for setting the direction of anistropy.\n"//&
1600 "\t 0 - Points along the grid i-direction.\n"//&
1601 "\t 1 - Points towards East.\n"//&
1602 "\t 2 - Points along the flow direction, U/|U|.", &
1604 select case (aniso_mode)
1606 call get_param(param_file, mdl,
"ANISO_GRID_DIR", aniso_grid_dir, &
1607 "The vector pointing in the direction of anistropy for "//&
1608 "horizont viscosity. n1,n2 are the i,j components relative "//&
1609 "to the grid.", units =
"nondim", fail_if_missing=.true.)
1611 call get_param(param_file, mdl,
"ANISO_GRID_DIR", aniso_grid_dir, &
1612 "The vector pointing in the direction of anistropy for "//&
1613 "horizont viscosity. n1,n2 are the i,j components relative "//&
1614 "to the spherical coordinates.", units =
"nondim", fail_if_missing=.true.)
1617 call get_param(param_file, mdl,
"BIHARMONIC", cs%biharmonic, &
1618 "If true, use a biharmonic horizontal viscosity. "//&
1619 "BIHARMONIC may be used with LAPLACIAN.", &
1621 if (cs%biharmonic .or. get_all)
then
1622 call get_param(param_file, mdl,
"AH", ah, &
1623 "The background biharmonic horizontal viscosity.", &
1624 units =
"m4 s-1", default=0.0, scale=us%m_to_L**4*us%T_to_s)
1625 call get_param(param_file, mdl,
"AH_VEL_SCALE", ah_vel_scale, &
1626 "The velocity scale which is multiplied by the cube of "//&
1627 "the grid spacing to calculate the biharmonic viscosity. "//&
1628 "The final viscosity is the largest of this scaled "//&
1629 "viscosity, the Smagorinsky and Leith viscosities, and AH.", &
1630 units=
"m s-1", default=0.0, scale=us%m_s_to_L_T)
1631 call get_param(param_file, mdl,
"AH_TIME_SCALE", ah_time_scale, &
1632 "A time scale whose inverse is multiplied by the fourth "//&
1633 "power of the grid spacing to calculate biharmonic viscosity. "//&
1634 "The final viscosity is the largest of all viscosity "//&
1635 "formulations in use. 0.0 means that it's not used.", &
1636 units=
"s", default=0.0, scale=us%s_to_T)
1637 call get_param(param_file, mdl,
"SMAGORINSKY_AH", cs%Smagorinsky_Ah, &
1638 "If true, use a biharmonic Smagorinsky nonlinear eddy "//&
1639 "viscosity.", default=.false.)
1640 call get_param(param_file, mdl,
"LEITH_AH", cs%Leith_Ah, &
1641 "If true, use a biharmonic Leith nonlinear eddy "//&
1642 "viscosity.", default=.false.)
1643 call get_param(param_file, mdl,
"BOUND_AH", cs%bound_Ah, &
1644 "If true, the biharmonic coefficient is locally limited "//&
1645 "to be stable.", default=.true.)
1646 call get_param(param_file, mdl,
"BETTER_BOUND_AH", cs%better_bound_Ah, &
1647 "If true, the biharmonic coefficient is locally limited "//&
1648 "to be stable with a better bounding than just BOUND_AH.", &
1649 default=cs%bound_Ah)
1650 call get_param(param_file, mdl,
"RE_AH", cs%Re_Ah, &
1651 "If nonzero, the biharmonic coefficient is scaled "//&
1652 "so that the biharmonic Reynolds number is equal to this.", &
1653 units=
"nondim", default=0.0)
1655 if (cs%Smagorinsky_Ah .or. get_all)
then
1656 call get_param(param_file, mdl,
"SMAG_BI_CONST",smag_bi_const, &
1657 "The nondimensional biharmonic Smagorinsky constant, "//&
1658 "typically 0.015 - 0.06.", units=
"nondim", default=0.0, &
1659 fail_if_missing = cs%Smagorinsky_Ah)
1660 call get_param(param_file, mdl,
"BOUND_CORIOLIS", bound_cor_def, default=.false.)
1661 call get_param(param_file, mdl,
"BOUND_CORIOLIS_BIHARM", cs%bound_Coriolis, &
1662 "If true use a viscosity that increases with the square "//&
1663 "of the velocity shears, so that the resulting viscous "//&
1664 "drag is of comparable magnitude to the Coriolis terms "//&
1665 "when the velocity differences between adjacent grid "//&
1666 "points is 0.5*BOUND_CORIOLIS_VEL. The default is the "//&
1667 "value of BOUND_CORIOLIS (or false).", default=bound_cor_def)
1668 if (cs%bound_Coriolis .or. get_all)
then
1669 call get_param(param_file, mdl,
"MAXVEL", maxvel, default=3.0e8)
1670 bound_cor_vel = maxvel
1671 call get_param(param_file, mdl,
"BOUND_CORIOLIS_VEL", bound_cor_vel, &
1672 "The velocity scale at which BOUND_CORIOLIS_BIHARM causes "//&
1673 "the biharmonic drag to have comparable magnitude to the "//&
1674 "Coriolis acceleration. The default is set by MAXVEL.", &
1675 units=
"m s-1", default=maxvel, scale=us%m_s_to_L_T)
1678 if (cs%Leith_Ah .or. get_all) &
1679 call get_param(param_file, mdl,
"LEITH_BI_CONST", leith_bi_const, &
1680 "The nondimensional biharmonic Leith constant, "//&
1681 "typical values are thus far undetermined.", units=
"nondim", default=0.0, &
1682 fail_if_missing = cs%Leith_Ah)
1684 call get_param(param_file, mdl,
"USE_LAND_MASK_FOR_HVISC", cs%use_land_mask, &
1685 "If true, use Use the land mask for the computation of thicknesses "//&
1686 "at velocity locations. This eliminates the dependence on arbitrary "//&
1687 "values over land or outside of the domain.", default=.true.)
1688 if (cs%better_bound_Ah .or. cs%better_bound_Kh .or. get_all) &
1689 call get_param(param_file, mdl,
"HORVISC_BOUND_COEF", cs%bound_coef, &
1690 "The nondimensional coefficient of the ratio of the "//&
1691 "viscosity bounds to the theoretical maximum for "//&
1692 "stability without considering other terms.", units=
"nondim", &
1694 call get_param(param_file, mdl,
"NOSLIP", cs%no_slip, &
1695 "If true, no slip boundary conditions are used; otherwise "//&
1696 "free slip boundary conditions are assumed. The "//&
1697 "implementation of the free slip BCs on a C-grid is much "//&
1698 "cleaner than the no slip BCs. The use of free slip BCs "//&
1699 "is strongly encouraged, and no slip BCs are not used with "//&
1700 "the biharmonic viscosity.", default=.false.)
1701 call get_param(param_file, mdl,
"USE_KH_BG_2D", cs%use_Kh_bg_2d, &
1702 "If true, read a file containing 2-d background harmonic "//&
1703 "viscosities. The final viscosity is the maximum of the other "//&
1704 "terms and this background value.", default=.false.)
1706 call get_param(param_file, mdl,
"USE_GME", cs%use_GME, &
1707 "If true, use the GM+E backscatter scheme in association \n"//&
1708 "with the Gent and McWilliams parameterization.", default=.false.)
1709 if (cs%use_GME)
then
1710 call get_param(param_file, mdl,
"SPLIT", split, &
1711 "Use the split time stepping if true.", default=.true., &
1713 if (.not. split)
call mom_error(fatal,
"ERROR: Currently, USE_GME = True "// &
1714 "cannot be used with SPLIT=False.")
1715 call get_param(param_file, mdl,
"GME_H0", cs%GME_h0, &
1716 "The strength of GME tapers quadratically to zero when the bathymetric "//&
1717 "depth is shallower than GME_H0.", units=
"m", scale=us%m_to_Z, &
1719 call get_param(param_file, mdl,
"GME_EFFICIENCY", cs%GME_efficiency, &
1720 "The nondimensional prefactor multiplying the GME coefficient.", &
1721 units=
"nondim", default=1.0)
1722 call get_param(param_file, mdl,
"GME_LIMITER", cs%GME_limiter, &
1723 "The absolute maximum value the GME coefficient is allowed to take.", &
1724 units=
"m2 s-1", scale=us%m_to_L**2*us%T_to_s, default=1.0e7)
1726 if (cs%Laplacian .or. cs%biharmonic)
then
1727 call get_param(param_file, mdl,
"DT", dt, &
1728 "The (baroclinic) dynamics time step.", units=
"s", scale=us%s_to_T, &
1729 fail_if_missing=.true.)
1732 if (cs%no_slip .and. cs%biharmonic) &
1733 call mom_error(fatal,
"ERROR: NOSLIP and BIHARMONIC cannot be defined "// &
1734 "at the same time in MOM.")
1735 if (.not.(cs%Laplacian .or. cs%biharmonic))
then
1737 if ( max(g%domain%niglobal, g%domain%njglobal)>2 )
call mom_error(warning, &
1738 "hor_visc_init: It is usually a very bad idea not to use either "//&
1739 "LAPLACIAN or BIHARMONIC viscosity.")
1742 deg2rad = atan(1.0) / 45.
1743 alloc_(cs%dx2h(isd:ied,jsd:jed)) ; cs%dx2h(:,:) = 0.0
1744 alloc_(cs%dy2h(isd:ied,jsd:jed)) ; cs%dy2h(:,:) = 0.0
1745 alloc_(cs%dx2q(isdb:iedb,jsdb:jedb)) ; cs%dx2q(:,:) = 0.0
1746 alloc_(cs%dy2q(isdb:iedb,jsdb:jedb)) ; cs%dy2q(:,:) = 0.0
1747 alloc_(cs%dx_dyT(isd:ied,jsd:jed)) ; cs%dx_dyT(:,:) = 0.0
1748 alloc_(cs%dy_dxT(isd:ied,jsd:jed)) ; cs%dy_dxT(:,:) = 0.0
1749 alloc_(cs%dx_dyBu(isdb:iedb,jsdb:jedb)) ; cs%dx_dyBu(:,:) = 0.0
1750 alloc_(cs%dy_dxBu(isdb:iedb,jsdb:jedb)) ; cs%dy_dxBu(:,:) = 0.0
1751 if (cs%Laplacian)
then
1752 alloc_(cs%grid_sp_h2(isd:ied,jsd:jed)) ; cs%grid_sp_h2(:,:) = 0.0
1753 alloc_(cs%Kh_bg_xx(isd:ied,jsd:jed)) ; cs%Kh_bg_xx(:,:) = 0.0
1754 alloc_(cs%Kh_bg_xy(isdb:iedb,jsdb:jedb)) ; cs%Kh_bg_xy(:,:) = 0.0
1755 if (cs%bound_Kh .or. cs%better_bound_Kh)
then
1756 alloc_(cs%Kh_Max_xx(isd:ied,jsd:jed)) ; cs%Kh_Max_xx(:,:) = 0.0
1757 alloc_(cs%Kh_Max_xy(isdb:iedb,jsdb:jedb)) ; cs%Kh_Max_xy(:,:) = 0.0
1759 if (cs%Smagorinsky_Kh)
then
1760 alloc_(cs%Laplac2_const_xx(isd:ied,jsd:jed)) ; cs%Laplac2_const_xx(:,:) = 0.0
1761 alloc_(cs%Laplac2_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Laplac2_const_xy(:,:) = 0.0
1763 if (cs%Leith_Kh)
then
1764 alloc_(cs%Laplac3_const_xx(isd:ied,jsd:jed)) ; cs%Laplac3_const_xx(:,:) = 0.0
1765 alloc_(cs%Laplac3_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Laplac3_const_xy(:,:) = 0.0
1768 alloc_(cs%reduction_xx(isd:ied,jsd:jed)) ; cs%reduction_xx(:,:) = 0.0
1769 alloc_(cs%reduction_xy(isdb:iedb,jsdb:jedb)) ; cs%reduction_xy(:,:) = 0.0
1770 if (cs%anisotropic)
then
1771 alloc_(cs%n1n2_h(isd:ied,jsd:jed)) ; cs%n1n2_h(:,:) = 0.0
1772 alloc_(cs%n1n1_m_n2n2_h(isd:ied,jsd:jed)) ; cs%n1n1_m_n2n2_h(:,:) = 0.0
1773 alloc_(cs%n1n2_q(isdb:iedb,jsdb:jedb)) ; cs%n1n2_q(:,:) = 0.0
1774 alloc_(cs%n1n1_m_n2n2_q(isdb:iedb,jsdb:jedb)) ; cs%n1n1_m_n2n2_q(:,:) = 0.0
1775 select case (aniso_mode)
1777 call align_aniso_tensor_to_grid(cs, aniso_grid_dir(1), aniso_grid_dir(2))
1781 cs%dynamic_aniso = .true.
1783 call mom_error(fatal,
"MOM_hor_visc: "//&
1784 "Runtime parameter ANISOTROPIC_MODE is out of range.")
1787 if (cs%use_Kh_bg_2d)
then
1788 alloc_(cs%Kh_bg_2d(isd:ied,jsd:jed)) ; cs%Kh_bg_2d(:,:) = 0.0
1789 call get_param(param_file, mdl,
"KH_BG_2D_FILENAME", filename, &
1790 'The filename containing a 2d map of "Kh".', &
1791 default=
'KH_background_2d.nc')
1792 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
1793 inputdir = slasher(inputdir)
1794 call mom_read_data(trim(inputdir)//trim(filename),
'Kh', cs%Kh_bg_2d, &
1795 g%domain, timelevel=1, scale=us%m_to_L**2*us%T_to_s)
1796 call pass_var(cs%Kh_bg_2d, g%domain)
1798 if (cs%biharmonic)
then
1799 alloc_(cs%Idx2dyCu(isdb:iedb,jsd:jed)) ; cs%Idx2dyCu(:,:) = 0.0
1800 alloc_(cs%Idx2dyCv(isd:ied,jsdb:jedb)) ; cs%Idx2dyCv(:,:) = 0.0
1801 alloc_(cs%Idxdy2u(isdb:iedb,jsd:jed)) ; cs%Idxdy2u(:,:) = 0.0
1802 alloc_(cs%Idxdy2v(isd:ied,jsdb:jedb)) ; cs%Idxdy2v(:,:) = 0.0
1803 alloc_(cs%Ah_bg_xx(isd:ied,jsd:jed)) ; cs%Ah_bg_xx(:,:) = 0.0
1804 alloc_(cs%Ah_bg_xy(isdb:iedb,jsdb:jedb)) ; cs%Ah_bg_xy(:,:) = 0.0
1805 alloc_(cs%grid_sp_h3(isd:ied,jsd:jed)) ; cs%grid_sp_h3(:,:) = 0.0
1806 if (cs%bound_Ah .or. cs%better_bound_Ah)
then
1807 alloc_(cs%Ah_Max_xx(isd:ied,jsd:jed)) ; cs%Ah_Max_xx(:,:) = 0.0
1808 alloc_(cs%Ah_Max_xy(isdb:iedb,jsdb:jedb)) ; cs%Ah_Max_xy(:,:) = 0.0
1810 if (cs%Smagorinsky_Ah)
then
1811 alloc_(cs%Biharm_const_xx(isd:ied,jsd:jed)) ; cs%Biharm_const_xx(:,:) = 0.0
1812 alloc_(cs%Biharm_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Biharm_const_xy(:,:) = 0.0
1813 if (cs%bound_Coriolis)
then
1814 alloc_(cs%Biharm_const2_xx(isd:ied,jsd:jed)) ; cs%Biharm_const2_xx(:,:) = 0.0
1815 alloc_(cs%Biharm_const2_xy(isdb:iedb,jsdb:jedb)) ; cs%Biharm_const2_xy(:,:) = 0.0
1818 if (cs%Leith_Ah)
then
1819 alloc_(cs%biharm6_const_xx(isd:ied,jsd:jed)) ; cs%biharm6_const_xx(:,:) = 0.0
1820 alloc_(cs%biharm6_const_xy(isdb:iedb,jsdb:jedb)) ; cs%biharm6_const_xy(:,:) = 0.0
1822 if (cs%Re_Ah > 0.0)
then
1823 alloc_(cs%Re_Ah_const_xx(isd:ied,jsd:jed)); cs%Re_Ah_const_xx(:,:) = 0.0
1824 alloc_(cs%Re_Ah_const_xy(isdb:iedb,jsdb:jedb)); cs%Re_Ah_const_xy(:,:) = 0.0
1827 do j=js-2,jeq+1 ;
do i=is-2,ieq+1
1828 cs%dx2q(i,j) = g%dxBu(i,j)*g%dxBu(i,j) ; cs%dy2q(i,j) = g%dyBu(i,j)*g%dyBu(i,j)
1829 cs%DX_dyBu(i,j) = g%dxBu(i,j)*g%IdyBu(i,j) ; cs%DY_dxBu(i,j) = g%dyBu(i,j)*g%IdxBu(i,j)
1831 do j=jsq-1,jeq+2 ;
do i=isq-1,ieq+2
1832 cs%dx2h(i,j) = g%dxT(i,j)*g%dxT(i,j) ; cs%dy2h(i,j) = g%dyT(i,j)*g%dyT(i,j)
1833 cs%DX_dyT(i,j) = g%dxT(i,j)*g%IdyT(i,j) ; cs%DY_dxT(i,j) = g%dyT(i,j)*g%IdxT(i,j)
1835 do j=jsq,jeq+1 ;
do i=isq,ieq+1
1836 cs%reduction_xx(i,j) = 1.0
1837 if ((g%dy_Cu(i,j) > 0.0) .and. (g%dy_Cu(i,j) < g%dyCu(i,j)) .and. &
1838 (g%dy_Cu(i,j) < g%dyCu(i,j) * cs%reduction_xx(i,j))) &
1839 cs%reduction_xx(i,j) = g%dy_Cu(i,j) / (g%dyCu(i,j))
1840 if ((g%dy_Cu(i-1,j) > 0.0) .and. (g%dy_Cu(i-1,j) < g%dyCu(i-1,j)) .and. &
1841 (g%dy_Cu(i-1,j) < g%dyCu(i-1,j) * cs%reduction_xx(i,j))) &
1842 cs%reduction_xx(i,j) = g%dy_Cu(i-1,j) / (g%dyCu(i-1,j))
1843 if ((g%dx_Cv(i,j) > 0.0) .and. (g%dx_Cv(i,j) < g%dxCv(i,j)) .and. &
1844 (g%dx_Cv(i,j) < g%dxCv(i,j) * cs%reduction_xx(i,j))) &
1845 cs%reduction_xx(i,j) = g%dx_Cv(i,j) / (g%dxCv(i,j))
1846 if ((g%dx_Cv(i,j-1) > 0.0) .and. (g%dx_Cv(i,j-1) < g%dxCv(i,j-1)) .and. &
1847 (g%dx_Cv(i,j-1) < g%dxCv(i,j-1) * cs%reduction_xx(i,j))) &
1848 cs%reduction_xx(i,j) = g%dx_Cv(i,j-1) / (g%dxCv(i,j-1))
1850 do j=js-1,jeq ;
do i=is-1,ieq
1851 cs%reduction_xy(i,j) = 1.0
1852 if ((g%dy_Cu(i,j) > 0.0) .and. (g%dy_Cu(i,j) < g%dyCu(i,j)) .and. &
1853 (g%dy_Cu(i,j) < g%dyCu(i,j) * cs%reduction_xy(i,j))) &
1854 cs%reduction_xy(i,j) = g%dy_Cu(i,j) / (g%dyCu(i,j))
1855 if ((g%dy_Cu(i,j+1) > 0.0) .and. (g%dy_Cu(i,j+1) < g%dyCu(i,j+1)) .and. &
1856 (g%dy_Cu(i,j+1) < g%dyCu(i,j+1) * cs%reduction_xy(i,j))) &
1857 cs%reduction_xy(i,j) = g%dy_Cu(i,j+1) / (g%dyCu(i,j+1))
1858 if ((g%dx_Cv(i,j) > 0.0) .and. (g%dx_Cv(i,j) < g%dxCv(i,j)) .and. &
1859 (g%dx_Cv(i,j) < g%dxCv(i,j) * cs%reduction_xy(i,j))) &
1860 cs%reduction_xy(i,j) = g%dx_Cv(i,j) / (g%dxCv(i,j))
1861 if ((g%dx_Cv(i+1,j) > 0.0) .and. (g%dx_Cv(i+1,j) < g%dxCv(i+1,j)) .and. &
1862 (g%dx_Cv(i+1,j) < g%dxCv(i+1,j) * cs%reduction_xy(i,j))) &
1863 cs%reduction_xy(i,j) = g%dx_Cv(i+1,j) / (g%dxCv(i+1,j))
1865 if (cs%Laplacian)
then
1868 if (cs%bound_Kh .or. cs%bound_Ah) kh_limit = 0.3 / (dt*4.0)
1871 min_grid_sp_h2 = huge(1.)
1872 do j=jsq,jeq+1 ;
do i=isq,ieq+1
1874 grid_sp_h2 = (2.0*cs%dx2h(i,j)*cs%dy2h(i,j)) / (cs%dx2h(i,j) + cs%dy2h(i,j))
1875 cs%grid_sp_h2(i,j) = grid_sp_h2
1876 grid_sp_h3 = grid_sp_h2*sqrt(grid_sp_h2)
1877 if (cs%Smagorinsky_Kh) cs%Laplac2_const_xx(i,j) = smag_lap_const * grid_sp_h2
1878 if (cs%Leith_Kh) cs%Laplac3_const_xx(i,j) = leith_lap_const * grid_sp_h3
1880 cs%Kh_bg_xx(i,j) = max(kh, kh_vel_scale * sqrt(grid_sp_h2))
1882 if (cs%use_Kh_bg_2d) cs%Kh_bg_xx(i,j) = max(cs%Kh_bg_2d(i,j), cs%Kh_bg_xx(i,j))
1884 if (kh_sin_lat>0.)
then
1885 slat_fn = abs( sin( deg2rad * g%geoLatT(i,j) ) ) ** kh_pwr_of_sine
1886 cs%Kh_bg_xx(i,j) = max(kh_sin_lat * slat_fn, cs%Kh_bg_xx(i,j))
1888 if (cs%bound_Kh .and. .not.cs%better_bound_Kh)
then
1890 cs%Kh_Max_xx(i,j) = kh_limit * grid_sp_h2
1891 cs%Kh_bg_xx(i,j) = min(cs%Kh_bg_xx(i,j), cs%Kh_Max_xx(i,j))
1893 min_grid_sp_h2 = min(grid_sp_h2, min_grid_sp_h2)
1895 call min_across_pes(min_grid_sp_h2)
1898 do j=js-1,jeq ;
do i=is-1,ieq
1900 grid_sp_q2 = (2.0*cs%dx2q(i,j)*cs%dy2q(i,j)) / (cs%dx2q(i,j) + cs%dy2q(i,j))
1901 grid_sp_q3 = grid_sp_q2*sqrt(grid_sp_q2)
1902 if (cs%Smagorinsky_Kh) cs%Laplac2_const_xy(i,j) = smag_lap_const * grid_sp_q2
1903 if (cs%Leith_Kh) cs%Laplac3_const_xy(i,j) = leith_lap_const * grid_sp_q3
1905 cs%Kh_bg_xy(i,j) = max(kh, kh_vel_scale * sqrt(grid_sp_q2))
1907 if (cs%use_Kh_bg_2d)
then
1908 cs%Kh_bg_xy(i,j) = max(cs%Kh_bg_xy(i,j), &
1909 0.25*((cs%Kh_bg_2d(i,j) + cs%Kh_bg_2d(i+1,j+1)) + &
1910 (cs%Kh_bg_2d(i+1,j) + cs%Kh_bg_2d(i,j+1))) )
1914 if (kh_sin_lat>0.)
then
1915 slat_fn = abs( sin( deg2rad * g%geoLatBu(i,j) ) ) ** kh_pwr_of_sine
1916 cs%Kh_bg_xy(i,j) = max(kh_sin_lat * slat_fn, cs%Kh_bg_xy(i,j))
1918 if (cs%bound_Kh .and. .not.cs%better_bound_Kh)
then
1920 cs%Kh_Max_xy(i,j) = kh_limit * grid_sp_q2
1921 cs%Kh_bg_xy(i,j) = min(cs%Kh_bg_xy(i,j), cs%Kh_Max_xy(i,j))
1925 if (cs%biharmonic)
then
1926 do j=js-1,jeq+1 ;
do i=isq-1,ieq+1
1927 cs%Idx2dyCu(i,j) = (g%IdxCu(i,j)*g%IdxCu(i,j)) * g%IdyCu(i,j)
1928 cs%Idxdy2u(i,j) = g%IdxCu(i,j) * (g%IdyCu(i,j)*g%IdyCu(i,j))
1930 do j=jsq-1,jeq+1 ;
do i=is-1,ieq+1
1931 cs%Idx2dyCv(i,j) = (g%IdxCv(i,j)*g%IdxCv(i,j)) * g%IdyCv(i,j)
1932 cs%Idxdy2v(i,j) = g%IdxCv(i,j) * (g%IdyCv(i,j)*g%IdyCv(i,j))
1934 cs%Ah_bg_xy(:,:) = 0.0
1937 if (cs%better_bound_Ah .or. cs%bound_Ah) ah_limit = 0.3 / (dt*64.0)
1938 if (cs%Smagorinsky_Ah .and. cs%bound_Coriolis) &
1939 boundcorconst = 1.0 / (5.0*(bound_cor_vel*bound_cor_vel))
1941 min_grid_sp_h4 = huge(1.)
1942 do j=jsq,jeq+1 ;
do i=isq,ieq+1
1943 grid_sp_h2 = (2.0*cs%dx2h(i,j)*cs%dy2h(i,j)) / (cs%dx2h(i,j)+cs%dy2h(i,j))
1944 grid_sp_h3 = grid_sp_h2*sqrt(grid_sp_h2)
1945 cs%grid_sp_h3(i,j) = grid_sp_h3
1946 if (cs%Smagorinsky_Ah)
then
1947 cs%Biharm_const_xx(i,j) = smag_bi_const * (grid_sp_h2 * grid_sp_h2)
1948 if (cs%bound_Coriolis)
then
1949 fmax = max(abs(g%CoriolisBu(i-1,j-1)), abs(g%CoriolisBu(i,j-1)), &
1950 abs(g%CoriolisBu(i-1,j)), abs(g%CoriolisBu(i,j)))
1951 cs%Biharm_const2_xx(i,j) = (grid_sp_h2 * grid_sp_h2 * grid_sp_h2) * &
1952 (fmax * boundcorconst)
1955 if (cs%Leith_Ah)
then
1956 cs%biharm6_const_xx(i,j) = leith_bi_const * (grid_sp_h3 * grid_sp_h3)
1958 cs%Ah_bg_xx(i,j) = max(ah, ah_vel_scale * grid_sp_h2 * sqrt(grid_sp_h2))
1959 if (cs%Re_Ah > 0.0) cs%Re_Ah_const_xx(i,j) = grid_sp_h3 / cs%Re_Ah
1960 if (ah_time_scale > 0.) cs%Ah_bg_xx(i,j) = &
1961 max(cs%Ah_bg_xx(i,j), (grid_sp_h2 * grid_sp_h2) / ah_time_scale)
1962 if (cs%bound_Ah .and. .not.cs%better_bound_Ah)
then
1963 cs%Ah_Max_xx(i,j) = ah_limit * (grid_sp_h2 * grid_sp_h2)
1964 cs%Ah_bg_xx(i,j) = min(cs%Ah_bg_xx(i,j), cs%Ah_Max_xx(i,j))
1966 min_grid_sp_h4 = min(grid_sp_h2**2, min_grid_sp_h4)
1968 call min_across_pes(min_grid_sp_h4)
1970 do j=js-1,jeq ;
do i=is-1,ieq
1971 grid_sp_q2 = (2.0*cs%dx2q(i,j)*cs%dy2q(i,j)) / (cs%dx2q(i,j)+cs%dy2q(i,j))
1972 grid_sp_q3 = grid_sp_q2*sqrt(grid_sp_q2)
1973 if (cs%Smagorinsky_Ah)
then
1974 cs%Biharm_const_xy(i,j) = smag_bi_const * (grid_sp_q2 * grid_sp_q2)
1975 if (cs%bound_Coriolis)
then
1976 cs%Biharm_const2_xy(i,j) = (grid_sp_q2 * grid_sp_q2 * grid_sp_q2) * &
1977 (abs(g%CoriolisBu(i,j)) * boundcorconst)
1980 if (cs%Leith_Ah)
then
1981 cs%biharm6_const_xy(i,j) = leith_bi_const * (grid_sp_q3 * grid_sp_q3)
1983 cs%Ah_bg_xy(i,j) = max(ah, ah_vel_scale * grid_sp_q2 * sqrt(grid_sp_q2))
1984 if (cs%Re_Ah > 0.0) cs%Re_Ah_const_xy(i,j) = grid_sp_q3 / cs%Re_Ah
1985 if (ah_time_scale > 0.) cs%Ah_bg_xy(i,j) = &
1986 max(cs%Ah_bg_xy(i,j), (grid_sp_q2 * grid_sp_q2) / ah_time_scale)
1987 if (cs%bound_Ah .and. .not.cs%better_bound_Ah)
then
1988 cs%Ah_Max_xy(i,j) = ah_limit * (grid_sp_q2 * grid_sp_q2)
1989 cs%Ah_bg_xy(i,j) = min(cs%Ah_bg_xy(i,j), cs%Ah_Max_xy(i,j))
1994 if (cs%Laplacian .and. cs%better_bound_Kh)
then
1995 do j=jsq,jeq+1 ;
do i=isq,ieq+1
1997 (cs%dy2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) * &
1998 max(g%IdyCu(i,j)*g%IareaCu(i,j), g%IdyCu(i-1,j)*g%IareaCu(i-1,j)) ), &
1999 (cs%dx2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) * &
2000 max(g%IdxCv(i,j)*g%IareaCv(i,j), g%IdxCv(i,j-1)*g%IareaCv(i,j-1)) ) )
2001 cs%Kh_Max_xx(i,j) = 0.0
2003 cs%Kh_Max_xx(i,j) = cs%bound_coef * 0.25 * idt / denom
2005 do j=js-1,jeq ;
do i=is-1,ieq
2007 (cs%dx2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) * &
2008 max(g%IdxCu(i,j)*g%IareaCu(i,j), g%IdxCu(i,j+1)*g%IareaCu(i,j+1)) ), &
2009 (cs%dy2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) * &
2010 max(g%IdyCv(i,j)*g%IareaCv(i,j), g%IdyCv(i+1,j)*g%IareaCv(i+1,j)) ) )
2011 cs%Kh_Max_xy(i,j) = 0.0
2013 cs%Kh_Max_xy(i,j) = cs%bound_coef * 0.25 * idt / denom
2016 call hchksum(cs%Kh_Max_xx,
"Kh_Max_xx", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
2017 call bchksum(cs%Kh_Max_xy,
"Kh_Max_xy", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
2022 if (cs%biharmonic .and. cs%better_bound_Ah)
then
2023 do j=js-1,jeq+1 ;
do i=isq-1,ieq+1
2024 u0u(i,j) = (cs%Idxdy2u(i,j)*(cs%dy2h(i+1,j)*cs%DY_dxT(i+1,j)*(g%IdyCu(i+1,j) + g%IdyCu(i,j)) + &
2025 cs%dy2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) ) + &
2026 cs%Idx2dyCu(i,j)*(cs%dx2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) + &
2027 cs%dx2q(i,j-1)*cs%DX_dyBu(i,j-1)*(g%IdxCu(i,j) + g%IdxCu(i,j-1)) ) )
2028 u0v(i,j) = (cs%Idxdy2u(i,j)*(cs%dy2h(i+1,j)*cs%DX_dyT(i+1,j)*(g%IdxCv(i+1,j) + g%IdxCv(i+1,j-1)) + &
2029 cs%dy2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) ) + &
2030 cs%Idx2dyCu(i,j)*(cs%dx2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) + &
2031 cs%dx2q(i,j-1)*cs%DY_dxBu(i,j-1)*(g%IdyCv(i+1,j-1) + g%IdyCv(i,j-1)) ) )
2033 do j=jsq-1,jeq+1 ;
do i=is-1,ieq+1
2034 v0u(i,j) = (cs%Idxdy2v(i,j)*(cs%dy2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) + &
2035 cs%dy2q(i-1,j)*cs%DX_dyBu(i-1,j)*(g%IdxCu(i-1,j+1) + g%IdxCu(i-1,j)) ) + &
2036 cs%Idx2dyCv(i,j)*(cs%dx2h(i,j+1)*cs%DY_dxT(i,j+1)*(g%IdyCu(i,j+1) + g%IdyCu(i-1,j+1)) + &
2037 cs%dx2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) ) )
2038 v0v(i,j) = (cs%Idxdy2v(i,j)*(cs%dy2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) + &
2039 cs%dy2q(i-1,j)*cs%DY_dxBu(i-1,j)*(g%IdyCv(i,j) + g%IdyCv(i-1,j)) ) + &
2040 cs%Idx2dyCv(i,j)*(cs%dx2h(i,j+1)*cs%DX_dyT(i,j+1)*(g%IdxCv(i,j+1) + g%IdxCv(i,j)) + &
2041 cs%dx2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) ) )
2043 do j=jsq,jeq+1 ;
do i=isq,ieq+1
2046 (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0u(i,j) + g%IdyCu(i-1,j)*u0u(i-1,j)) + &
2047 cs%DX_dyT(i,j)*(g%IdxCv(i,j)*v0u(i,j) + g%IdxCv(i,j-1)*v0u(i,j-1))) * &
2048 max(g%IdyCu(i,j)*g%IareaCu(i,j), g%IdyCu(i-1,j)*g%IareaCu(i-1,j)) ), &
2050 (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0v(i,j) + g%IdyCu(i-1,j)*u0v(i-1,j)) + &
2051 cs%DX_dyT(i,j)*(g%IdxCv(i,j)*v0v(i,j) + g%IdxCv(i,j-1)*v0v(i,j-1))) * &
2052 max(g%IdxCv(i,j)*g%IareaCv(i,j), g%IdxCv(i,j-1)*g%IareaCv(i,j-1)) ) )
2053 cs%Ah_Max_xx(i,j) = 0.0
2055 cs%Ah_Max_xx(i,j) = cs%bound_coef * 0.5 * idt / denom
2057 do j=js-1,jeq ;
do i=is-1,ieq
2060 (cs%DX_dyBu(i,j)*(u0u(i,j+1)*g%IdxCu(i,j+1) + u0u(i,j)*g%IdxCu(i,j)) + &
2061 cs%DY_dxBu(i,j)*(v0u(i+1,j)*g%IdyCv(i+1,j) + v0u(i,j)*g%IdyCv(i,j))) * &
2062 max(g%IdxCu(i,j)*g%IareaCu(i,j), g%IdxCu(i,j+1)*g%IareaCu(i,j+1)) ), &
2064 (cs%DX_dyBu(i,j)*(u0v(i,j+1)*g%IdxCu(i,j+1) + u0v(i,j)*g%IdxCu(i,j)) + &
2065 cs%DY_dxBu(i,j)*(v0v(i+1,j)*g%IdyCv(i+1,j) + v0v(i,j)*g%IdyCv(i,j))) * &
2066 max(g%IdyCv(i,j)*g%IareaCv(i,j), g%IdyCv(i+1,j)*g%IareaCv(i+1,j)) ) )
2067 cs%Ah_Max_xy(i,j) = 0.0
2069 cs%Ah_Max_xy(i,j) = cs%bound_coef * 0.5 * idt / denom
2072 call hchksum(cs%Ah_Max_xx,
"Ah_Max_xx", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
2073 call bchksum(cs%Ah_Max_xy,
"Ah_Max_xy", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
2077 cs%id_diffu = register_diag_field(
'ocean_model',
'diffu', diag%axesCuL, time, &
2078 'Zonal Acceleration from Horizontal Viscosity',
'm s-2', conversion=us%L_T2_to_m_s2)
2079 cs%id_diffv = register_diag_field(
'ocean_model',
'diffv', diag%axesCvL, time, &
2080 'Meridional Acceleration from Horizontal Viscosity',
'm s-2', conversion=us%L_T2_to_m_s2)
2098 cs%id_hf_diffu_2d = register_diag_field(
'ocean_model',
'hf_diffu_2d', diag%axesCu1, time, &
2099 'Depth-sum Fractional Thickness-weighted Zonal Acceleration from Horizontal Viscosity',
'm s-2', &
2100 conversion=us%L_T2_to_m_s2)
2101 if ((cs%id_hf_diffu_2d > 0) .and. (
present(adp)))
then
2102 call safe_alloc_ptr(adp%diag_hfrac_u,g%IsdB,g%IedB,g%jsd,g%jed,g%ke)
2105 cs%id_hf_diffv_2d = register_diag_field(
'ocean_model',
'hf_diffv_2d', diag%axesCv1, time, &
2106 'Depth-sum Fractional Thickness-weighted Meridional Acceleration from Horizontal Viscosity',
'm s-2', &
2107 conversion=us%L_T2_to_m_s2)
2108 if ((cs%id_hf_diffv_2d > 0) .and. (
present(adp)))
then
2109 call safe_alloc_ptr(adp%diag_hfrac_v,g%isd,g%ied,g%JsdB,g%JedB,g%ke)
2112 if (cs%biharmonic)
then
2113 cs%id_Ah_h = register_diag_field(
'ocean_model',
'Ahh', diag%axesTL, time, &
2114 'Biharmonic Horizontal Viscosity at h Points',
'm4 s-1', conversion=us%L_to_m**4*us%s_to_T, &
2115 cmor_field_name=
'difmxybo', &
2116 cmor_long_name=
'Ocean lateral biharmonic viscosity', &
2117 cmor_standard_name=
'ocean_momentum_xy_biharmonic_diffusivity')
2118 cs%id_Ah_q = register_diag_field(
'ocean_model',
'Ahq', diag%axesBL, time, &
2119 'Biharmonic Horizontal Viscosity at q Points',
'm4 s-1', conversion=us%L_to_m**4*us%s_to_T)
2120 cs%id_grid_Re_Ah = register_diag_field(
'ocean_model',
'grid_Re_Ah', diag%axesTL, time, &
2121 'Grid Reynolds number for the Biharmonic horizontal viscosity at h points',
'nondim')
2123 if (cs%id_grid_Re_Ah > 0) &
2126 cs%min_grid_Ah = spacing(1.) * min_grid_sp_h4 * idt
2128 if (cs%Laplacian)
then
2129 cs%id_Kh_h = register_diag_field(
'ocean_model',
'Khh', diag%axesTL, time, &
2130 'Laplacian Horizontal Viscosity at h Points',
'm2 s-1', conversion=us%L_to_m**2*us%s_to_T, &
2131 cmor_field_name=
'difmxylo', &
2132 cmor_long_name=
'Ocean lateral Laplacian viscosity', &
2133 cmor_standard_name=
'ocean_momentum_xy_laplacian_diffusivity')
2134 cs%id_Kh_q = register_diag_field(
'ocean_model',
'Khq', diag%axesBL, time, &
2135 'Laplacian Horizontal Viscosity at q Points',
'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
2136 cs%id_grid_Re_Kh = register_diag_field(
'ocean_model',
'grid_Re_Kh', diag%axesTL, time, &
2137 'Grid Reynolds number for the Laplacian horizontal viscosity at h points',
'nondim')
2138 cs%id_vort_xy_q = register_diag_field(
'ocean_model',
'vort_xy_q', diag%axesBL, time, &
2139 'Vertical vorticity at q Points',
's-1', conversion=us%s_to_T)
2140 cs%id_div_xx_h = register_diag_field(
'ocean_model',
'div_xx_h', diag%axesTL, time, &
2141 'Horizontal divergence at h Points',
's-1', conversion=us%s_to_T)
2142 cs%id_sh_xy_q = register_diag_field(
'ocean_model',
'sh_xy_q', diag%axesBL, time, &
2143 'Shearing strain at q Points',
's-1', conversion=us%s_to_T)
2144 cs%id_sh_xx_h = register_diag_field(
'ocean_model',
'sh_xx_h', diag%axesTL, time, &
2145 'Horizontal tension at h Points',
's-1', conversion=us%s_to_T)
2147 if (cs%id_grid_Re_Kh > 0) &
2150 cs%min_grid_Kh = spacing(1.) * min_grid_sp_h2 * idt
2152 if (cs%use_GME)
then
2153 cs%id_GME_coeff_h = register_diag_field(
'ocean_model',
'GME_coeff_h', diag%axesTL, time, &
2154 'GME coefficient at h Points',
'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
2155 cs%id_GME_coeff_q = register_diag_field(
'ocean_model',
'GME_coeff_q', diag%axesBL, time, &
2156 'GME coefficient at q Points',
'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
2157 cs%id_FrictWork_GME = register_diag_field(
'ocean_model',
'FrictWork_GME',diag%axesTL,time,&
2158 'Integral work done by lateral friction terms in GME (excluding diffusion of energy)', &
2159 'W m-2', conversion=us%RZ3_T3_to_W_m2*us%L_to_Z**2)
2161 cs%id_FrictWork = register_diag_field(
'ocean_model',
'FrictWork',diag%axesTL,time,&
2162 'Integral work done by lateral friction terms', &
2163 'W m-2', conversion=us%RZ3_T3_to_W_m2*us%L_to_Z**2)
2164 cs%id_FrictWorkIntz = register_diag_field(
'ocean_model',
'FrictWorkIntz',diag%axesT1,time, &
2165 'Depth integrated work done by lateral friction', &
2166 'W m-2', conversion=us%RZ3_T3_to_W_m2*us%L_to_Z**2, &
2167 cmor_field_name=
'dispkexyfo', &
2168 cmor_long_name=
'Depth integrated ocean kinetic energy dissipation due to lateral friction',&
2169 cmor_standard_name=
'ocean_kinetic_energy_dissipation_per_unit_area_due_to_xy_friction')
2170 if (cs%Laplacian .or. get_all)
then
2172 end subroutine hor_visc_init
2175 subroutine align_aniso_tensor_to_grid(CS, n1, n2)
2177 real,
intent(in) :: n1
2178 real,
intent(in) :: n2
2180 real :: recip_n2_norm
2182 recip_n2_norm = n1**2 + n2**2
2183 if (recip_n2_norm > 0.) recip_n2_norm = 1./recip_n2_norm
2184 cs%n1n2_h(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2185 cs%n1n2_q(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2186 cs%n1n1_m_n2n2_h(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm
2187 cs%n1n1_m_n2n2_q(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm
2188 end subroutine align_aniso_tensor_to_grid
2191 subroutine smooth_gme(CS,G,GME_flux_h,GME_flux_q)
2195 real,
dimension(SZI_(G),SZJ_(G)),
optional,
intent(inout) :: GME_flux_h
2197 real,
dimension(SZIB_(G),SZJB_(G)),
optional,
intent(inout) :: GME_flux_q
2200 real,
dimension(SZI_(G),SZJ_(G)) :: GME_flux_h_original
2201 real,
dimension(SZIB_(G),SZJB_(G)) :: GME_flux_q_original
2202 real :: wc, ww, we, wn, ws
2203 integer :: i, j, k, s
2206 if (
present(gme_flux_h))
then
2207 call pass_var(gme_flux_h, g%Domain)
2208 gme_flux_h_original = gme_flux_h
2213 if (g%mask2dT(i,j)==0.) cycle
2215 ww = 0.125 * g%mask2dT(i-1,j)
2216 we = 0.125 * g%mask2dT(i+1,j)
2217 ws = 0.125 * g%mask2dT(i,j-1)
2218 wn = 0.125 * g%mask2dT(i,j+1)
2219 wc = 1.0 - (ww+we+wn+ws)
2220 gme_flux_h(i,j) = wc * gme_flux_h_original(i,j) &
2221 + ww * gme_flux_h_original(i-1,j) &
2222 + we * gme_flux_h_original(i+1,j) &
2223 + ws * gme_flux_h_original(i,j-1) &
2224 + wn * gme_flux_h_original(i,j+1)
2228 if (
present(gme_flux_q))
then
2229 call pass_var(gme_flux_q, g%Domain, position=corner, complete=.true.)
2230 gme_flux_q_original = gme_flux_q
2232 do j = g%JscB, g%JecB
2233 do i = g%IscB, g%IecB
2235 if (g%mask2dBu(i,j)==0.) cycle
2237 ww = 0.125 * g%mask2dBu(i-1,j)
2238 we = 0.125 * g%mask2dBu(i+1,j)
2239 ws = 0.125 * g%mask2dBu(i,j-1)
2240 wn = 0.125 * g%mask2dBu(i,j+1)
2241 wc = 1.0 - (ww+we+wn+ws)
2242 gme_flux_q(i,j) = wc * gme_flux_q_original(i,j) &
2243 + ww * gme_flux_q_original(i-1,j) &
2244 + we * gme_flux_q_original(i+1,j) &
2245 + ws * gme_flux_q_original(i,j-1) &
2246 + wn * gme_flux_q_original(i,j+1)
2250 end subroutine smooth_gme
2252 subroutine hor_visc_end(CS)
2255 if (cs%Laplacian .or. cs%biharmonic)
then
2256 dealloc_(cs%dx2h) ; dealloc_(cs%dx2q) ; dealloc_(cs%dy2h) ; dealloc_(cs%dy2q)
2257 dealloc_(cs%dx_dyT) ; dealloc_(cs%dy_dxT) ; dealloc_(cs%dx_dyBu) ; dealloc_(cs%dy_dxBu)
2258 dealloc_(cs%reduction_xx) ; dealloc_(cs%reduction_xy)
2260 if (cs%Laplacian)
then
2261 dealloc_(cs%Kh_bg_xx) ; dealloc_(cs%Kh_bg_xy)
2262 dealloc_(cs%grid_sp_h2)
2263 if (cs%bound_Kh)
then
2264 dealloc_(cs%Kh_Max_xx) ; dealloc_(cs%Kh_Max_xy)
2266 if (cs%Smagorinsky_Kh)
then
2267 dealloc_(cs%Laplac2_const_xx) ; dealloc_(cs%Laplac2_const_xy)
2269 if (cs%Leith_Kh)
then
2270 dealloc_(cs%Laplac3_const_xx) ; dealloc_(cs%Laplac3_const_xy)
2273 if (cs%biharmonic)
then
2274 dealloc_(cs%grid_sp_h3)
2275 dealloc_(cs%Idx2dyCu) ; dealloc_(cs%Idx2dyCv)
2276 dealloc_(cs%Idxdy2u) ; dealloc_(cs%Idxdy2v)
2277 dealloc_(cs%Ah_bg_xx) ; dealloc_(cs%Ah_bg_xy)
2278 if (cs%bound_Ah)
then
2279 dealloc_(cs%Ah_Max_xx) ; dealloc_(cs%Ah_Max_xy)
2281 if (cs%Smagorinsky_Ah)
then
2282 dealloc_(cs%Biharm6_const_xx) ; dealloc_(cs%Biharm6_const_xy)
2284 if (cs%Leith_Ah)
then
2285 dealloc_(cs%Biharm_const_xx) ; dealloc_(cs%Biharm_const_xy)
2287 if (cs%Re_Ah > 0.0)
then
2288 dealloc_(cs%Re_Ah_const_xx) ; dealloc_(cs%Re_Ah_const_xy)
2291 if (cs%anisotropic)
then
2294 dealloc_(cs%n1n1_m_n2n2_h)
2295 dealloc_(cs%n1n1_m_n2n2_q)
2298 end subroutine hor_visc_end