6 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
7 use mom_cpu_clock,
only : clock_module_driver, clock_module, clock_routine
19 implicit none ;
private
21 #include <MOM_memory.h>
23 public calculate_kappa_shear, calc_kappa_shear_vertex, kappa_shear_init
24 public kappa_shear_is_used, kappa_shear_at_vertex
64 integer :: max_rino_it
68 logical :: dkdq_iteration_bug
71 logical :: ks_at_vertex
73 logical :: eliminate_massless
78 real :: kappa_src_max_chg
86 logical :: all_layer_tke_bug
90 logical :: restrictive_tolerance_check
96 logical :: debug = .false.
100 integer :: id_kd_shear = -1, id_tke = -1, id_ild2 = -1, id_dz_int = -1
109 subroutine calculate_kappa_shear(u_in, v_in, h, tv, p_surf, kappa_io, tke_io, &
110 kv_io, dt, G, GV, US, CS, initialize_all)
114 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
116 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
118 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
123 real,
dimension(:,:),
pointer :: p_surf
124 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), &
125 intent(inout) :: kappa_io
129 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), &
130 intent(out) :: tke_io
132 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), &
133 intent(inout) :: kv_io
137 real,
intent(in) :: dt
140 logical,
optional,
intent(in) :: initialize_all
144 real,
dimension(SZI_(G),SZK_(GV)) :: &
148 real,
dimension(SZI_(G),SZK_(GV)+1) :: &
151 real,
dimension(SZK_(GV)) :: &
158 real,
dimension(SZK_(GV)+1) :: &
169 logical :: use_temperature
171 logical :: new_kappa = .true.
174 integer,
dimension(SZK_(GV)+1) :: kc
177 real,
dimension(SZK_(GV)+1) :: kf
179 integer :: is, ie, js, je, i, j, k, nz, nzc
181 is = g%isc ; ie = g%iec; js = g%jsc ; je = g%jec ; nz = gv%ke
183 use_temperature = .false. ;
if (
associated(tv%T)) use_temperature = .true.
184 new_kappa = .true. ;
if (
present(initialize_all)) new_kappa = initialize_all
187 dz_massless = 0.1*sqrt(k0dt)
192 do k=1,nz ;
do i=is,ie
193 h_2d(i,k) = h(i,j,k)*gv%H_to_Z
194 u_2d(i,k) = u_in(i,j,k) ; v_2d(i,k) = v_in(i,j,k)
196 if (use_temperature)
then ;
do k=1,nz ;
do i=is,ie
197 t_2d(i,k) = tv%T(i,j,k) ; s_2d(i,k) = tv%S(i,j,k)
198 enddo ;
enddo ;
else ;
do k=1,nz ;
do i=is,ie
199 rho_2d(i,k) = gv%Rlay(k)
200 enddo ;
enddo ;
endif
201 if (.not.new_kappa)
then ;
do k=1,nz+1 ;
do i=is,ie
202 kappa_2d(i,k) = kappa_io(i,j,k)
203 enddo ;
enddo ;
endif
208 do i=is,ie ;
if (g%mask2dT(i,j) > 0.5)
then
212 if (cs%eliminate_massless)
then
216 dz(k) = 0.0 ; u0xdz(k) = 0.0 ; v0xdz(k) = 0.0
217 t0xdz(k) = 0.0 ; s0xdz(k) = 0.0
221 if ((k>cs%nkml) .and. (dz(nzc) > 0.0) .and. &
222 (h_2d(i,k) > dz_massless)) nzc = nzc+1
229 dz(nzc) = dz(nzc) + h_2d(i,k)
230 u0xdz(nzc) = u0xdz(nzc) + u_2d(i,k)*h_2d(i,k)
231 v0xdz(nzc) = v0xdz(nzc) + v_2d(i,k)*h_2d(i,k)
232 if (use_temperature)
then
233 t0xdz(nzc) = t0xdz(nzc) + t_2d(i,k)*h_2d(i,k)
234 s0xdz(nzc) = s0xdz(nzc) + s_2d(i,k)*h_2d(i,k)
236 t0xdz(nzc) = t0xdz(nzc) + rho_2d(i,k)*h_2d(i,k)
237 s0xdz(nzc) = s0xdz(nzc) + rho_2d(i,k)*h_2d(i,k)
243 do k=1,nzc ; idz(k) = 1.0 / dz(k) ;
enddo
247 kf(1) = 0.0 ; dz_in_lay = h_2d(i,1)
249 if (kc(k) > kc(k-1))
then
250 kf(k) = 0.0 ; dz_in_lay = h_2d(i,k)
252 kf(k) = dz_in_lay*idz(kc(k)) ; dz_in_lay = dz_in_lay + h_2d(i,k)
259 u0xdz(k) = u_2d(i,k)*dz(k) ; v0xdz(k) = v_2d(i,k)*dz(k)
261 if (use_temperature)
then
263 t0xdz(k) = t_2d(i,k)*dz(k) ; s0xdz(k) = s_2d(i,k)*dz(k)
267 t0xdz(k) = rho_2d(i,k)*dz(k) ; s0xdz(k) = rho_2d(i,k)*dz(k)
271 do k=1,nzc+1 ; kc(k) = k ; kf(k) = 0.0 ;
enddo
273 f2 = 0.25 * ((g%CoriolisBu(i,j)**2 + g%CoriolisBu(i-1,j-1)**2) + &
274 (g%CoriolisBu(i,j-1)**2 + g%CoriolisBu(i-1,j)**2))
275 surface_pres = 0.0 ;
if (
associated(p_surf)) surface_pres = p_surf(i,j)
282 do k=1,nzc+1 ; kappa(k) = us%m2_s_to_Z2_T*1.0 ;
enddo
284 do k=1,nzc+1 ; kappa(k) = kappa_2d(i,k) ;
enddo
287 call kappa_shear_column(kappa, tke, dt, nzc, f2, surface_pres, &
288 dz, u0xdz, v0xdz, t0xdz, s0xdz, kappa_avg, &
289 tke_avg, tv, cs, gv, us)
295 kappa_2d(i,k) = kappa_avg(k)
296 if (cs%all_layer_TKE_bug)
then
299 tke_2d(i,k) = tke_avg(k)
304 if (kf(k) == 0.0)
then
305 kappa_2d(i,k) = kappa_avg(kc(k))
306 tke_2d(i,k) = tke_avg(kc(k))
308 kappa_2d(i,k) = (1.0-kf(k)) * kappa_avg(kc(k)) + &
309 kf(k) * kappa_avg(kc(k)+1)
310 tke_2d(i,k) = (1.0-kf(k)) * tke_avg(kc(k)) + &
311 kf(k) * tke_avg(kc(k)+1)
318 kappa_2d(i,k) = 0.0 ; tke_2d(i,k) = 0.0
322 do k=1,nz+1 ;
do i=is,ie
323 kappa_io(i,j,k) = g%mask2dT(i,j) * kappa_2d(i,k)
324 tke_io(i,j,k) = g%mask2dT(i,j) * tke_2d(i,k)
325 kv_io(i,j,k) = ( g%mask2dT(i,j) * kappa_2d(i,k) ) * cs%Prandtl_turb
331 call hchksum(kappa_io,
"kappa", g%HI, scale=us%Z2_T_to_m2_s)
332 call hchksum(tke_io,
"tke", g%HI, scale=us%Z_to_m**2*us%s_to_T**2)
335 if (cs%id_Kd_shear > 0)
call post_data(cs%id_Kd_shear, kappa_io, cs%diag)
336 if (cs%id_TKE > 0)
call post_data(cs%id_TKE, tke_io, cs%diag)
338 end subroutine calculate_kappa_shear
342 subroutine calc_kappa_shear_vertex(u_in, v_in, h, T_in, S_in, tv, p_surf, kappa_io, tke_io, &
343 kv_io, dt, G, GV, US, CS, initialize_all)
347 real,
dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
349 real,
dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
351 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
353 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
355 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
360 real,
dimension(:,:),
pointer :: p_surf
362 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), &
363 intent(out) :: kappa_io
365 real,
dimension(SZIB_(G),SZJB_(G),SZK_(GV)+1), &
366 intent(out) :: tke_io
368 real,
dimension(SZIB_(G),SZJB_(G),SZK_(GV)+1), &
369 intent(inout) :: kv_io
373 real,
intent(in) :: dt
376 logical,
optional,
intent(in) :: initialize_all
380 real,
dimension(SZIB_(G),SZK_(GV)) :: &
384 real,
dimension(SZIB_(G),SZK_(GV)+1,2) :: &
386 real,
dimension(SZIB_(G),SZK_(GV)+1) :: &
388 real,
dimension(SZK_(GV)) :: &
395 real,
dimension(SZK_(GV)+1) :: &
408 logical :: use_temperature
410 logical :: new_kappa = .true.
414 integer,
dimension(SZK_(GV)+1) :: kc
417 real,
dimension(SZK_(GV)+1) :: kf
419 integer :: isb, ieb, jsb, jeb, i, j, k, nz, nzc, j2, j2m1
422 isb = g%isc-1 ; ieb = g%iecB ; jsb = g%jsc-1 ; jeb = g%jecB ; nz = gv%ke
424 use_temperature = .false. ;
if (
associated(tv%T)) use_temperature = .true.
425 new_kappa = .true. ;
if (
present(initialize_all)) new_kappa = initialize_all
428 dz_massless = 0.1*sqrt(k0dt)
429 i_prandtl = 0.0 ;
if (cs%Prandtl_turb > 0.0) i_prandtl = 1.0 / cs%Prandtl_turb
434 j2 = mod(j,2)+1 ; j2m1 = 3-j2
437 do k=1,nz ;
do i=isb,ieb
438 u_2d(i,k) = (u_in(i,j,k) * (g%mask2dCu(i,j) * (h(i,j,k) + h(i+1,j,k))) + &
439 u_in(i,j+1,k) * (g%mask2dCu(i,j+1) * (h(i,j+1,k) + h(i+1,j+1,k))) ) / &
440 ((g%mask2dCu(i,j) * (h(i,j,k) + h(i+1,j,k)) + &
441 g%mask2dCu(i,j+1) * (h(i,j+1,k) + h(i+1,j+1,k))) + gv%H_subroundoff)
442 v_2d(i,k) = (v_in(i,j,k) * (g%mask2dCv(i,j) * (h(i,j,k) + h(i,j+1,k))) + &
443 v_in(i+1,j,k) * (g%mask2dCv(i+1,j) * (h(i+1,j,k) + h(i+1,j+1,k))) ) / &
444 ((g%mask2dCv(i,j) * (h(i,j,k) + h(i,j+1,k)) + &
445 g%mask2dCv(i+1,j) * (h(i+1,j,k) + h(i+1,j+1,k))) + gv%H_subroundoff)
446 i_hwt = 1.0 / (((g%mask2dT(i,j) * h(i,j,k) + g%mask2dT(i+1,j+1) * h(i+1,j+1,k)) + &
447 (g%mask2dT(i+1,j) * h(i+1,j,k) + g%mask2dT(i,j+1) * h(i,j+1,k))) + &
449 if (use_temperature)
then
450 t_2d(i,k) = ( ((g%mask2dT(i,j) * h(i,j,k)) * t_in(i,j,k) + &
451 (g%mask2dT(i+1,j+1) * h(i+1,j+1,k)) * t_in(i+1,j+1,k)) + &
452 ((g%mask2dT(i+1,j) * h(i+1,j,k)) * t_in(i+1,j,k) + &
453 (g%mask2dT(i,j+1) * h(i,j+1,k)) * t_in(i,j+1,k)) ) * i_hwt
454 s_2d(i,k) = ( ((g%mask2dT(i,j) * h(i,j,k)) * s_in(i,j,k) + &
455 (g%mask2dT(i+1,j+1) * h(i+1,j+1,k)) * s_in(i+1,j+1,k)) + &
456 ((g%mask2dT(i+1,j) * h(i+1,j,k)) * s_in(i+1,j,k) + &
457 (g%mask2dT(i,j+1) * h(i,j+1,k)) * s_in(i,j+1,k)) ) * i_hwt
459 h_2d(i,k) = gv%H_to_Z * ((g%mask2dT(i,j) * h(i,j,k) + g%mask2dT(i+1,j+1) * h(i+1,j+1,k)) + &
460 (g%mask2dT(i+1,j) * h(i+1,j,k) + g%mask2dT(i,j+1) * h(i,j+1,k)) ) / &
461 ((g%mask2dT(i,j) + g%mask2dT(i+1,j+1)) + &
462 (g%mask2dT(i+1,j) + g%mask2dT(i,j+1)) + 1.0e-36 )
467 if (.not.use_temperature)
then ;
do k=1,nz ;
do i=isb,ieb
468 rho_2d(i,k) = gv%Rlay(k)
469 enddo ;
enddo ;
endif
470 if (.not.new_kappa)
then ;
do k=1,nz+1 ;
do i=isb,ieb
471 kappa_2d(i,k,j2) = kv_io(i,j,k) * i_prandtl
472 enddo ;
enddo ;
endif
477 do i=isb,ieb ;
if ((g%mask2dCu(i,j) + g%mask2dCu(i,j+1)) + &
478 (g%mask2dCv(i,j) + g%mask2dCv(i+1,j)) > 0.0)
then
482 if (cs%eliminate_massless)
then
486 dz(k) = 0.0 ; u0xdz(k) = 0.0 ; v0xdz(k) = 0.0
487 t0xdz(k) = 0.0 ; s0xdz(k) = 0.0
491 if ((k>cs%nkml) .and. (dz(nzc) > 0.0) .and. &
492 (h_2d(i,k) > dz_massless)) nzc = nzc+1
499 dz(nzc) = dz(nzc) + h_2d(i,k)
500 u0xdz(nzc) = u0xdz(nzc) + u_2d(i,k)*h_2d(i,k)
501 v0xdz(nzc) = v0xdz(nzc) + v_2d(i,k)*h_2d(i,k)
502 if (use_temperature)
then
503 t0xdz(nzc) = t0xdz(nzc) + t_2d(i,k)*h_2d(i,k)
504 s0xdz(nzc) = s0xdz(nzc) + s_2d(i,k)*h_2d(i,k)
506 t0xdz(nzc) = t0xdz(nzc) + rho_2d(i,k)*h_2d(i,k)
507 s0xdz(nzc) = s0xdz(nzc) + rho_2d(i,k)*h_2d(i,k)
513 do k=1,nzc ; idz(k) = 1.0 / dz(k) ;
enddo
517 kf(1) = 0.0 ; dz_in_lay = h_2d(i,1)
519 if (kc(k) > kc(k-1))
then
520 kf(k) = 0.0 ; dz_in_lay = h_2d(i,k)
522 kf(k) = dz_in_lay*idz(kc(k)) ; dz_in_lay = dz_in_lay + h_2d(i,k)
529 u0xdz(k) = u_2d(i,k)*dz(k) ; v0xdz(k) = v_2d(i,k)*dz(k)
531 if (use_temperature)
then
533 t0xdz(k) = t_2d(i,k)*dz(k) ; s0xdz(k) = s_2d(i,k)*dz(k)
537 t0xdz(k) = rho_2d(i,k)*dz(k) ; s0xdz(k) = rho_2d(i,k)*dz(k)
541 do k=1,nzc+1 ; kc(k) = k ; kf(k) = 0.0 ;
enddo
543 f2 = g%CoriolisBu(i,j)**2
545 if (
associated(p_surf))
then
546 if (cs%psurf_bug)
then
548 surface_pres = 0.25 * ((p_surf(i,j) + p_surf(i+1,j+1)) + &
549 (p_surf(i+1,j) + p_surf(i,j+1)))
551 surface_pres = ((g%mask2dT(i,j) * p_surf(i,j) + g%mask2dT(i+1,j+1) * p_surf(i+1,j+1)) + &
552 (g%mask2dT(i+1,j) * p_surf(i+1,j) + g%mask2dT(i,j+1) * p_surf(i,j+1)) ) / &
553 ((g%mask2dT(i,j) + g%mask2dT(i+1,j+1)) + &
554 (g%mask2dT(i+1,j) + g%mask2dT(i,j+1)) + 1.0e-36 )
562 do k=1,nzc+1 ; kappa(k) = us%m2_s_to_Z2_T*1.0 ;
enddo
564 do k=1,nzc+1 ; kappa(k) = kappa_2d(i,k,j2) ;
enddo
567 call kappa_shear_column(kappa, tke, dt, nzc, f2, surface_pres, &
568 dz, u0xdz, v0xdz, t0xdz, s0xdz, kappa_avg, &
569 tke_avg, tv, cs, gv, us)
574 kappa_2d(i,k,j2) = kappa_avg(k)
575 if (cs%all_layer_TKE_bug)
then
578 tke_2d(i,k) = tke_avg(k)
583 if (kf(k) == 0.0)
then
584 kappa_2d(i,k,j2) = kappa_avg(kc(k))
585 tke_2d(i,k) = tke_avg(kc(k))
587 kappa_2d(i,k,j2) = (1.0-kf(k)) * kappa_avg(kc(k)) + kf(k) * kappa_avg(kc(k)+1)
588 tke_2d(i,k) = (1.0-kf(k)) * tke_avg(kc(k)) + kf(k) * tke_avg(kc(k)+1)
595 kappa_2d(i,k,j2) = 0.0 ; tke_2d(i,k) = 0.0
599 do k=1,nz+1 ;
do i=isb,ieb
600 tke_io(i,j,k) = g%mask2dBu(i,j) * tke_2d(i,k)
601 kv_io(i,j,k) = ( g%mask2dBu(i,j) * kappa_2d(i,k,j2) ) * cs%Prandtl_turb
603 if (j>=g%jsc)
then ;
do k=1,nz+1 ;
do i=g%isc,g%iec
605 kappa_io(i,j,k) = g%mask2dT(i,j) * 0.25 * &
606 ((kappa_2d(i-1,k,j2m1) + kappa_2d(i,k,j2)) + &
607 (kappa_2d(i-1,k,j2) + kappa_2d(i,k,j2m1)))
608 enddo ;
enddo ;
endif
613 call hchksum(kappa_io,
"kappa", g%HI, scale=us%Z2_T_to_m2_s)
614 call bchksum(tke_io,
"tke", g%HI, scale=us%Z_to_m**2*us%s_to_T**2)
617 if (cs%id_Kd_shear > 0)
call post_data(cs%id_Kd_shear, kappa_io, cs%diag)
618 if (cs%id_TKE > 0)
call post_data(cs%id_TKE, tke_io, cs%diag)
620 end subroutine calc_kappa_shear_vertex
624 subroutine kappa_shear_column(kappa, tke, dt, nzc, f2, surface_pres, &
625 dz, u0xdz, v0xdz, T0xdz, S0xdz, kappa_avg, &
626 tke_avg, tv, CS, GV, US, I_Ld2_1d, dz_Int_1d)
628 real,
dimension(SZK_(GV)+1), &
629 intent(inout) :: kappa
630 real,
dimension(SZK_(GV)+1), &
633 integer,
intent(in) :: nzc
634 real,
intent(in) :: f2
635 real,
intent(in) :: surface_pres
636 real,
dimension(SZK_(GV)), &
638 real,
dimension(SZK_(GV)), &
640 real,
dimension(SZK_(GV)), &
642 real,
dimension(SZK_(GV)), &
644 real,
dimension(SZK_(GV)), &
646 real,
dimension(SZK_(GV)+1), &
647 intent(out) :: kappa_avg
648 real,
dimension(SZK_(GV)+1), &
649 intent(out) :: tke_avg
650 real,
intent(in) :: dt
657 real,
dimension(SZK_(GV)+1), &
658 optional,
intent(out) :: I_Ld2_1d
659 real,
dimension(SZK_(GV)+1), &
660 optional,
intent(out) :: dz_Int_1d
664 real,
dimension(nzc) :: &
673 real,
dimension(nzc+1) :: &
708 real :: dist_from_bot
720 real :: tol_dksrc_low
734 logical :: use_temperature
736 integer :: ks_kappa, ke_kappa
737 integer :: dt_halvings
740 integer :: dt_refinements
743 integer :: k, itt, itt_dt
749 ri_crit = cs%Rino_crit
750 gr0 = gv%Rho0 * gv%g_Earth
751 g_r0 = (us%L_to_Z**2 * gv%g_Earth) / (gv%Rho0)
754 tol_dksrc = cs%kappa_src_max_chg
755 if (tol_dksrc == 10.0)
then
759 tol_dksrc_low = (tol_dksrc - 0.5)/tol_dksrc
761 tol2 = 2.0*cs%kappa_tol_err
763 use_temperature = .false. ;
if (
associated(tv%T)) use_temperature = .true.
767 do k=1,nzc ; idz(k) = 1.0 / dz(k) ;
enddo
770 i_dz_int(1) = 2.0 / dz(1)
771 dist_from_top(1) = 0.0
773 i_dz_int(k) = 2.0 / (dz(k-1) + dz(k))
774 dist_from_top(k) = dist_from_top(k-1) + dz(k-1)
776 i_dz_int(nzc+1) = 2.0 / dz(nzc)
781 a1(2) = k0dt*i_dz_int(2)
782 b1 = 1.0 / (dz(1) + a1(2))
783 u(1) = b1 * u0xdz(1) ; v(1) = b1 * v0xdz(1)
784 t(1) = b1 * t0xdz(1) ; sal(1) = b1 * s0xdz(1)
785 c1(2) = a1(2) * b1 ; d1 = dz(1) * b1
787 bd1 = dz(k) + d1*a1(k)
788 a1(k+1) = k0dt*i_dz_int(k+1)
789 b1 = 1.0 / (bd1 + a1(k+1))
790 u(k) = b1 * (u0xdz(k) + a1(k)*u(k-1))
791 v(k) = b1 * (v0xdz(k) + a1(k)*v(k-1))
792 t(k) = b1 * (t0xdz(k) + a1(k)*t(k-1))
793 sal(k) = b1 * (s0xdz(k) + a1(k)*sal(k-1))
794 c1(k+1) = a1(k+1) * b1 ; d1 = bd1 * b1
799 b1 = 1.0 / ((dz(nzc) + d1*a1(nzc)) + k0dt*i_dz_int(nzc+1))
800 u(nzc) = b1 * (u0xdz(nzc) + a1(nzc)*u(nzc-1))
801 v(nzc) = b1 * (v0xdz(nzc) + a1(nzc)*v(nzc-1))
803 b1 = 1.0 / (dz(nzc) + d1*a1(nzc))
804 t(nzc) = b1 * (t0xdz(nzc) + a1(nzc)*t(nzc-1))
805 sal(nzc) = b1 * (s0xdz(nzc) + a1(nzc)*sal(nzc-1))
807 u(k) = u(k) + c1(k+1)*u(k+1) ; v(k) = v(k) + c1(k+1)*v(k+1)
808 t(k) = t(k) + c1(k+1)*t(k+1) ; sal(k) = sal(k) + c1(k+1)*sal(k+1)
812 b1 = 1.0 / (dz(1) + k0dt*i_dz_int(2))
813 u(1) = b1 * u0xdz(1) ; v(1) = b1 * v0xdz(1)
815 t(1) = b1 * t0xdz(1) ; sal(1) = b1 * s0xdz(1)
824 dz_int(1) = 0.0 ; dz_int(2) = dz(1)
826 norm = 1.0 / (dz(k)*(dz(k-1)+dz(k+1)) + 2.0*dz(k-1)*dz(k+1))
827 dz_int(k) = dz_int(k) + dz(k) * ( ((dz(k)+dz(k+1)) * dz(k-1)) * norm)
828 dz_int(k+1) = dz(k) * ( ((dz(k-1)+dz(k)) * dz(k+1)) * norm)
830 dz_int(nzc) = dz_int(nzc) + dz(nzc) ; dz_int(nzc+1) = 0.0
834 dist_from_bot = dist_from_bot + dz(k)
835 i_l2_bdry(k) = (dist_from_top(k) + dist_from_bot)**2 / &
836 (dist_from_top(k) * dist_from_bot)**2
840 if (use_temperature)
then
841 pressure(1) = surface_pres
843 pressure(k) = pressure(k-1) + gr0*dz(k-1)
844 t_int(k) = 0.5*(t(k-1) + t(k))
845 sal_int(k) = 0.5*(sal(k-1) + sal(k))
848 tv%eqn_of_state, (/2,nzc/), scale=-g_r0 )
850 do k=1,nzc+1 ; dbuoy_dt(k) = -g_r0 ; dbuoy_ds(k) = 0.0 ;
enddo
861 call calculate_projected_state(kappa, u, v, t, sal, 0.0, nzc, dz, i_dz_int, &
862 dbuoy_dt, dbuoy_ds, u, v, t, sal, gv, us, &
863 n2=n2, s2=s2, vel_underflow=cs%vel_underflow)
870 kappa_avg(k) = 0.0 ; tke_avg(k) = 0.0
871 local_src_avg(k) = 0.0
873 if ( dz_int(k) > 0.0 ) &
874 local_src_avg(k) = 0.1 * k0dt * i_dz_int(k) / dz_int(k)
880 do itt=1,cs%max_KS_it
887 call find_kappa_tke(n2, s2, kappa, idz, dz_int, i_l2_bdry, f2, &
888 nzc, cs, gv, us, k_q, tke, kappa_out, kappa_src, local_src)
893 ks_kappa = gv%ke+1 ; ke_kappa = 0
894 do k=2,nzc ;
if (kappa_out(k) > 0.0)
then
897 do k=nzc,ks_kappa,-1 ;
if (kappa_out(k) > 0.0)
then
900 if (ke_kappa == nzc) kappa_out(nzc+1) = 0.0
906 if ((ke_kappa < ks_kappa) .or. (itt==cs%max_KS_it))
then
912 tol_max(k) = kappa_src(k) + tol_dksrc * local_src(k)
913 tol_min(k) = kappa_src(k) - tol_dksrc_low * local_src(k)
914 tol_chg(k) = tol2 * local_src_avg(k)
917 do itt_dt=1,(cs%max_KS_it+1-itt)/2
925 call calculate_projected_state(kappa_out, u, v, t, sal, 0.5*dt_test, nzc, dz, i_dz_int, &
926 dbuoy_dt, dbuoy_ds, u_test, v_test, t_test, s_test, &
927 gv, us, n2, s2, ks_int=ks_kappa, ke_int=ke_kappa, &
928 vel_underflow=cs%vel_underflow)
931 do k=max(ks_kappa-1,2),min(ke_kappa+1,nzc)
932 if (n2(k) < ri_crit * s2(k))
then
933 k_src(k) = (2.0 * cs%Shearmix_rate * sqrt(s2(k))) * &
934 ((ri_crit*s2(k) - n2(k)) / (ri_crit*s2(k) + cs%FRi_curvature*n2(k)))
935 if (cs%restrictive_tolerance_check)
then
936 if ((k_src(k) > min(tol_max(k), kappa_src(k) + idtt*tol_chg(k))) .or. &
937 (k_src(k) < max(tol_min(k), kappa_src(k) - idtt*tol_chg(k))))
then
938 valid_dt = .false. ;
exit
941 if ((k_src(k) > max(tol_max(k), kappa_src(k) + idtt*tol_chg(k))) .or. &
942 (k_src(k) < min(tol_min(k), kappa_src(k) - idtt*tol_chg(k))))
then
943 valid_dt = .false. ;
exit
947 if (0.0 < min(tol_min(k), kappa_src(k) - idtt*tol_chg(k)))
then
948 valid_dt = .false. ; k_src(k) = 0.0 ;
exit
954 dt_test = 0.5*dt_test
956 if ((dt_test < dt_rem) .and. valid_dt)
then
958 do itt_dt=1,dt_refinements
959 call calculate_projected_state(kappa_out, u, v, t, sal, 0.5*(dt_test+dt_inc), &
960 nzc, dz, i_dz_int, dbuoy_dt, dbuoy_ds, u_test, v_test, t_test, s_test, &
961 gv, us, n2, s2, ks_int=ks_kappa, ke_int=ke_kappa, vel_underflow=cs%vel_underflow)
963 idtt = 1.0 / (dt_test+dt_inc)
964 do k=max(ks_kappa-1,2),min(ke_kappa+1,nzc)
965 if (n2(k) < ri_crit * s2(k))
then
966 k_src(k) = (2.0 * cs%Shearmix_rate * sqrt(s2(k))) * &
967 ((ri_crit*s2(k) - n2(k)) / (ri_crit*s2(k) + cs%FRi_curvature*n2(k)))
968 if ((k_src(k) > max(tol_max(k), kappa_src(k) + idtt*tol_chg(k))) .or. &
969 (k_src(k) < min(tol_min(k), kappa_src(k) - idtt*tol_chg(k))))
then
970 valid_dt = .false. ;
exit
973 if (0.0 < min(tol_min(k), kappa_src(k) - idtt*tol_chg(k)))
then
974 valid_dt = .false. ; k_src(k) = 0.0 ;
exit
979 if (valid_dt) dt_test = dt_test + dt_inc
986 dt_now = min(dt_test*(1.0+cs%kappa_tol_err)+dt_inc, dt_rem)
988 local_src_avg(k) = local_src_avg(k) + dt_now * local_src(k)
995 if (ke_kappa < ks_kappa)
then
998 dt_wt = dt_rem / dt ; dt_rem = 0.0
1003 tke_avg(k) = tke_avg(k) + dt_wt*tke(k)
1008 call calculate_projected_state(kappa_out, u, v, t, sal, dt_now, nzc, dz, i_dz_int, &
1009 dbuoy_dt, dbuoy_ds, u_test, v_test, t_test, s_test, &
1010 gv, us, n2=n2, s2=s2, ks_int=ks_kappa, ke_int=ke_kappa, &
1011 vel_underflow=cs%vel_underflow)
1015 do k=1,nzc+1 ; k_q_tmp(k) = k_q(k) ;
enddo
1016 call find_kappa_tke(n2, s2, kappa_out, idz, dz_int, i_l2_bdry, f2, &
1017 nzc, cs, gv, us, k_q_tmp, tke_pred, kappa_pred)
1020 ks_kappa = gv%ke+1 ; ke_kappa = 0
1022 kappa_mid(k) = 0.5*(kappa_out(k) + kappa_pred(k))
1023 if ((kappa_mid(k) > 0.0) .and. (k<ks_kappa)) ks_kappa = k
1024 if (kappa_mid(k) > 0.0) ke_kappa = k
1028 call calculate_projected_state(kappa_mid, u, v, t, sal, dt_now, nzc, dz, i_dz_int, &
1029 dbuoy_dt, dbuoy_ds, u_test, v_test, t_test, s_test, &
1030 gv, us, n2=n2, s2=s2, ks_int=ks_kappa, ke_int=ke_kappa, &
1031 vel_underflow=cs%vel_underflow)
1035 call find_kappa_tke(n2, s2, kappa_out, idz, dz_int, i_l2_bdry, f2, &
1036 nzc, cs, gv, us, k_q, tke_pred, kappa_pred)
1040 dt_wt = dt_now / dt ; dt_rem = dt_rem - dt_now
1042 kappa_mid(k) = 0.5*(kappa_out(k) + kappa_pred(k))
1043 kappa_avg(k) = kappa_avg(k) + kappa_mid(k)*dt_wt
1044 tke_avg(k) = tke_avg(k) + dt_wt*0.5*(tke_pred(k) + tke(k))
1045 kappa(k) = kappa_pred(k)
1050 if (dt_rem > 0.0)
then
1053 call calculate_projected_state(kappa_mid, u, v, t, sal, dt_now, nzc, &
1054 dz, i_dz_int, dbuoy_dt, dbuoy_ds, u, v, t, sal, &
1055 gv, us, n2, s2, vel_underflow=cs%vel_underflow)
1059 if (dt_rem <= 0.0)
exit
1063 if (
present(i_ld2_1d))
then
1064 do k=1,gv%ke+1 ; i_ld2_1d(k) = 0.0 ;
enddo
1065 do k=2,nzc ;
if (tke(k) > 0.0) &
1066 i_ld2_1d(k) = i_l2_bdry(k) + (n2(k) / cs%lambda**2 + f2) / tke(k)
1069 if (
present(dz_int_1d))
then
1070 do k=1,nzc+1 ; dz_int_1d(k) = dz_int(k) ;
enddo
1071 do k=nzc+2,gv%ke ; dz_int_1d(k) = 0.0 ;
enddo
1074 end subroutine kappa_shear_column
1079 subroutine calculate_projected_state(kappa, u0, v0, T0, S0, dt, nz, &
1080 dz, I_dz_int, dbuoy_dT, dbuoy_dS, &
1081 u, v, T, Sal, GV, US, N2, S2, ks_int, ke_int, vel_underflow)
1082 integer,
intent(in) :: nz
1084 real,
dimension(nz+1),
intent(in) :: kappa
1086 real,
dimension(nz),
intent(in) :: u0
1087 real,
dimension(nz),
intent(in) :: v0
1088 real,
dimension(nz),
intent(in) :: T0
1089 real,
dimension(nz),
intent(in) :: S0
1090 real,
dimension(nz),
intent(in) :: dz
1091 real,
dimension(nz+1),
intent(in) :: I_dz_int
1093 real,
dimension(nz+1),
intent(in) :: dbuoy_dT
1095 real,
dimension(nz+1),
intent(in) :: dbuoy_dS
1097 real,
intent(in) :: dt
1098 real,
dimension(nz),
intent(inout) :: u
1099 real,
dimension(nz),
intent(inout) :: v
1100 real,
dimension(nz),
intent(inout) :: T
1101 real,
dimension(nz),
intent(inout) :: Sal
1104 real,
dimension(nz+1),
optional, &
1106 real,
dimension(nz+1),
optional, &
1108 integer,
optional,
intent(in) :: ks_int
1109 integer,
optional,
intent(in) :: ke_int
1111 real,
optional,
intent(in) :: vel_underflow
1116 real,
dimension(nz+1) :: c1
1119 real :: underflow_vel
1120 real :: a_a, a_b, b1, d1, bd1, b1nz_0
1121 integer :: k, ks, ke
1124 if (
present(ks_int)) ks = max(ks_int-1,1)
1125 if (
present(ke_int)) ke = min(ke_int,nz)
1126 underflow_vel = 0.0 ;
if (
present(vel_underflow)) underflow_vel = vel_underflow
1131 a_b = dt*(kappa(ks+1)*i_dz_int(ks+1))
1132 b1 = 1.0 / (dz(ks) + a_b)
1133 c1(ks+1) = a_b * b1 ; d1 = dz(ks) * b1
1135 u(ks) = (b1 * dz(ks))*u0(ks) ; v(ks) = (b1 * dz(ks))*v0(ks)
1136 t(ks) = (b1 * dz(ks))*t0(ks) ; sal(ks) = (b1 * dz(ks))*s0(ks)
1139 a_b = dt*(kappa(k+1)*i_dz_int(k+1))
1140 bd1 = dz(k) + d1*a_a
1141 b1 = 1.0 / (bd1 + a_b)
1142 c1(k+1) = a_b * b1 ; d1 = bd1 * b1
1144 u(k) = b1 * (dz(k)*u0(k) + a_a*u(k-1))
1145 v(k) = b1 * (dz(k)*v0(k) + a_a*v(k-1))
1146 t(k) = b1 * (dz(k)*t0(k) + a_a*t(k-1))
1147 sal(k) = b1 * (dz(k)*s0(k) + a_a*sal(k-1))
1154 b1 = 1.0 / (dz(ke) + d1*a_a)
1155 t(ke) = b1 * (dz(ke)*t0(ke) + a_a*t(ke-1))
1156 sal(ke) = b1 * (dz(ke)*s0(ke) + a_a*sal(ke-1))
1162 a_b = dt*(kappa(nz+1)*i_dz_int(nz+1))
1163 b1nz_0 = 1.0 / ((dz(nz) + d1*a_a) + a_b)
1167 u(ke) = b1nz_0 * (dz(ke)*u0(ke) + a_a*u(ke-1))
1168 v(ke) = b1nz_0 * (dz(ke)*v0(ke) + a_a*v(ke-1))
1169 if (abs(u(ke)) < underflow_vel) u(ke) = 0.0
1170 if (abs(v(ke)) < underflow_vel) v(ke) = 0.0
1173 u(k) = u(k) + c1(k+1)*u(k+1)
1174 v(k) = v(k) + c1(k+1)*v(k+1)
1175 if (abs(u(k)) < underflow_vel) u(k) = 0.0
1176 if (abs(v(k)) < underflow_vel) v(k) = 0.0
1177 t(k) = t(k) + c1(k+1)*t(k+1)
1178 sal(k) = sal(k) + c1(k+1)*sal(k+1)
1182 u(k) = u0(k) ; v(k) = v0(k) ; t(k) = t0(k) ; sal(k) = s0(k)
1183 if (abs(u(k)) < underflow_vel) u(k) = 0.0
1184 if (abs(v(k)) < underflow_vel) v(k) = 0.0
1188 if (
present(s2))
then
1190 l2_to_z2 = us%L_to_Z**2
1191 s2(1) = 0.0 ; s2(nz+1) = 0.0
1193 s2(ks) = ((u(ks)-u0(ks-1))**2 + (v(ks)-v0(ks-1))**2) * (l2_to_z2*i_dz_int(ks)**2)
1195 s2(k) = ((u(k)-u(k-1))**2 + (v(k)-v(k-1))**2) * (l2_to_z2*i_dz_int(k)**2)
1198 s2(ke+1) = ((u0(ke+1)-u(ke))**2 + (v0(ke+1)-v(ke))**2) * (l2_to_z2*i_dz_int(ke+1)**2)
1201 if (
present(n2))
then
1202 n2(1) = 0.0 ; n2(nz+1) = 0.0
1204 n2(ks) = max(0.0, i_dz_int(ks) * &
1205 (dbuoy_dt(ks) * (t0(ks-1)-t(ks)) + dbuoy_ds(ks) * (s0(ks-1)-sal(ks))))
1207 n2(k) = max(0.0, i_dz_int(k) * &
1208 (dbuoy_dt(k) * (t(k-1)-t(k)) + dbuoy_ds(k) * (sal(k-1)-sal(k))))
1211 n2(ke+1) = max(0.0, i_dz_int(ke+1) * &
1212 (dbuoy_dt(ke+1) * (t(ke)-t0(ke+1)) + dbuoy_ds(ke+1) * (sal(ke)-s0(ke+1))))
1215 end subroutine calculate_projected_state
1218 subroutine find_kappa_tke(N2, S2, kappa_in, Idz, dz_Int, I_L2_bdry, f2, &
1219 nz, CS, GV, US, K_Q, tke, kappa, kappa_src, local_src)
1220 integer,
intent(in) :: nz
1221 real,
dimension(nz+1),
intent(in) :: N2
1222 real,
dimension(nz+1),
intent(in) :: S2
1223 real,
dimension(nz+1),
intent(in) :: kappa_in
1225 real,
dimension(nz+1),
intent(in) :: dz_Int
1227 real,
dimension(nz+1),
intent(in) :: I_L2_bdry
1229 real,
dimension(nz),
intent(in) :: Idz
1230 real,
intent(in) :: f2
1234 real,
dimension(nz+1),
intent(inout) :: K_Q
1237 real,
dimension(nz+1),
intent(out) :: tke
1239 real,
dimension(nz+1),
intent(out) :: kappa
1241 real,
dimension(nz+1),
optional, &
1242 intent(out) :: kappa_src
1243 real,
dimension(nz+1),
optional, &
1244 intent(out) :: local_src
1249 real,
dimension(nz) :: &
1252 real,
dimension(nz+1) :: &
1273 real :: cQcomp, cKcomp
1287 real :: eden1, eden2, I_eden, ome
1288 real :: diffusive_src
1296 real :: decay_term_k
1297 real :: decay_term_Q
1307 real,
parameter :: roundoff = 1.0e-16
1310 logical :: tke_noflux_bottom_BC = .false.
1311 logical :: tke_noflux_top_BC = .false.
1312 logical :: do_Newton
1313 logical :: abort_Newton
1315 logical :: was_Newton
1316 logical :: within_tolerance
1318 integer :: ks_src, ke_src
1319 integer :: ks_kappa, ke_kappa, ke_tke
1320 integer :: ks_kappa_prev, ke_kappa_prev
1321 integer :: itt, k, k2
1324 logical,
parameter :: debug_soln = .false.
1325 real :: K_err_lin, Q_err_lin, TKE_src_norm
1326 real,
dimension(nz+1) :: &
1331 c_n2 = cs%C_N**2 ; c_s2 = cs%C_S**2
1332 q0 = cs%TKE_bg ; kappa0 = cs%kappa_0
1333 tke_min = max(cs%TKE_bg, 1.0e-20*us%m_to_Z**2*us%T_to_s**2)
1334 ri_crit = cs%Rino_crit
1335 ilambda2 = 1.0 / cs%lambda**2
1336 kappa_trunc = cs%kappa_trunc
1337 do_newton = .false. ; abort_newton = .false.
1338 tol_err = cs%kappa_tol_err
1341 ks_kappa = 2 ; ke_kappa = nz ; ks_kappa_prev = 2 ; ke_kappa_prev = nz
1343 ke_src = 0 ; ks_src = nz+1
1345 if (n2(k) < ri_crit * s2(k))
then
1349 k_src(k) = (2.0 * cs%Shearmix_rate * sqrt(s2(k))) * &
1350 ((ri_crit*s2(k) - n2(k)) / (ri_crit*s2(k) + cs%FRi_curvature*n2(k)))
1352 if (ks_src > k) ks_src = k
1359 if (ks_src > ke_src)
then
1361 kappa(k) = 0.0 ; k_q(k) = 0.0 ; tke(k) = tke_min
1363 if (
present(kappa_src))
then ;
do k=1,nz+1 ; kappa_src(k) = 0.0 ;
enddo ;
endif
1364 if (
present(local_src))
then ;
do k=1,nz+1 ; local_src(k) = 0.0 ;
enddo ;
endif
1369 kappa(k) = kappa_in(k)
1371 tke_decay(k) = sqrt(c_n2*n2(k) + c_s2*s2(k))
1372 if ((kappa(k) > 0.0) .and. (k_q(k) > 0.0))
then
1373 tke(k) = kappa(k) / k_q(k)
1379 kappa(1) = 0.0 ; kappa(nz+1) = 0.0
1384 eden2 = kappa0 * idz(nz)
1385 if (tke_noflux_bottom_bc)
then
1386 eden1 = dz_int(nz+1)*tke_decay(nz+1)
1387 i_eden = 1.0 / (eden2 + eden1)
1388 e1(nz+1) = eden2 * i_eden ; ome = eden1 * i_eden
1390 e1(nz+1) = 0.0 ; ome = 1.0
1393 eden1 = dz_int(k)*tke_decay(k) + ome * eden2
1394 eden2 = kappa0 * idz(k-1)
1395 i_eden = 1.0 / (eden2 + eden1)
1396 e1(k) = eden2 * i_eden ; ome = eden1 * i_eden
1402 do itt=1,cs%max_RiNo_it
1408 if (debug_soln)
then ;
do k=1,nz+1 ; kappa_prev(k) = kappa(k) ; tke_prev(k) = tke(k) ;
enddo ;
endif
1410 if (.not.do_newton)
then
1416 ke_tke = max(ke_kappa,ke_kappa_prev)+1
1418 do k=1,min(ke_tke,nz)
1419 aq(k) = (0.5*(kappa(k)+kappa(k+1)) + kappa0) * idz(k)
1422 if (tke_noflux_top_bc)
then
1423 tke_src = kappa0*s2(1) + q0 * tke_decay(1)
1424 bqd1 = dz_int(1) * tke_decay(1)
1425 bq = 1.0 / (bqd1 + aq(1))
1426 tke(1) = bq * (dz_int(1)*tke_src)
1427 cq(2) = aq(1) * bq ; cqcomp = bqd1 * bq
1429 tke(1) = q0 ; cq(2) = 0.0 ; cqcomp = 1.0
1433 tke_src = (kappa(k) + kappa0)*s2(k) + q0*tke_decay(k)
1434 bqd1 = dz_int(k)*(tke_decay(k) + n2(k)*k_q(k)) + cqcomp*aq(k-1)
1435 bq = 1.0 / (bqd1 + aq(k))
1436 tke(k) = bq * (dz_int(k)*tke_src + aq(k-1)*tke(k-1))
1437 cq(k+1) = aq(k) * bq ; cqcomp = bqd1 * bq
1439 if ((ke_tke == nz+1) .and. .not.(tke_noflux_bottom_bc))
then
1444 tke_src = kappa0*s2(k) + q0*tke_decay(k)
1447 bq = 1.0 / (dz_int(k)*tke_decay(k) + cqcomp*aq(k-1))
1448 tke(k) = max(tke_min, bq * (dz_int(k)*tke_src + aq(k-1)*tke(k-1)))
1449 dq(k) = tke(k) + dq(k)
1451 bq = 1.0 / ((dz_int(k)*tke_decay(k) + cqcomp*aq(k-1)) + aq(k))
1452 cq(k+1) = aq(k) * bq
1455 tke(k) = max((bq * (dz_int(k)*tke_src + aq(k-1)*tke(k-1)) + &
1456 cq(k+1)*(tke(k+1) - e1(k+1)*tke(k))) / (1.0 - cq(k+1)*e1(k+1)), tke_min)
1457 dq(k) = tke(k) + dq(k)
1462 dq(k) = e1(k)*dq(k-1)
1463 tke(k) = max(tke(k) + dq(k), tke_min)
1464 if (abs(dq(k)) < roundoff*tke(k))
exit
1467 if (dq(k2) == 0.0)
exit
1473 tke(k) = max(tke(k) + cq(k+1)*tke(k+1), tke_min)
1474 dq(k) = tke(k) + dq(k)
1481 ke_kappa_prev = ke_kappa ; ks_kappa_prev = ks_kappa
1484 ck(2) = 0.0 ; ckcomp = 1.0
1485 if (itt == 1)
then ;
dO k=2,nz
1486 i_ld2(k) = (n2(k)*ilambda2 + f2) / tke(k) + i_l2_bdry(k)
1491 i_ld2(k) = (n2(k)*ilambda2 + f2) / tke(k) + i_l2_bdry(k)
1492 bkd1 = dz_int(k)*i_ld2(k) + ckcomp*idz(k-1)
1493 bk = 1.0 / (bkd1 + idz(k))
1495 kappa(k) = bk * (idz(k-1)*kappa(k-1) + dz_int(k) * k_src(k))
1496 ck(k+1) = idz(k) * bk ; ckcomp = bkd1 * bk
1499 if (kappa(k) < ckcomp*kappa_trunc)
then
1501 if (k > ke_src)
then ; ke_kappa = k-1 ; k_q(k) = 0.0 ;
exit ;
endif
1502 elseif (kappa(k) < 2.0*ckcomp*kappa_trunc)
then
1503 kappa(k) = 2.0 * (kappa(k) - ckcomp*kappa_trunc)
1506 k_q(ke_kappa) = kappa(ke_kappa) / tke(ke_kappa)
1507 dk(ke_kappa) = dk(ke_kappa) + kappa(ke_kappa)
1508 do k=ke_kappa+2,ke_kappa_prev
1509 dk(k) = -kappa(k) ; kappa(k) = 0.0 ; k_q(k) = 0.0
1511 do k=ke_kappa-1,2,-1
1512 kappa(k) = kappa(k) + ck(k+1)*kappa(k+1)
1514 if (kappa(k) <= kappa_trunc)
then
1516 if (k < ks_src)
then ; ks_kappa = k+1 ; k_q(k) = 0.0 ;
exit ;
endif
1517 elseif (kappa(k) < 2.0*kappa_trunc)
then
1518 kappa(k) = 2.0 * (kappa(k) - kappa_trunc)
1521 dk(k) = dk(k) + kappa(k)
1522 k_q(k) = kappa(k) / tke(k)
1524 do k=ks_kappa_prev,ks_kappa-2 ; kappa(k) = 0.0 ; k_q(k) = 0.0 ;
enddo
1529 ks_kappa_prev = ks_kappa ; ke_kappa_prev = ke_kappa ; ke_kappa = nz
1531 dk(1) = 0.0 ; ck(2) = 0.0 ; ckcomp = 1.0 ; dkdq(1) = 0.0
1532 aq(1) = (0.5*(kappa(1)+kappa(2))+kappa0) * idz(1)
1533 dqdz(1) = 0.5*(tke(1) - tke(2))*idz(1)
1534 if (tke_noflux_top_bc)
then
1535 tke_src = dz_int(1) * (kappa0*s2(1) - (tke(1) - q0)*tke_decay(1)) - &
1536 aq(1) * (tke(1) - tke(2))
1538 bq = 1.0 / (aq(1) + dz_int(1)*tke_decay(1))
1540 cqcomp = (dz_int(1)*tke_decay(1)) * bq
1541 dqmdk(2) = -dqdz(1) * bq
1542 dq(1) = bq * tke_src
1544 dq(1) = 0.0 ; cq(2) = 0.0 ; cqcomp = 1.0 ; dqmdk(2) = 0.0
1548 i_ld2(k) = (n2(k)*ilambda2 + f2) * i_q + i_l2_bdry(k)
1550 kap_src = dz_int(k) * (k_src(k) - i_ld2(k)*kappa(k)) + &
1551 idz(k-1)*(kappa(k-1)-kappa(k)) - idz(k)*(kappa(k)-kappa(k+1))
1555 decay_term_k = -idz(k-1)*dqmdk(k)*dkdq(k-1) + dz_int(k)*i_ld2(k)
1556 if (decay_term_k < 0.0)
then ; abort_newton = .true. ;
exit ;
endif
1557 bk = 1.0 / (idz(k) + idz(k-1)*ckcomp + decay_term_k)
1559 ck(k+1) = bk * idz(k)
1560 ckcomp = bk * (idz(k-1)*ckcomp + decay_term_k)
1561 if (cs%dKdQ_iteration_bug)
then
1562 dkdq(k) = bk * (idz(k-1)*dkdq(k-1)*cq(k) + &
1563 us%m_to_Z*(n2(k)*ilambda2 + f2) * i_q**2 * kappa(k) )
1565 dkdq(k) = bk * (idz(k-1)*dkdq(k-1)*cq(k) + &
1566 dz_int(k)*(n2(k)*ilambda2 + f2) * i_q**2 * kappa(k) )
1568 dk(k) = bk * (kap_src + idz(k-1)*dk(k-1) + idz(k-1)*dkdq(k-1)*dq(k-1))
1571 if (dk(k) <= ckcomp*(kappa_trunc - kappa(k)))
then
1572 dk(k) = -ckcomp*kappa(k)
1574 elseif (dk(k) < ckcomp*(2.0*kappa_trunc - kappa(k)))
then
1575 dk(k) = 2.0 * dk(k) - ckcomp*(2.0*kappa_trunc - kappa(k))
1579 aq(k) = (0.5*(kappa(k)+kappa(k+1))+kappa0) * idz(k)
1580 dqdz(k) = 0.5*(tke(k) - tke(k+1))*idz(k)
1581 tke_src = dz_int(k) * (((kappa(k) + kappa0)*s2(k) - kappa(k)*n2(k)) - &
1582 (tke(k) - q0)*tke_decay(k)) - &
1583 (aq(k) * (tke(k) - tke(k+1)) - aq(k-1) * (tke(k-1) - tke(k)))
1584 v1 = aq(k-1) + dqdz(k-1)*dkdq(k-1)
1585 v2 = (v1*dqmdk(k) + dqdz(k-1)*ck(k)) + &
1586 ((dqdz(k-1) - dqdz(k)) + dz_int(k)*(s2(k) - n2(k)))
1590 decay_term_q = dz_int(k)*tke_decay(k) - dqdz(k-1)*dkdq(k-1)*cq(k) - v2*dkdq(k)
1591 if (decay_term_q < 0.0)
then ; abort_newton = .true. ;
exit ;
endif
1592 bq = 1.0 / (aq(k) + (cqcomp*aq(k-1) + decay_term_q))
1594 cq(k+1) = aq(k) * bq
1595 cqcomp = (cqcomp*aq(k-1) + decay_term_q) * bq
1596 dqmdk(k+1) = (v2 * ck(k+1) - dqdz(k)) * bq
1599 dq(k) = max(bq * ((v1 * dq(k-1) + dqdz(k-1)*dk(k-1)) + &
1600 (v2 * dk(k) + tke_src)), cqcomp*(-0.5*tke(k)))
1603 if ((itt > 1) .and. (k > ke_src) .and. (dk(k) == 0.0) .and. &
1604 ((kappa(k) + kappa(k+1)) == 0.0))
then
1606 ke_kappa = k-1 ;
exit
1609 if ((ke_kappa == nz) .and. (.not. abort_newton))
then
1610 dk(nz+1) = 0.0 ; dkdq(nz+1) = 0.0
1611 if (tke_noflux_bottom_bc)
then
1613 tke_src = dz_int(k) * (kappa0*s2(k) - (tke(k) - q0)*tke_decay(k)) + &
1614 aq(k-1) * (tke(k-1) - tke(k))
1616 v1 = aq(k-1) + dqdz(k-1)*dkdq(k-1)
1617 decay_term_q = max(0.0, dz_int(k)*tke_decay(k) - dqdz(k-1)*dkdq(k-1)*cq(k))
1618 if (decay_term_q < 0.0)
then
1619 abort_newton = .true.
1621 bq = 1.0 / (aq(k) + (cqcomp*aq(k-1) + decay_term_q))
1623 dq(k) = max(bq * ((v1 * dq(k-1) + dqdz(k-1)*dk(k-1)) + tke_src), -0.5*tke(k))
1624 tke(k) = max(tke(k) + dq(k), tke_min)
1629 elseif (.not. abort_newton)
then
1631 dq(ke_kappa+1) = dq(ke_kappa+1) / (1.0 - cq(ke_kappa+2)*e1(ke_kappa+2))
1632 tke(ke_kappa+1) = max(tke(ke_kappa+1) + dq(ke_kappa+1), tke_min)
1633 do k=ke_kappa+2,nz+1
1634 if (debug_soln .and. (k < nz+1))
then
1636 aq(k) = (0.5*(kappa(k)+kappa(k+1))+kappa0) * idz(k)
1643 dq(k) = max(e1(k)*dq(k-1),-0.5*tke(k))
1644 tke(k) = max(tke(k) + dq(k), tke_min)
1645 if (abs(dq(k)) < roundoff*tke(k))
exit
1647 if (debug_soln)
then ;
do k2=k+1,nz+1 ; dq(k2) = 0.0 ; dk(k2) = 0.0 ;
enddo ;
endif
1649 if (.not. abort_newton)
then
1652 dq(k) = max(dq(k) + (cq(k+1)*dq(k+1) + dqmdk(k+1) * dk(k+1)), -0.5*tke(k))
1653 tke(k) = max(tke(k) + dq(k), tke_min)
1654 dk(k) = dk(k) + (ck(k+1)*dk(k+1) + dkdq(k) * dq(k))
1656 if (dk(k) <= kappa_trunc - kappa(k))
then
1659 if ((k < ks_src) .and. (k+1 > ks_kappa)) ks_kappa = k+1
1660 elseif (dk(k) < 2.0*kappa_trunc - kappa(k))
then
1661 dk(k) = 2.0*dk(k) - (2.0*kappa_trunc - kappa(k))
1662 kappa(k) = max(kappa(k) + dk(k), 0.0)
1663 if (k<=ks_kappa) ks_kappa = 2
1665 kappa(k) = kappa(k) + dk(k)
1666 if (k<=ks_kappa) ks_kappa = 2
1669 dq(1) = max(dq(1) + cq(2)*dq(2) + dqmdk(2) * dk(2), tke_min - tke(1))
1670 tke(1) = max(tke(1) + dq(1), tke_min)
1676 if (debug_soln)
then ;
do k=2,nz
1683 i_ld2_debug(k) = (n2(k)*ilambda2 + f2) * i_q + i_l2_bdry(k)
1685 kap_src = dz_int(k) * (k_src(k) - i_ld2(k)*kappa_prev(k)) + &
1686 (idz(k-1)*(kappa_prev(k-1)-kappa_prev(k)) - &
1687 idz(k)*(kappa_prev(k)-kappa_prev(k+1)))
1688 k_err_lin = -idz(k-1)*(dk(k-1)-dk(k)) + idz(k)*(dk(k)-dk(k+1)) + &
1689 dz_int(k)*i_ld2_debug(k)*dk(k) - kap_src - &
1690 dz_int(k)*(n2(k)*ilambda2 + f2)*i_q**2*kappa_prev(k) * dq(k)
1692 tke_src = dz_int(k) * ((kappa_prev(k) + kappa0)*s2(k) - &
1693 kappa_prev(k)*n2(k) - (tke_prev(k) - q0)*tke_decay(k)) - &
1694 (aq(k) * (tke_prev(k) - tke_prev(k+1)) - aq(k-1) * (tke_prev(k-1) - tke_prev(k)))
1695 q_err_lin = tke_src + (aq(k-1) * (dq(k-1)-dq(k)) - aq(k) * (dq(k)-dq(k+1))) - &
1696 0.5*(tke_prev(k)-tke_prev(k+1))*idz(k) * (dk(k) + dk(k+1)) - &
1697 0.5*(tke_prev(k)-tke_prev(k-1))*idz(k-1)* (dk(k-1) + dk(k)) + &
1698 dz_int(k) * (dk(k) * (s2(k) - n2(k)) - dq(k)*tke_decay(k))
1704 if ((tol_err < newton_err) .and. (.not.abort_newton))
then
1706 newton_test = newton_err ;
if (do_newton) newton_test = 2.0*newton_err
1707 was_newton = do_newton
1708 within_tolerance = .true. ; do_newton = .true.
1709 do k=min(ks_kappa,ks_kappa_prev),max(ke_kappa,ke_kappa_prev)
1710 kappa_mean = kappa0 + (kappa(k) - 0.5*dk(k))
1711 if (abs(dk(k)) > newton_test * kappa_mean)
then
1712 if (do_newton) abort_newton = .true.
1713 within_tolerance = .false. ; do_newton = .false. ;
exit
1714 elseif (abs(dk(k)) > tol_err * kappa_mean)
then
1715 within_tolerance = .false. ;
if (.not.do_newton)
exit
1717 if (abs(dq(k)) > newton_test*(tke(k) - 0.5*dq(k)))
then
1718 if (do_newton) abort_newton = .true.
1719 do_newton = .false. ;
if (.not.within_tolerance)
exit
1724 within_tolerance = .true.
1725 do k=min(ks_kappa,ks_kappa_prev),max(ke_kappa,ke_kappa_prev)
1726 if (abs(dk(k)) > tol_err * (kappa0 + (kappa(k) - 0.5*dk(k))))
then
1727 within_tolerance = .false. ;
exit
1732 if (abort_newton)
then
1733 do_newton = .false. ; abort_newton = .false.
1735 newton_err = 0.5*newton_err
1737 do k=2,nz ; k_q(k) = kappa(k) / max(tke(k), tke_min) ;
enddo
1740 if (within_tolerance)
exit
1745 do k=1,ks_kappa-1 ; k_q(k) = 0.0 ;
enddo
1746 do k=ks_kappa,ke_kappa ; k_q(k) = kappa(k) / tke(k) ;
enddo
1747 do k=ke_kappa+1,nz+1 ; k_q(k) = 0.0 ;
enddo
1750 if (
present(local_src))
then
1751 local_src(1) = 0.0 ; local_src(nz+1) = 0.0
1753 diffusive_src = idz(k-1)*(kappa(k-1)-kappa(k)) + idz(k)*(kappa(k+1)-kappa(k))
1754 chg_by_k0 = kappa0 * ((idz(k-1)+idz(k)) / dz_int(k) + i_ld2(k))
1755 if (diffusive_src <= 0.0)
then
1756 local_src(k) = k_src(k) + chg_by_k0
1758 local_src(k) = (k_src(k) + chg_by_k0) + diffusive_src / dz_int(k)
1762 if (
present(kappa_src))
then
1763 kappa_src(1) = 0.0 ; kappa_src(nz+1) = 0.0
1765 kappa_src(k) = k_src(k)
1769 end subroutine find_kappa_tke
1772 function kappa_shear_init(Time, G, GV, US, param_file, diag, CS)
1773 type(time_type),
intent(in) :: time
1779 type(
diag_ctrl),
target,
intent(inout) :: diag
1783 logical :: kappa_shear_init
1786 logical :: merge_mixedlayer
1787 logical :: debug_shear
1788 logical :: just_read
1790 # include "version_variable.h"
1791 character(len=40) :: mdl =
"MOM_kappa_shear"
1792 real :: kappa_0_unscaled
1796 if (
associated(cs))
then
1797 call mom_error(warning,
"kappa_shear_init called with an associated "// &
1798 "control structure.")
1813 call get_param(param_file, mdl,
"USE_JACKSON_PARAM", kappa_shear_init, default=.false., do_not_log=.true.)
1815 "Parameterization of shear-driven turbulence following Jackson, Hallberg and Legg, JPO 2008", &
1816 log_to_all=.true., debugging=kappa_shear_init, all_default=.not.kappa_shear_init)
1817 call get_param(param_file, mdl,
"USE_JACKSON_PARAM", kappa_shear_init, &
1818 "If true, use the Jackson-Hallberg-Legg (JPO 2008) "//&
1819 "shear mixing parameterization.", default=.false.)
1820 just_read = .not.kappa_shear_init
1821 call get_param(param_file, mdl,
"VERTEX_SHEAR", cs%KS_at_vertex, &
1822 "If true, do the calculations of the shear-driven mixing "//&
1823 "at the cell vertices (i.e., the vorticity points).", &
1824 default=.false., do_not_log=just_read)
1825 call get_param(param_file, mdl,
"RINO_CRIT", cs%RiNo_crit, &
1826 "The critical Richardson number for shear mixing.", &
1827 units=
"nondim", default=0.25, do_not_log=just_read)
1828 call get_param(param_file, mdl,
"SHEARMIX_RATE", cs%Shearmix_rate, &
1829 "A nondimensional rate scale for shear-driven entrainment. "//&
1830 "Jackson et al find values in the range of 0.085-0.089.", &
1831 units=
"nondim", default=0.089, do_not_log=just_read)
1832 call get_param(param_file, mdl,
"MAX_RINO_IT", cs%max_RiNo_it, &
1833 "The maximum number of iterations that may be used to "//&
1834 "estimate the Richardson number driven mixing.", &
1835 units=
"nondim", default=50, do_not_log=just_read)
1836 call get_param(param_file, mdl,
"KD", kd_normal, default=0.0, do_not_log=.true.)
1837 call get_param(param_file, mdl,
"KD_KAPPA_SHEAR_0", cs%kappa_0, &
1838 "The background diffusivity that is used to smooth the "//&
1839 "density and shear profiles before solving for the "//&
1840 "diffusivities. The default is the greater of KD and 1e-7 m2 s-1.", &
1841 units=
"m2 s-1", default=max(kd_normal, 1.0e-7), scale=us%m2_s_to_Z2_T, &
1842 unscaled=kappa_0_unscaled, do_not_log=just_read)
1843 call get_param(param_file, mdl,
"KD_TRUNC_KAPPA_SHEAR", cs%kappa_trunc, &
1844 "The value of shear-driven diffusivity that is considered negligible "//&
1845 "and is rounded down to 0. The default is 1% of KD_KAPPA_SHEAR_0.", &
1846 units=
"m2 s-1", default=0.01*kappa_0_unscaled, scale=us%m2_s_to_Z2_T, do_not_log=just_read)
1847 call get_param(param_file, mdl,
"FRI_CURVATURE", cs%FRi_curvature, &
1848 "The nondimensional curvature of the function of the "//&
1849 "Richardson number in the kappa source term in the "//&
1850 "Jackson et al. scheme.", units=
"nondim", default=-0.97, do_not_log=just_read)
1851 call get_param(param_file, mdl,
"TKE_N_DECAY_CONST", cs%C_N, &
1852 "The coefficient for the decay of TKE due to "//&
1853 "stratification (i.e. proportional to N*tke). "//&
1854 "The values found by Jackson et al. are 0.24-0.28.", &
1855 units=
"nondim", default=0.24, do_not_log=just_read)
1858 call get_param(param_file, mdl,
"TKE_SHEAR_DECAY_CONST", cs%C_S, &
1859 "The coefficient for the decay of TKE due to shear (i.e. "//&
1860 "proportional to |S|*tke). The values found by Jackson "//&
1861 "et al. are 0.14-0.12.", units=
"nondim", default=0.14, do_not_log=just_read)
1862 call get_param(param_file, mdl,
"KAPPA_BUOY_SCALE_COEF", cs%lambda, &
1863 "The coefficient for the buoyancy length scale in the "//&
1864 "kappa equation. The values found by Jackson et al. are "//&
1865 "in the range of 0.81-0.86.", units=
"nondim", default=0.82, do_not_log=just_read)
1866 call get_param(param_file, mdl,
"KAPPA_N_OVER_S_SCALE_COEF2", cs%lambda2_N_S, &
1867 "The square of the ratio of the coefficients of the "//&
1868 "buoyancy and shear scales in the diffusivity equation, "//&
1869 "Set this to 0 (the default) to eliminate the shear scale. "//&
1870 "This is only used if USE_JACKSON_PARAM is true.", &
1871 units=
"nondim", default=0.0, do_not_log=just_read)
1872 call get_param(param_file, mdl,
"KAPPA_SHEAR_TOL_ERR", cs%kappa_tol_err, &
1873 "The fractional error in kappa that is tolerated. "//&
1874 "Iteration stops when changes between subsequent "//&
1875 "iterations are smaller than this everywhere in a "//&
1876 "column. The peak diffusivities usually converge most "//&
1877 "rapidly, and have much smaller errors than this.", &
1878 units=
"nondim", default=0.1, do_not_log=just_read)
1879 call get_param(param_file, mdl,
"TKE_BACKGROUND", cs%TKE_bg, &
1880 "A background level of TKE used in the first iteration "//&
1881 "of the kappa equation. TKE_BACKGROUND could be 0.", &
1882 units=
"m2 s-2", default=0.0, scale=us%m_to_Z**2*us%T_to_s**2)
1883 call get_param(param_file, mdl,
"KAPPA_SHEAR_ELIM_MASSLESS", cs%eliminate_massless, &
1884 "If true, massless layers are merged with neighboring "//&
1885 "massive layers in this calculation. The default is "//&
1886 "true and I can think of no good reason why it should "//&
1887 "be false. This is only used if USE_JACKSON_PARAM is true.", &
1888 default=.true., do_not_log=just_read)
1889 call get_param(param_file, mdl,
"MAX_KAPPA_SHEAR_IT", cs%max_KS_it, &
1890 "The maximum number of iterations that may be used to "//&
1891 "estimate the time-averaged diffusivity.", units=
"nondim", &
1892 default=13, do_not_log=just_read)
1893 call get_param(param_file, mdl,
"PRANDTL_TURB", cs%Prandtl_turb, &
1894 "The turbulent Prandtl number applied to shear instability.", &
1895 units=
"nondim", default=1.0, do_not_log=.true.)
1896 call get_param(param_file, mdl,
"VEL_UNDERFLOW", cs%vel_underflow, &
1897 "A negligibly small velocity magnitude below which velocity components are set "//&
1898 "to 0. A reasonable value might be 1e-30 m/s, which is less than an "//&
1899 "Angstrom divided by the age of the universe.", &
1900 units=
"m s-1", default=0.0, scale=us%m_s_to_L_T, do_not_log=just_read)
1901 call get_param(param_file, mdl,
"KAPPA_SHEAR_MAX_KAP_SRC_CHG", cs%kappa_src_max_chg, &
1902 "The maximum permitted increase in the kappa source within an iteration relative "//&
1903 "to the local source; this must be greater than 1. The lower limit for the "//&
1904 "permitted fractional decrease is (1 - 0.5/kappa_src_max_chg). These limits "//&
1905 "could perhaps be made dynamic with an improved iterative solver.", &
1906 default=10.0, units=
"nondim", do_not_log=just_read)
1908 call get_param(param_file, mdl,
"DEBUG", cs%debug, &
1909 "If true, write out verbose debugging data.", &
1910 default=.false., debuggingparam=.true., do_not_log=just_read)
1911 call get_param(param_file, mdl,
"DEBUG_KAPPA_SHEAR", debug_shear, &
1912 "If true, write debugging data for the kappa-shear code.", &
1913 default=.false., debuggingparam=.true., do_not_log=.true.)
1914 if (debug_shear) cs%debug = .true.
1915 call get_param(param_file, mdl,
"KAPPA_SHEAR_VERTEX_PSURF_BUG", cs%psurf_bug, &
1916 "If true, do a simple average of the cell surface pressures to get a pressure "//&
1917 "at the corner if VERTEX_SHEAR=True. Otherwise mask out any land points in "//&
1918 "the average.", default=.true., do_not_log=(just_read .or. (.not.cs%KS_at_vertex)))
1920 call get_param(param_file, mdl,
"KAPPA_SHEAR_ITER_BUG", cs%dKdQ_iteration_bug, &
1921 "If true, use an older, dimensionally inconsistent estimate of the "//&
1922 "derivative of diffusivity with energy in the Newton's method iteration. "//&
1923 "The bug causes undercorrections when dz > 1 m.", default=.false., do_not_log=just_read)
1924 call get_param(param_file, mdl,
"KAPPA_SHEAR_ALL_LAYER_TKE_BUG", cs%all_layer_TKE_bug, &
1925 "If true, report back the latest estimate of TKE instead of the time average "//&
1926 "TKE when there is mass in all layers. Otherwise always report the time "//&
1927 "averaged TKE, as is currently done when there are some massless layers.", &
1928 default=.false., do_not_log=just_read)
1929 call get_param(param_file, mdl,
"USE_RESTRICTIVE_TOLERANCE_CHECK", cs%restrictive_tolerance_check, &
1930 "If true, uses the more restrictive tolerance check to determine if a timestep "//&
1931 "is acceptable for the KS_it outer iteration loop. False uses the original less "//&
1932 "restrictive check.", default=.false., do_not_log=just_read)
1940 call get_param(param_file, mdl,
"KAPPA_SHEAR_MERGE_ML",merge_mixedlayer, &
1941 "If true, combine the mixed layers together before solving the "//&
1942 "kappa-shear equations.", default=.true., do_not_log=just_read)
1943 if (merge_mixedlayer) cs%nkml = gv%nkml
1947 if (.not. kappa_shear_init)
return
1951 cs%id_Kd_shear = register_diag_field(
'ocean_model',
'Kd_shear', diag%axesTi, time, &
1952 'Shear-driven Diapycnal Diffusivity',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
1953 cs%id_TKE = register_diag_field(
'ocean_model',
'TKE_shear', diag%axesTi, time, &
1954 'Shear-driven Turbulent Kinetic Energy',
'm2 s-2', conversion=us%Z_to_m**2*us%s_to_T**2)
1956 end function kappa_shear_init
1960 logical function kappa_shear_is_used(param_file)
1964 character(len=40) :: mdl =
"MOM_kappa_shear"
1967 call get_param(param_file, mdl,
"USE_JACKSON_PARAM", kappa_shear_is_used, &
1968 default=.false., do_not_log=.true.)
1969 end function kappa_shear_is_used
1974 logical function kappa_shear_at_vertex(param_file)
1978 character(len=40) :: mdl =
"MOM_kappa_shear"
1979 logical :: do_kappa_shear
1982 kappa_shear_at_vertex = .false.
1984 call get_param(param_file, mdl,
"USE_JACKSON_PARAM", do_kappa_shear, &
1985 default=.false., do_not_log=.true.)
1986 if (do_kappa_shear) &
1987 call get_param(param_file, mdl,
"VERTEX_SHEAR", kappa_shear_at_vertex, &
1988 "If true, do the calculations of the shear-driven mixing "//&
1989 "at the cell vertices (i.e., the vorticity points).", &
1990 default=.false., do_not_log=.true.)
1992 end function kappa_shear_at_vertex