6 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
10 use mom_domains,
only : sum_across_pes, max_across_pes
33 implicit none ;
private 35 #include <MOM_memory.h> 37 public tracer_hordiff, tracer_hor_diff_init, tracer_hor_diff_end
42 real :: khtr_slope_cff
45 real :: khtr_passivity_coeff
48 real :: khtr_passivity_min
55 logical :: diffuse_ml_interior
57 logical :: check_diffusive_cfl
60 logical :: use_neutral_diffusion
62 logical :: use_lateral_boundary_diffusion
64 logical :: recalc_neutral_surf
72 logical :: show_call_tree
73 logical :: first_call = .true.
75 integer :: id_khtr_u = -1
76 integer :: id_khtr_v = -1
77 integer :: id_khtr_h = -1
78 integer :: id_cfl = -1
79 integer :: id_khdt_x = -1
80 integer :: id_khdt_y = -1
83 type(group_pass_type) :: pass_t
89 real,
dimension(:,:),
pointer :: p => null()
93 integer,
dimension(:,:),
pointer :: p => null()
97 integer :: id_clock_diffuse, id_clock_epimix, id_clock_pass, id_clock_sync
106 subroutine tracer_hordiff(h, dt, MEKE, VarMix, G, GV, US, CS, Reg, tv, do_online_flag, read_khdt_x, read_khdt_y)
108 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
110 real,
intent(in) :: dt
124 logical,
optional,
intent(in) :: do_online_flag
126 real,
dimension(SZIB_(G),SZJ_(G)), &
127 optional,
intent(in) :: read_khdt_x
129 real,
dimension(SZI_(G),SZJB_(G)), &
130 optional,
intent(in) :: read_khdt_y
134 real,
dimension(SZI_(G),SZJ_(G)) :: &
141 real,
dimension(SZIB_(G),SZJ_(G)) :: &
147 real,
dimension(SZI_(G),SZJB_(G)) :: &
156 logical :: use_varmix, resoln_scaled, do_online, use_eady
157 integer :: s_idx, t_idx
158 integer :: i, j, k, m, is, ie, js, je, nz, ntr, itt, num_itts
170 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
173 if (
present(do_online_flag)) do_online = do_online_flag
175 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_tracer_hor_diff: "// &
176 "register_tracer must be called before tracer_hordiff.")
177 if (loc(reg)==0)
call mom_error(fatal,
"MOM_tracer_hor_diff: "// &
178 "register_tracer must be called before tracer_hordiff.")
179 if ((reg%ntr==0) .or. ((cs%KhTr <= 0.0) .and. .not.
associated(varmix)) )
return 181 if (cs%show_call_tree)
call calltree_enter(
"tracer_hordiff(), MOM_tracer_hor_diff.F90")
183 call cpu_clock_begin(id_clock_diffuse)
187 h_neglect = gv%H_subroundoff
189 if (cs%Diffuse_ML_interior .and. cs%first_call)
then ;
if (is_root_pe())
then 190 do m=1,ntr ;
if (
associated(reg%Tr(m)%df_x) .or.
associated(reg%Tr(m)%df_y)) &
191 call mom_error(warning,
"tracer_hordiff: Tracer "//trim(reg%Tr(m)%name)// &
192 " has associated 3-d diffusive flux diagnostics. These are not "//&
193 "valid when DIFFUSE_ML_TO_INTERIOR is defined. Use 2-d tracer "//&
194 "diffusion diagnostics instead to get accurate total fluxes.")
197 cs%first_call = .false.
199 if (cs%debug)
call mom_tracer_chksum(
"Before tracer diffusion ", reg%Tr, ntr, g)
201 use_varmix = .false. ; resoln_scaled = .false. ; use_eady = .false.
202 if (
Associated(varmix))
then 203 use_varmix = varmix%use_variable_mixing
204 resoln_scaled = varmix%Resoln_scaled_KhTr
205 use_eady = cs%KhTr_Slope_Cff > 0.
208 call cpu_clock_begin(id_clock_pass)
212 call cpu_clock_end(id_clock_pass)
214 if (cs%show_call_tree)
call calltree_waypoint(
"Calculating diffusivity (tracer_hordiff)")
219 do j=js,je ;
do i=is-1,ie
221 if (use_eady) kh_loc = kh_loc + cs%KhTr_Slope_Cff*varmix%L2u(i,j)*varmix%SN_u(i,j)
222 if (
associated(meke%Kh)) &
223 kh_loc = kh_loc + meke%KhTr_fac*sqrt(meke%Kh(i,j)*meke%Kh(i+1,j))
224 if (cs%KhTr_max > 0.) kh_loc = min(kh_loc, cs%KhTr_max)
226 kh_loc = kh_loc * 0.5*(varmix%Res_fn_h(i,j) + varmix%Res_fn_h(i+1,j))
227 kh_u(i,j) = max(kh_loc, cs%KhTr_min)
228 if (cs%KhTr_passivity_coeff>0.)
then 229 rd_dx=0.5*( varmix%Rd_dx_h(i,j)+varmix%Rd_dx_h(i+1,j) )
230 kh_loc = kh_u(i,j)*max( cs%KhTr_passivity_min, cs%KhTr_passivity_coeff*rd_dx )
231 if (cs%KhTr_max > 0.) kh_loc = min(kh_loc, cs%KhTr_max)
232 kh_u(i,j) = max(kh_loc, cs%KhTr_min)
236 do j=js-1,je ;
do i=is,ie
238 if (use_eady) kh_loc = kh_loc + cs%KhTr_Slope_Cff*varmix%L2v(i,j)*varmix%SN_v(i,j)
239 if (
associated(meke%Kh)) &
240 kh_loc = kh_loc + meke%KhTr_fac*sqrt(meke%Kh(i,j)*meke%Kh(i,j+1))
241 if (cs%KhTr_max > 0.) kh_loc = min(kh_loc, cs%KhTr_max)
243 kh_loc = kh_loc * 0.5*(varmix%Res_fn_h(i,j) + varmix%Res_fn_h(i,j+1))
244 kh_v(i,j) = max(kh_loc, cs%KhTr_min)
245 if (cs%KhTr_passivity_coeff>0.)
then 246 rd_dx = 0.5*( varmix%Rd_dx_h(i,j)+varmix%Rd_dx_h(i,j+1) )
247 kh_loc = kh_v(i,j)*max( cs%KhTr_passivity_min, cs%KhTr_passivity_coeff*rd_dx )
248 if (cs%KhTr_max > 0.) kh_loc = min(kh_loc, cs%KhTr_max)
249 kh_v(i,j) = max(kh_loc, cs%KhTr_min)
254 do j=js,je ;
do i=is-1,ie
255 khdt_x(i,j) = dt*(kh_u(i,j)*(g%dy_Cu(i,j)*g%IdxCu(i,j)))
258 do j=js-1,je ;
do i=is,ie
259 khdt_y(i,j) = dt*(kh_v(i,j)*(g%dx_Cv(i,j)*g%IdyCv(i,j)))
261 elseif (resoln_scaled)
then 263 do j=js,je ;
do i=is-1,ie
264 res_fn = 0.5 * (varmix%Res_fn_h(i,j) + varmix%Res_fn_h(i+1,j))
265 kh_u(i,j) = max(cs%KhTr * res_fn, cs%KhTr_min)
266 khdt_x(i,j) = dt*(cs%KhTr*(g%dy_Cu(i,j)*g%IdxCu(i,j))) * res_fn
269 do j=js-1,je ;
do i=is,ie
270 res_fn = 0.5*(varmix%Res_fn_h(i,j) + varmix%Res_fn_h(i,j+1))
271 kh_v(i,j) = max(cs%KhTr * res_fn, cs%KhTr_min)
272 khdt_y(i,j) = dt*(cs%KhTr*(g%dx_Cv(i,j)*g%IdyCv(i,j))) * res_fn
275 if (cs%id_KhTr_u > 0)
then 277 do j=js,je ;
do i=is-1,ie
279 khdt_x(i,j) = dt*(cs%KhTr*(g%dy_Cu(i,j)*g%IdxCu(i,j)))
283 do j=js,je ;
do i=is-1,ie
284 khdt_x(i,j) = dt*(cs%KhTr*(g%dy_Cu(i,j)*g%IdxCu(i,j)))
287 if (cs%id_KhTr_v > 0)
then 289 do j=js-1,je ;
do i=is,ie
291 khdt_y(i,j) = dt*(cs%KhTr*(g%dx_Cv(i,j)*g%IdyCv(i,j)))
295 do j=js-1,je ;
do i=is,ie
296 khdt_y(i,j) = dt*(cs%KhTr*(g%dx_Cv(i,j)*g%IdyCv(i,j)))
301 if (cs%max_diff_CFL > 0.0)
then 302 if ((cs%id_KhTr_u > 0) .or. (cs%id_KhTr_h > 0))
then 304 do j=js,je ;
do i=is-1,ie
305 khdt_max = 0.125*cs%max_diff_CFL * min(g%areaT(i,j), g%areaT(i+1,j))
306 if (khdt_x(i,j) > khdt_max)
then 307 khdt_x(i,j) = khdt_max
308 if (dt*(g%dy_Cu(i,j)*g%IdxCu(i,j)) > 0.0) &
309 kh_u(i,j) = khdt_x(i,j) / (dt*(g%dy_Cu(i,j)*g%IdxCu(i,j)))
314 do j=js,je ;
do i=is-1,ie
315 khdt_max = 0.125*cs%max_diff_CFL * min(g%areaT(i,j), g%areaT(i+1,j))
316 khdt_x(i,j) = min(khdt_x(i,j), khdt_max)
319 if ((cs%id_KhTr_v > 0) .or. (cs%id_KhTr_h > 0))
then 321 do j=js-1,je ;
do i=is,ie
322 khdt_max = 0.125*cs%max_diff_CFL * min(g%areaT(i,j), g%areaT(i,j+1))
323 if (khdt_y(i,j) > khdt_max)
then 324 khdt_y(i,j) = khdt_max
325 if (dt*(g%dx_Cv(i,j)*g%IdyCv(i,j)) > 0.0) &
326 kh_v(i,j) = khdt_y(i,j) / (dt*(g%dx_Cv(i,j)*g%IdyCv(i,j)))
331 do j=js-1,je ;
do i=is,ie
332 khdt_max = 0.125*cs%max_diff_CFL * min(g%areaT(i,j), g%areaT(i,j+1))
333 khdt_y(i,j) = min(khdt_y(i,j), khdt_max)
340 do j=js,je ;
do i=is-1,ie
341 khdt_x(i,j) = us%m_to_L**2*read_khdt_x(i,j)
344 do j=js-1,je ;
do i=is,ie
345 khdt_y(i,j) = us%m_to_L**2*read_khdt_y(i,j)
350 if (cs%check_diffusive_CFL)
then 351 if (cs%show_call_tree)
call calltree_waypoint(
"Checking diffusive CFL (tracer_hordiff)")
353 do j=js,je ;
do i=is,ie
354 cfl(i,j) = 2.0*((khdt_x(i-1,j) + khdt_x(i,j)) + &
355 (khdt_y(i,j-1) + khdt_y(i,j))) * g%IareaT(i,j)
356 if (max_cfl < cfl(i,j)) max_cfl = cfl(i,j)
358 call cpu_clock_begin(id_clock_sync)
359 call max_across_pes(max_cfl)
360 call cpu_clock_end(id_clock_sync)
361 num_itts = max(1, ceiling(max_cfl - 4.0*epsilon(max_cfl)))
362 i_numitts = 1.0 / (
real(num_itts))
363 if (cs%id_CFL > 0)
call post_data(cs%id_CFL, cfl, cs%diag, mask=g%mask2dT)
364 elseif (cs%max_diff_CFL > 0.0)
then 365 num_itts = max(1, ceiling(cs%max_diff_CFL - 4.0*epsilon(cs%max_diff_CFL)))
366 i_numitts = 1.0 / (
real(num_itts))
368 num_itts = 1 ; i_numitts = 1.0
372 if (
associated(reg%Tr(m)%df_x))
then 373 do k=1,nz ;
do j=js,je ;
do i=is-1,ie
374 reg%Tr(m)%df_x(i,j,k) = 0.0
375 enddo ;
enddo ;
enddo 377 if (
associated(reg%Tr(m)%df_y))
then 378 do k=1,nz ;
do j=js-1,je ;
do i=is,ie
379 reg%Tr(m)%df_y(i,j,k) = 0.0
380 enddo ;
enddo ;
enddo 382 if (
associated(reg%Tr(m)%df2d_x))
then 383 do j=js,je ;
do i=is-1,ie ; reg%Tr(m)%df2d_x(i,j) = 0.0 ;
enddo ;
enddo 385 if (
associated(reg%Tr(m)%df2d_y))
then 386 do j=js-1,je ;
do i=is,ie ; reg%Tr(m)%df2d_y(i,j) = 0.0 ;
enddo ;
enddo 390 if (cs%use_lateral_boundary_diffusion)
then 392 if (cs%show_call_tree)
call calltree_waypoint(
"Calling lateral boundary mixing (tracer_hordiff)")
394 call do_group_pass(cs%pass_t, g%Domain, clock=id_clock_pass)
396 do j=js-1,je ;
do i=is,ie
397 coef_y(i,j) = i_numitts * khdt_y(i,j)
401 coef_x(i,j) = i_numitts * khdt_x(i,j)
406 if (cs%show_call_tree)
call calltree_waypoint(
"Calling lateral boundary diffusion (tracer_hordiff)",itt)
408 call do_group_pass(cs%pass_t, g%Domain, clock=id_clock_pass)
410 call lateral_boundary_diffusion(g, gv, us, h, coef_x, coef_y, i_numitts*dt, reg, &
411 cs%lateral_boundary_diffusion_CSp)
415 if (cs%use_neutral_diffusion)
then 417 if (cs%show_call_tree)
call calltree_waypoint(
"Calling neutral diffusion coeffs (tracer_hordiff)")
419 call do_group_pass(cs%pass_t, g%Domain, clock=id_clock_pass)
424 if (
associated(tv%p_surf))
then 425 call neutral_diffusion_calc_coeffs(g, gv, us, h, tv%T, tv%S, cs%neutral_diffusion_CSp, p_surf=tv%p_surf)
427 call neutral_diffusion_calc_coeffs(g, gv, us, h, tv%T, tv%S, cs%neutral_diffusion_CSp)
429 do j=js-1,je ;
do i=is,ie
430 coef_y(i,j) = i_numitts * khdt_y(i,j)
434 coef_x(i,j) = i_numitts * khdt_x(i,j)
439 if (cs%show_call_tree)
call calltree_waypoint(
"Calling neutral diffusion (tracer_hordiff)",itt)
441 call do_group_pass(cs%pass_t, g%Domain, clock=id_clock_pass)
442 if (cs%recalc_neutral_surf)
then 443 if (
associated(tv%p_surf))
then 444 call neutral_diffusion_calc_coeffs(g, gv, us, h, tv%T, tv%S, cs%neutral_diffusion_CSp, p_surf=tv%p_surf)
446 call neutral_diffusion_calc_coeffs(g, gv, us, h, tv%T, tv%S, cs%neutral_diffusion_CSp)
450 call neutral_diffusion(g, gv, h, coef_x, coef_y, i_numitts*dt, reg, us, cs%neutral_diffusion_CSp)
455 if (cs%show_call_tree)
call calltree_waypoint(
"Calculating horizontal diffusion (tracer_hordiff)")
457 call do_group_pass(cs%pass_t, g%Domain, clock=id_clock_pass)
461 if (cs%Diffuse_ML_interior)
then 463 if (cs%ML_KhTr_scale <= 0.0) cycle
464 scale = i_numitts * cs%ML_KhTr_scale
466 if ((k>gv%nkml) .and. (k<=gv%nk_rho_varies)) cycle
469 do j=js-1,je ;
do i=is,ie
470 coef_y(i,j) = ((scale * khdt_y(i,j))*2.0*(h(i,j,k)*h(i,j+1,k))) / &
471 (h(i,j,k)+h(i,j+1,k)+h_neglect)
476 coef_x(i,j) = ((scale * khdt_x(i,j))*2.0*(h(i,j,k)*h(i+1,j,k))) / &
477 (h(i,j,k)+h(i+1,j,k)+h_neglect)
481 ihdxdy(i,j) = g%IareaT(i,j) / (h(i,j,k)+h_neglect)
486 do j=js,je ;
do i=is,ie
487 dtr(i,j) = ihdxdy(i,j) * &
488 ((coef_x(i-1,j) * (reg%Tr(m)%t(i-1,j,k) - reg%Tr(m)%t(i,j,k)) - &
489 coef_x(i,j) * (reg%Tr(m)%t(i,j,k) - reg%Tr(m)%t(i+1,j,k))) + &
490 (coef_y(i,j-1) * (reg%Tr(m)%t(i,j-1,k) - reg%Tr(m)%t(i,j,k)) - &
491 coef_y(i,j) * (reg%Tr(m)%t(i,j,k) - reg%Tr(m)%t(i,j+1,k))))
493 if (
associated(reg%Tr(m)%df_x))
then ;
do j=js,je ;
do i=g%IscB,g%IecB
494 reg%Tr(m)%df_x(i,j,k) = reg%Tr(m)%df_x(i,j,k) + coef_x(i,j) &
495 * (reg%Tr(m)%t(i,j,k) - reg%Tr(m)%t(i+1,j,k)) * idt
496 enddo ;
enddo ;
endif 497 if (
associated(reg%Tr(m)%df_y))
then ;
do j=g%JscB,g%JecB ;
do i=is,ie
498 reg%Tr(m)%df_y(i,j,k) = reg%Tr(m)%df_y(i,j,k) + coef_y(i,j) &
499 * (reg%Tr(m)%t(i,j,k) - reg%Tr(m)%t(i,j+1,k)) * idt
500 enddo ;
enddo ;
endif 501 if (
associated(reg%Tr(m)%df2d_x))
then ;
do j=js,je ;
do i=g%IscB,g%IecB
502 reg%Tr(m)%df2d_x(i,j) = reg%Tr(m)%df2d_x(i,j) + coef_x(i,j) &
503 * (reg%Tr(m)%t(i,j,k) - reg%Tr(m)%t(i+1,j,k)) * idt
504 enddo ;
enddo ;
endif 505 if (
associated(reg%Tr(m)%df2d_y))
then ;
do j=g%JscB,g%JecB ;
do i=is,ie
506 reg%Tr(m)%df2d_y(i,j) = reg%Tr(m)%df2d_y(i,j) + coef_y(i,j) &
507 * (reg%Tr(m)%t(i,j,k) - reg%Tr(m)%t(i,j+1,k)) * idt
508 enddo ;
enddo ;
endif 509 do j=js,je ;
do i=is,ie
510 reg%Tr(m)%t(i,j,k) = reg%Tr(m)%t(i,j,k) + dtr(i,j)
519 call cpu_clock_end(id_clock_diffuse)
522 if (cs%Diffuse_ML_interior)
then 523 if (cs%show_call_tree)
call calltree_waypoint(
"Calling epipycnal_ML_diff (tracer_hordiff)")
524 if (cs%debug)
call mom_tracer_chksum(
"Before epipycnal diff ", reg%Tr, ntr, g)
526 call cpu_clock_begin(id_clock_epimix)
527 call tracer_epipycnal_ml_diff(h, dt, reg%Tr, ntr, khdt_x, khdt_y, g, gv, us, &
529 call cpu_clock_end(id_clock_epimix)
532 if (cs%debug)
call mom_tracer_chksum(
"After tracer diffusion ", reg%Tr, ntr, g)
535 if (cs%id_KhTr_u > 0)
then 536 do j=js,je ;
do i=is-1,ie
537 kh_u(i,j) = g%mask2dCu(i,j)*kh_u(i,j)
539 call post_data(cs%id_KhTr_u, kh_u, cs%diag, mask=g%mask2dCu)
541 if (cs%id_KhTr_v > 0)
then 542 do j=js-1,je ;
do i=is,ie
543 kh_v(i,j) = g%mask2dCv(i,j)*kh_v(i,j)
545 call post_data(cs%id_KhTr_v, kh_v, cs%diag, mask=g%mask2dCv)
547 if (cs%id_KhTr_h > 0)
then 549 do j=js,je ;
do i=is-1,ie
550 kh_u(i,j) = g%mask2dCu(i,j)*kh_u(i,j)
552 do j=js-1,je ;
do i=is,ie
553 kh_v(i,j) = g%mask2dCv(i,j)*kh_v(i,j)
555 do j=js,je ;
do i=is,ie
556 normalize = 1.0 / ((g%mask2dCu(i-1,j)+g%mask2dCu(i,j)) + &
557 (g%mask2dCv(i,j-1)+g%mask2dCv(i,j)) + gv%H_subroundoff)
558 kh_h(i,j) = normalize*g%mask2dT(i,j)*((kh_u(i-1,j)+kh_u(i,j)) + &
559 (kh_v(i,j-1)+kh_v(i,j)))
561 call post_data(cs%id_KhTr_h, kh_h, cs%diag, mask=g%mask2dT)
566 call uvchksum(
"After tracer diffusion khdt_[xy]", khdt_x, khdt_y, &
567 g%HI, haloshift=0, symmetric=.true., scale=us%L_to_m**2, &
569 if (cs%use_neutral_diffusion)
then 570 call uvchksum(
"After tracer diffusion Coef_[xy]", coef_x, coef_y, &
571 g%HI, haloshift=0, symmetric=.true., scale=us%L_to_m**2, &
576 if (cs%id_khdt_x > 0)
call post_data(cs%id_khdt_x, khdt_x, cs%diag)
577 if (cs%id_khdt_y > 0)
call post_data(cs%id_khdt_y, khdt_y, cs%diag)
579 if (cs%show_call_tree)
call calltree_leave(
"tracer_hordiff()")
581 end subroutine tracer_hordiff
587 subroutine tracer_epipycnal_ml_diff(h, dt, Tr, ntr, khdt_epi_x, khdt_epi_y, G, &
588 GV, US, CS, tv, num_itts)
589 type(ocean_grid_type),
intent(inout) :: G
590 type(verticalGrid_type),
intent(in) :: GV
591 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
592 real,
intent(in) :: dt
593 type(tracer_type),
intent(inout) :: Tr(:)
594 integer,
intent(in) :: ntr
595 real,
dimension(SZIB_(G),SZJ_(G)),
intent(in) :: khdt_epi_x
598 real,
dimension(SZI_(G),SZJB_(G)),
intent(in) :: khdt_epi_y
601 type(unit_scale_type),
intent(in) :: US
602 type(tracer_hor_diff_CS),
intent(inout) :: CS
603 type(thermo_var_ptrs),
intent(in) :: tv
604 integer,
intent(in) :: num_itts
607 real,
dimension(SZI_(G), SZJ_(G)) :: &
609 real,
dimension(SZI_(G), SZJ_(G), max(1,GV%nk_rho_varies)) :: &
614 type(p2d),
dimension(SZJ_(G)) :: &
615 deep_wt_Lu, deep_wt_Ru, &
618 type(p2d),
dimension(SZJB_(G)) :: &
619 deep_wt_Lv, deep_wt_Rv, &
622 type(p2di),
dimension(SZJ_(G)) :: &
625 type(p2di),
dimension(SZJB_(G)) :: &
629 real,
dimension(SZI_(G), SZJ_(G), SZK_(G)) :: &
631 real,
dimension(SZI_(G), SZJ_(G), SZK_(G)) :: Tr_flux_3d, Tr_adj_vert_L, Tr_adj_vert_R
633 real,
dimension(SZI_(G), SZK_(G), SZJ_(G)) :: &
636 integer,
dimension(SZI_(G), SZK_(G), SZJ_(G)) :: &
640 real,
dimension(SZK_(G)) :: &
647 integer,
dimension(SZI_(G), SZJ_(G)) :: &
653 integer,
dimension(SZJ_(G)) :: &
655 integer,
dimension(SZIB_(G), SZJ_(G)) :: &
657 integer,
dimension(SZI_(G), SZJB_(G)) :: &
663 real :: rho_pair, rho_a, rho_b
679 integer,
dimension(SZK_(G)*2) :: &
682 logical,
dimension(SZK_(G)*2) :: &
687 real :: p_ref_cv(SZI_(G))
689 integer,
dimension(2) :: EOSdom
690 integer :: k_max, k_min, k_test, itmp
691 integer :: i, j, k, k2, m, is, ie, js, je, nz, nkmb
692 integer :: isd, ied, jsd, jed, IsdB, IedB, k_size
693 integer :: kL, kR, kLa, kLb, kRa, kRb, nP, itt, ns, max_itt
694 integer :: PEmax_kRho
695 integer :: isv, iev, jsv, jev
697 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
698 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
699 isdb = g%IsdB ; iedb = g%IedB
701 nkmb = gv%nk_rho_varies
703 if (num_itts <= 1)
then 704 max_itt = 1 ; i_maxitt = 1.0
706 max_itt = num_itts ; i_maxitt = 1.0 / (
real(max_itt))
709 do i=is-2,ie+2 ; p_ref_cv(i) = tv%P_Ref ;
enddo 710 eosdom(:) = eos_domain(g%HI,halo=2)
712 call do_group_pass(cs%pass_t, g%Domain, clock=id_clock_pass)
715 do k=1,nkmb ;
do j=js-2,je+2
716 call calculate_density(tv%T(:,j,k),tv%S(:,j,k), p_ref_cv, rho_coord(:,j,k), &
717 tv%eqn_of_state, eosdom)
720 do j=js-2,je+2 ;
do i=is-2,ie+2
721 rml_max(i,j) = rho_coord(i,j,1)
722 num_srt(i,j) = 0 ; max_krho(i,j) = 0
724 do k=2,nkmb ;
do j=js-2,je+2 ;
do i=is-2,ie+2
725 if (rml_max(i,j) < rho_coord(i,j,k)) rml_max(i,j) = rho_coord(i,j,k)
726 enddo ;
enddo ;
enddo 731 do j=js-2,je+2 ;
do i=is-2,ie+2 ;
if (g%mask2dT(i,j) > 0.5)
then 732 if ((rml_max(i,j) > gv%Rlay(nz)) .or. (nkmb+1 > nz))
then ; max_krho(i,j) = nz+1
733 elseif ((rml_max(i,j) <= gv%Rlay(nkmb+1)) .or. (nkmb+2 > nz))
then ; max_krho(i,j) = nkmb+1
735 k_min = nkmb+2 ; k_max = nz
737 k_test = (k_min + k_max) / 2
738 if (rml_max(i,j) <= gv%Rlay(k_test-1))
then ; k_max = k_test-1
739 elseif (gv%Rlay(k_test) < rml_max(i,j))
then ; k_min = k_test+1
740 else ; max_krho(i,j) = k_test ;
exit ;
endif 742 if (k_min == k_max)
then ; max_krho(i,j) = k_max ;
exit ;
endif 745 endif ;
enddo ;
enddo 748 do j=js-1,je+1 ;
do i=is-1,ie+1
749 k_end_srt(i,j) = max(max_krho(i,j), max_krho(i-1,j), max_krho(i+1,j), &
750 max_krho(i,j-1), max_krho(i,j+1))
751 if (pemax_krho < k_end_srt(i,j)) pemax_krho = k_end_srt(i,j)
753 if (pemax_krho > nz) pemax_krho = nz
755 h_exclude = 10.0*(gv%Angstrom_H + gv%H_subroundoff)
759 do k=1,nkmb ;
do i=is-1,ie+1 ;
if (g%mask2dT(i,j) > 0.5)
then 760 if (h(i,j,k) > h_exclude)
then 761 num_srt(i,j) = num_srt(i,j) + 1 ; ns = num_srt(i,j)
763 rho_srt(i,ns,j) = rho_coord(i,j,k)
764 h_srt(i,ns,j) = h(i,j,k)
766 endif ;
enddo ;
enddo 767 do k=nkmb+1,pemax_krho ;
do i=is-1,ie+1 ;
if (g%mask2dT(i,j) > 0.5)
then 768 if ((k<=k_end_srt(i,j)) .and. (h(i,j,k) > h_exclude))
then 769 num_srt(i,j) = num_srt(i,j) + 1 ; ns = num_srt(i,j)
771 rho_srt(i,ns,j) = gv%Rlay(k)
772 h_srt(i,ns,j) = h(i,j,k)
774 endif ;
enddo ;
enddo 779 do j=js-1,je+1;
do i=is-1,ie+1
780 do k=2,num_srt(i,j) ;
if (rho_srt(i,k,j) < rho_srt(i,k-1,j))
then 782 do k2 = k,2,-1 ;
if (rho_srt(i,k2,j) >= rho_srt(i,k2-1,j))
exit 783 itmp = k0_srt(i,k2-1,j) ; k0_srt(i,k2-1,j) = k0_srt(i,k2,j) ; k0_srt(i,k2,j) = itmp
784 tmp = rho_srt(i,k2-1,j) ; rho_srt(i,k2-1,j) = rho_srt(i,k2,j) ; rho_srt(i,k2,j) = tmp
785 tmp = h_srt(i,k2-1,j) ; h_srt(i,k2-1,j) = h_srt(i,k2,j) ; h_srt(i,k2,j) = tmp
792 do i=is-1,ie+1 ; max_srt(j) = max(max_srt(j), num_srt(i,j)) ;
enddo 797 k_size = max(2*max_srt(j),1)
798 allocate(deep_wt_lu(j)%p(isdb:iedb,k_size))
799 allocate(deep_wt_ru(j)%p(isdb:iedb,k_size))
800 allocate(hp_lu(j)%p(isdb:iedb,k_size))
801 allocate(hp_ru(j)%p(isdb:iedb,k_size))
802 allocate(k0a_lu(j)%p(isdb:iedb,k_size))
803 allocate(k0a_ru(j)%p(isdb:iedb,k_size))
804 allocate(k0b_lu(j)%p(isdb:iedb,k_size))
805 allocate(k0b_ru(j)%p(isdb:iedb,k_size))
815 do j=js,je ;
do i=is-1,ie ;
if (g%mask2dCu(i,j) > 0.5)
then 818 do k=1,num_srt(i,j) ; h_demand_l(k) = 0.0 ; h_used_l(k) = 0.0 ;
enddo 819 do k=1,num_srt(i+1,j) ; h_demand_r(k) = 0.0 ; h_used_r(k) = 0.0 ;
enddo 826 if (rho_srt(i,1,j) < rho_srt(i+1,1,j))
then 828 do kl=2,num_srt(i,j) ;
if (rho_srt(i,kl,j) >= rho_srt(i+1,1,j))
exit ;
enddo 829 elseif (rho_srt(i+1,1,j) < rho_srt(i,1,j))
then 831 do kr=2,num_srt(i+1,j) ;
if (rho_srt(i+1,kr,j) >= rho_srt(i,1,j))
exit ;
enddo 837 if ((kl > num_srt(i,j)) .or. (kr > num_srt(i+1,j)))
exit 839 if (rho_srt(i,kl,j) > rho_srt(i+1,kr,j))
then 842 rho_pair = rho_srt(i+1,kr,j)
844 k0b_lu(j)%p(i,k) = k0_srt(i,kl,j) ; k0b_ru(j)%p(i,k) = k0_srt(i+1,kr,j)
845 k0a_lu(j)%p(i,k) = k0_srt(i,kl-1,j) ; k0a_ru(j)%p(i,k) = k0b_ru(j)%p(i,k)
846 kbs_lp(k) = kl ; kbs_rp(k) = kr
848 rho_a = rho_srt(i,kl-1,j) ; rho_b = rho_srt(i,kl,j)
849 wt_b = 1.0 ;
if (abs(rho_a - rho_b) > abs(rho_pair - rho_a)) &
850 wt_b = (rho_pair - rho_a) / (rho_b - rho_a)
851 deep_wt_lu(j)%p(i,k) = wt_b ; deep_wt_ru(j)%p(i,k) = 1.0
853 h_demand_l(kl) = h_demand_l(kl) + 0.5*h_srt(i+1,kr,j) * wt_b
854 h_demand_l(kl-1) = h_demand_l(kl-1) + 0.5*h_srt(i+1,kr,j) * (1.0-wt_b)
856 kr = kr+1 ; left_set(k) = .false. ; right_set(k) = .true.
857 elseif (rho_srt(i,kl,j) < rho_srt(i+1,kr,j))
then 860 rho_pair = rho_srt(i,kl,j)
861 k0b_lu(j)%p(i,k) = k0_srt(i,kl,j) ; k0b_ru(j)%p(i,k) = k0_srt(i+1,kr,j)
862 k0a_lu(j)%p(i,k) = k0b_lu(j)%p(i,k) ; k0a_ru(j)%p(i,k) = k0_srt(i+1,kr-1,j)
864 kbs_lp(k) = kl ; kbs_rp(k) = kr
866 rho_a = rho_srt(i+1,kr-1,j) ; rho_b = rho_srt(i+1,kr,j)
867 wt_b = 1.0 ;
if (abs(rho_a - rho_b) > abs(rho_pair - rho_a)) &
868 wt_b = (rho_pair - rho_a) / (rho_b - rho_a)
869 deep_wt_lu(j)%p(i,k) = 1.0 ; deep_wt_ru(j)%p(i,k) = wt_b
871 h_demand_r(kr) = h_demand_r(kr) + 0.5*h_srt(i,kl,j) * wt_b
872 h_demand_r(kr-1) = h_demand_r(kr-1) + 0.5*h_srt(i,kl,j) * (1.0-wt_b)
874 kl = kl+1 ; left_set(k) = .true. ; right_set(k) = .false.
875 elseif ((k0_srt(i,kl,j) <= nkmb) .or. (k0_srt(i+1,kr,j) <= nkmb))
then 878 k0b_lu(j)%p(i,k) = k0_srt(i,kl,j) ; k0b_ru(j)%p(i,k) = k0_srt(i+1,kr,j)
879 k0a_lu(j)%p(i,k) = k0b_lu(j)%p(i,k) ; k0a_ru(j)%p(i,k) = k0b_ru(j)%p(i,k)
880 kbs_lp(k) = kl ; kbs_rp(k) = kr
881 deep_wt_lu(j)%p(i,k) = 1.0 ; deep_wt_ru(j)%p(i,k) = 1.0
883 h_demand_l(kl) = h_demand_l(kl) + 0.5*h_srt(i+1,kr,j)
884 h_demand_r(kr) = h_demand_r(kr) + 0.5*h_srt(i,kl,j)
886 kl = kl+1 ; kr = kr+1 ; left_set(k) = .true. ; right_set(k) = .true.
890 h_demand_l(kl) = h_demand_l(kl) + 0.5*h_srt(i+1,kr,j)
891 h_demand_r(kr) = h_demand_r(kr) + 0.5*h_srt(i,kl,j)
893 kl = kl+1 ; kr = kr+1
899 do k=1,num_srt(i+1,j)
900 h_supply_frac_r(k) = 1.0
901 if (h_demand_r(k) > 0.5*h_srt(i+1,k,j)) &
902 h_supply_frac_r(k) = 0.5*h_srt(i+1,k,j) / h_demand_r(k)
905 h_supply_frac_l(k) = 1.0
906 if (h_demand_l(k) > 0.5*h_srt(i,k,j)) &
907 h_supply_frac_l(k) = 0.5*h_srt(i,k,j) / h_demand_l(k)
912 kl = kbs_lp(k) ; kr = kbs_rp(k)
913 hp_lu(j)%p(i,k) = 0.0 ; hp_ru(j)%p(i,k) = 0.0
914 if (left_set(k))
then 915 if (deep_wt_ru(j)%p(i,k) < 1.0)
then 916 hp_ru(j)%p(i,k) = 0.5*h_srt(i,kl,j) * min(h_supply_frac_r(kr), h_supply_frac_r(kr-1))
917 wt_b = deep_wt_ru(j)%p(i,k)
918 h_used_r(kr-1) = h_used_r(kr-1) + (1.0 - wt_b)*hp_ru(j)%p(i,k)
919 h_used_r(kr) = h_used_r(kr) + wt_b*hp_ru(j)%p(i,k)
921 hp_ru(j)%p(i,k) = 0.5*h_srt(i,kl,j) * h_supply_frac_r(kr)
922 h_used_r(kr) = h_used_r(kr) + hp_ru(j)%p(i,k)
925 if (right_set(k))
then 926 if (deep_wt_lu(j)%p(i,k) < 1.0)
then 927 hp_lu(j)%p(i,k) = 0.5*h_srt(i+1,kr,j) * min(h_supply_frac_l(kl), h_supply_frac_l(kl-1))
928 wt_b = deep_wt_lu(j)%p(i,k)
929 h_used_l(kl-1) = h_used_l(kl-1) + (1.0 - wt_b)*hp_lu(j)%p(i,k)
930 h_used_l(kl) = h_used_l(kl) + wt_b*hp_lu(j)%p(i,k)
932 hp_lu(j)%p(i,k) = 0.5*h_srt(i+1,kr,j) * h_supply_frac_l(kl)
933 h_used_l(kl) = h_used_l(kl) + hp_lu(j)%p(i,k)
941 if (left_set(k)) hp_lu(j)%p(i,k) = hp_lu(j)%p(i,k) + &
942 (h_srt(i,kbs_lp(k),j) - h_used_l(kbs_lp(k)))
943 if (right_set(k)) hp_ru(j)%p(i,k) = hp_ru(j)%p(i,k) + &
944 (h_srt(i+1,kbs_rp(k),j) - h_used_r(kbs_rp(k)))
947 endif ;
enddo ;
enddo 950 k_size = max(max_srt(j)+max_srt(j+1),1)
951 allocate(deep_wt_lv(j)%p(isd:ied,k_size))
952 allocate(deep_wt_rv(j)%p(isd:ied,k_size))
953 allocate(hp_lv(j)%p(isd:ied,k_size))
954 allocate(hp_rv(j)%p(isd:ied,k_size))
955 allocate(k0a_lv(j)%p(isd:ied,k_size))
956 allocate(k0a_rv(j)%p(isd:ied,k_size))
957 allocate(k0b_lv(j)%p(isd:ied,k_size))
958 allocate(k0b_rv(j)%p(isd:ied,k_size))
968 do j=js-1,je ;
do i=is,ie ;
if (g%mask2dCv(i,j) > 0.5)
then 971 do k=1,num_srt(i,j) ; h_demand_l(k) = 0.0 ; h_used_l(k) = 0.0 ;
enddo 972 do k=1,num_srt(i,j+1) ; h_demand_r(k) = 0.0 ; h_used_r(k) = 0.0 ;
enddo 979 if (rho_srt(i,1,j) < rho_srt(i,1,j+1))
then 981 do kl=2,num_srt(i,j) ;
if (rho_srt(i,kl,j) >= rho_srt(i,1,j+1))
exit ;
enddo 982 elseif (rho_srt(i,1,j+1) < rho_srt(i,1,j))
then 984 do kr=2,num_srt(i,j+1) ;
if (rho_srt(i,kr,j+1) >= rho_srt(i,1,j))
exit ;
enddo 990 if ((kl > num_srt(i,j)) .or. (kr > num_srt(i,j+1)))
exit 992 if (rho_srt(i,kl,j) > rho_srt(i,kr,j+1))
then 995 rho_pair = rho_srt(i,kr,j+1)
997 k0b_lv(j)%p(i,k) = k0_srt(i,kl,j) ; k0b_rv(j)%p(i,k) = k0_srt(i,kr,j+1)
998 k0a_lv(j)%p(i,k) = k0_srt(i,kl-1,j) ; k0a_rv(j)%p(i,k) = k0b_rv(j)%p(i,k)
999 kbs_lp(k) = kl ; kbs_rp(k) = kr
1001 rho_a = rho_srt(i,kl-1,j) ; rho_b = rho_srt(i,kl,j)
1002 wt_b = 1.0 ;
if (abs(rho_a - rho_b) > abs(rho_pair - rho_a)) &
1003 wt_b = (rho_pair - rho_a) / (rho_b - rho_a)
1004 deep_wt_lv(j)%p(i,k) = wt_b ; deep_wt_rv(j)%p(i,k) = 1.0
1006 h_demand_l(kl) = h_demand_l(kl) + 0.5*h_srt(i,kr,j+1) * wt_b
1007 h_demand_l(kl-1) = h_demand_l(kl-1) + 0.5*h_srt(i,kr,j+1) * (1.0-wt_b)
1009 kr = kr+1 ; left_set(k) = .false. ; right_set(k) = .true.
1010 elseif (rho_srt(i,kl,j) < rho_srt(i,kr,j+1))
then 1013 rho_pair = rho_srt(i,kl,j)
1014 k0b_lv(j)%p(i,k) = k0_srt(i,kl,j) ; k0b_rv(j)%p(i,k) = k0_srt(i,kr,j+1)
1015 k0a_lv(j)%p(i,k) = k0b_lv(j)%p(i,k) ; k0a_rv(j)%p(i,k) = k0_srt(i,kr-1,j+1)
1017 kbs_lp(k) = kl ; kbs_rp(k) = kr
1019 rho_a = rho_srt(i,kr-1,j+1) ; rho_b = rho_srt(i,kr,j+1)
1020 wt_b = 1.0 ;
if (abs(rho_a - rho_b) > abs(rho_pair - rho_a)) &
1021 wt_b = (rho_pair - rho_a) / (rho_b - rho_a)
1022 deep_wt_lv(j)%p(i,k) = 1.0 ; deep_wt_rv(j)%p(i,k) = wt_b
1024 h_demand_r(kr) = h_demand_r(kr) + 0.5*h_srt(i,kl,j) * wt_b
1025 h_demand_r(kr-1) = h_demand_r(kr-1) + 0.5*h_srt(i,kl,j) * (1.0-wt_b)
1027 kl = kl+1 ; left_set(k) = .true. ; right_set(k) = .false.
1028 elseif ((k0_srt(i,kl,j) <= nkmb) .or. (k0_srt(i,kr,j+1) <= nkmb))
then 1031 k0b_lv(j)%p(i,k) = k0_srt(i,kl,j) ; k0b_rv(j)%p(i,k) = k0_srt(i,kr,j+1)
1032 k0a_lv(j)%p(i,k) = k0b_lv(j)%p(i,k) ; k0a_rv(j)%p(i,k) = k0b_rv(j)%p(i,k)
1033 kbs_lp(k) = kl ; kbs_rp(k) = kr
1034 deep_wt_lv(j)%p(i,k) = 1.0 ; deep_wt_rv(j)%p(i,k) = 1.0
1036 h_demand_l(kl) = h_demand_l(kl) + 0.5*h_srt(i,kr,j+1)
1037 h_demand_r(kr) = h_demand_r(kr) + 0.5*h_srt(i,kl,j)
1039 kl = kl+1 ; kr = kr+1 ; left_set(k) = .true. ; right_set(k) = .true.
1043 h_demand_l(kl) = h_demand_l(kl) + 0.5*h_srt(i,kr,j+1)
1044 h_demand_r(kr) = h_demand_r(kr) + 0.5*h_srt(i,kl,j)
1046 kl = kl+1 ; kr = kr+1
1052 do k=1,num_srt(i,j+1)
1053 h_supply_frac_r(k) = 1.0
1054 if (h_demand_r(k) > 0.5*h_srt(i,k,j+1)) &
1055 h_supply_frac_r(k) = 0.5*h_srt(i,k,j+1) / h_demand_r(k)
1058 h_supply_frac_l(k) = 1.0
1059 if (h_demand_l(k) > 0.5*h_srt(i,k,j)) &
1060 h_supply_frac_l(k) = 0.5*h_srt(i,k,j) / h_demand_l(k)
1065 kl = kbs_lp(k) ; kr = kbs_rp(k)
1066 hp_lv(j)%p(i,k) = 0.0 ; hp_rv(j)%p(i,k) = 0.0
1067 if (left_set(k))
then 1068 if (deep_wt_rv(j)%p(i,k) < 1.0)
then 1069 hp_rv(j)%p(i,k) = 0.5*h_srt(i,kl,j) * min(h_supply_frac_r(kr), h_supply_frac_r(kr-1))
1070 wt_b = deep_wt_rv(j)%p(i,k)
1071 h_used_r(kr-1) = h_used_r(kr-1) + (1.0 - wt_b) * hp_rv(j)%p(i,k)
1072 h_used_r(kr) = h_used_r(kr) + wt_b * hp_rv(j)%p(i,k)
1074 hp_rv(j)%p(i,k) = 0.5*h_srt(i,kl,j) * h_supply_frac_r(kr)
1075 h_used_r(kr) = h_used_r(kr) + hp_rv(j)%p(i,k)
1078 if (right_set(k))
then 1079 if (deep_wt_lv(j)%p(i,k) < 1.0)
then 1080 hp_lv(j)%p(i,k) = 0.5*h_srt(i,kr,j+1) * min(h_supply_frac_l(kl), h_supply_frac_l(kl-1))
1081 wt_b = deep_wt_lv(j)%p(i,k)
1082 h_used_l(kl-1) = h_used_l(kl-1) + (1.0 - wt_b) * hp_lv(j)%p(i,k)
1083 h_used_l(kl) = h_used_l(kl) + wt_b * hp_lv(j)%p(i,k)
1085 hp_lv(j)%p(i,k) = 0.5*h_srt(i,kr,j+1) * h_supply_frac_l(kl)
1086 h_used_l(kl) = h_used_l(kl) + hp_lv(j)%p(i,k)
1094 if (left_set(k)) hp_lv(j)%p(i,k) = hp_lv(j)%p(i,k) + &
1095 (h_srt(i,kbs_lp(k),j) - h_used_l(kbs_lp(k)))
1096 if (right_set(k)) hp_rv(j)%p(i,k) = hp_rv(j)%p(i,k) + &
1097 (h_srt(i,kbs_rp(k),j+1) - h_used_r(kbs_rp(k)))
1101 endif ;
enddo ;
enddo 1106 do k=1,pemax_krho ;
do j=js-1,je+1 ;
do i=is-1,ie+1
1107 tr_flux_conv(i,j,k) = 0.0
1108 enddo ;
enddo ;
enddo 1113 call do_group_pass(cs%pass_t, g%Domain, clock=id_clock_pass)
1123 do j=js,je ;
do i=is-1,ie ;
if (g%mask2dCu(i,j) > 0.5)
then 1127 if (npu(i,j) >= 1)
then 1128 tr_min_face = min(tr(m)%t(i,j,1), tr(m)%t(i+1,j,1))
1129 tr_max_face = max(tr(m)%t(i,j,1), tr(m)%t(i+1,j,1))
1131 tr_min_face = min(tr_min_face, tr(m)%t(i,j,k), tr(m)%t(i+1,j,k))
1132 tr_max_face = max(tr_max_face, tr(m)%t(i,j,k), tr(m)%t(i+1,j,k))
1136 kla = nkmb+1 ;
if (max_krho(i,j) < nz+1) kla = max_krho(i,j)
1137 klb = kla ;
if (max_krho(i,j) < nz) klb = max_krho(i,j)+1
1138 kra = nkmb+1 ;
if (max_krho(i+1,j) < nz+1) kra = max_krho(i+1,j)
1139 krb = kra ;
if (max_krho(i+1,j) < nz) krb = max_krho(i+1,j)+1
1140 tr_la = tr_min_face ; tr_lb = tr_la ; tr_ra = tr_la ; tr_rb = tr_la
1141 if (h(i,j,kla) > h_exclude) tr_la = tr(m)%t(i,j,kla)
1142 if (h(i,j,klb) > h_exclude) tr_la = tr(m)%t(i,j,klb)
1143 if (h(i+1,j,kra) > h_exclude) tr_ra = tr(m)%t(i+1,j,kra)
1144 if (h(i+1,j,krb) > h_exclude) tr_rb = tr(m)%t(i+1,j,krb)
1145 tr_min_face = min(tr_min_face, tr_la, tr_lb, tr_ra, tr_rb)
1146 tr_max_face = max(tr_max_face, tr_la, tr_lb, tr_ra, tr_rb)
1150 tr_lb = tr(m)%t(i,j,k0b_lu(j)%p(i,k))
1151 tr_rb = tr(m)%t(i+1,j,k0b_ru(j)%p(i,k))
1152 tr_la = tr_lb ; tr_ra = tr_rb
1153 if (deep_wt_lu(j)%p(i,k) < 1.0) tr_la = tr(m)%t(i,j,k0a_lu(j)%p(i,k))
1154 if (deep_wt_ru(j)%p(i,k) < 1.0) tr_ra = tr(m)%t(i+1,j,k0a_ru(j)%p(i,k))
1155 tr_min_face = min(tr_min_face, tr_la, tr_lb, tr_ra, tr_rb)
1156 tr_max_face = max(tr_max_face, tr_la, tr_lb, tr_ra, tr_rb)
1161 klb = k0b_lu(j)%p(i,k) ; tr_lb = tr(m)%t(i,j,klb) ; tr_av_l = tr_lb
1162 if (deep_wt_lu(j)%p(i,k) < 1.0)
then 1163 kla = k0a_lu(j)%p(i,k) ; tr_la = tr(m)%t(i,j,kla)
1164 wt_b = deep_wt_lu(j)%p(i,k)
1165 tr_av_l = wt_b*tr_lb + (1.0-wt_b)*tr_la
1168 krb = k0b_ru(j)%p(i,k) ; tr_rb = tr(m)%t(i+1,j,krb) ; tr_av_r = tr_rb
1169 if (deep_wt_ru(j)%p(i,k) < 1.0)
then 1170 kra = k0a_ru(j)%p(i,k) ; tr_ra = tr(m)%t(i+1,j,kra)
1171 wt_b = deep_wt_ru(j)%p(i,k)
1172 tr_av_r = wt_b*tr_rb + (1.0-wt_b)*tr_ra
1175 h_l = hp_lu(j)%p(i,k) ; h_r = hp_ru(j)%p(i,k)
1176 tr_flux = i_maxitt * khdt_epi_x(i,j) * (tr_av_l - tr_av_r) * &
1177 ((2.0 * h_l * h_r) / (h_l + h_r))
1180 if (deep_wt_lu(j)%p(i,k) >= 1.0)
then 1181 tr_flux_conv(i,j,klb) = tr_flux_conv(i,j,klb) - tr_flux
1184 wt_b = deep_wt_lu(j)%p(i,k) ; wt_a = 1.0 - wt_b
1185 vol = hp_lu(j)%p(i,k) * g%areaT(i,j)
1193 if (tr_flux > 0.0)
then 1194 if (tr_la < tr_lb)
then ;
if (vol*(tr_la-tr_min_face) < tr_flux) &
1195 tr_adj_vert = -wt_a * min(tr_flux - vol * (tr_la-tr_min_face), &
1196 (vol*wt_b) * (tr_lb - tr_la))
1197 else ;
if (vol*(tr_lb-tr_min_face) < tr_flux) &
1198 tr_adj_vert = wt_b * min(tr_flux - vol * (tr_lb-tr_min_face), &
1199 (vol*wt_a) * (tr_la - tr_lb))
1201 elseif (tr_flux < 0.0)
then 1202 if (tr_la > tr_lb)
then ;
if (vol * (tr_max_face-tr_la) < -tr_flux) &
1203 tr_adj_vert = wt_a * min(-tr_flux - vol * (tr_max_face-tr_la), &
1204 (vol*wt_b) * (tr_la - tr_lb))
1205 else ;
if (vol*(tr_max_face-tr_lb) < -tr_flux) &
1206 tr_adj_vert = -wt_b * min(-tr_flux - vol * (tr_max_face-tr_lb), &
1207 (vol*wt_a)*(tr_lb - tr_la))
1211 tr_flux_conv(i,j,kla) = tr_flux_conv(i,j,kla) - (wt_a*tr_flux + tr_adj_vert)
1212 tr_flux_conv(i,j,klb) = tr_flux_conv(i,j,klb) - (wt_b*tr_flux - tr_adj_vert)
1215 if (deep_wt_ru(j)%p(i,k) >= 1.0)
then 1216 tr_flux_conv(i+1,j,krb) = tr_flux_conv(i+1,j,krb) + tr_flux
1219 wt_b = deep_wt_ru(j)%p(i,k) ; wt_a = 1.0 - wt_b
1220 vol = hp_ru(j)%p(i,k) * g%areaT(i+1,j)
1228 if (tr_flux < 0.0)
then 1229 if (tr_ra < tr_rb)
then ;
if (vol * (tr_ra-tr_min_face) < -tr_flux) &
1230 tr_adj_vert = -wt_a * min(-tr_flux - vol * (tr_ra-tr_min_face), &
1231 (vol*wt_b) * (tr_rb - tr_ra))
1232 else ;
if (vol*(tr_rb-tr_min_face) < (-tr_flux)) &
1233 tr_adj_vert = wt_b * min(-tr_flux - vol * (tr_rb-tr_min_face), &
1234 (vol*wt_a) * (tr_ra - tr_rb))
1236 elseif (tr_flux > 0.0)
then 1237 if (tr_ra > tr_rb)
then ;
if (vol * (tr_max_face-tr_ra) < tr_flux) &
1238 tr_adj_vert = wt_a * min(tr_flux - vol * (tr_max_face-tr_ra), &
1239 (vol*wt_b) * (tr_ra - tr_rb))
1240 else ;
if (vol*(tr_max_face-tr_rb) < tr_flux) &
1241 tr_adj_vert = -wt_b * min(tr_flux - vol * (tr_max_face-tr_rb), &
1242 (vol*wt_a)*(tr_rb - tr_ra))
1246 tr_flux_conv(i+1,j,kra) = tr_flux_conv(i+1,j,kra) + &
1247 (wt_a*tr_flux - tr_adj_vert)
1248 tr_flux_conv(i+1,j,krb) = tr_flux_conv(i+1,j,krb) + &
1249 (wt_b*tr_flux + tr_adj_vert)
1251 if (
associated(tr(m)%df2d_x)) &
1252 tr(m)%df2d_x(i,j) = tr(m)%df2d_x(i,j) + tr_flux * idt
1254 endif ;
enddo ;
enddo 1263 do j=js-1,je ;
do i=is,ie ;
if (g%mask2dCv(i,j) > 0.5)
then 1267 if (npv(i,j) >= 1)
then 1268 tr_min_face = min(tr(m)%t(i,j,1), tr(m)%t(i,j+1,1))
1269 tr_max_face = max(tr(m)%t(i,j,1), tr(m)%t(i,j+1,1))
1271 tr_min_face = min(tr_min_face, tr(m)%t(i,j,k), tr(m)%t(i,j+1,k))
1272 tr_max_face = max(tr_max_face, tr(m)%t(i,j,k), tr(m)%t(i,j+1,k))
1276 kla = nkmb+1 ;
if (max_krho(i,j) < nz+1) kla = max_krho(i,j)
1277 klb = kla ;
if (max_krho(i,j) < nz) klb = max_krho(i,j)+1
1278 kra = nkmb+1 ;
if (max_krho(i,j+1) < nz+1) kra = max_krho(i,j+1)
1279 krb = kra ;
if (max_krho(i,j+1) < nz) krb = max_krho(i,j+1)+1
1280 tr_la = tr_min_face ; tr_lb = tr_la ; tr_ra = tr_la ; tr_rb = tr_la
1281 if (h(i,j,kla) > h_exclude) tr_la = tr(m)%t(i,j,kla)
1282 if (h(i,j,klb) > h_exclude) tr_la = tr(m)%t(i,j,klb)
1283 if (h(i,j+1,kra) > h_exclude) tr_ra = tr(m)%t(i,j+1,kra)
1284 if (h(i,j+1,krb) > h_exclude) tr_rb = tr(m)%t(i,j+1,krb)
1285 tr_min_face = min(tr_min_face, tr_la, tr_lb, tr_ra, tr_rb)
1286 tr_max_face = max(tr_max_face, tr_la, tr_lb, tr_ra, tr_rb)
1290 tr_lb = tr(m)%t(i,j,k0b_lv(j)%p(i,k)) ; tr_rb = tr(m)%t(i,j+1,k0b_rv(j)%p(i,k))
1291 tr_la = tr_lb ; tr_ra = tr_rb
1292 if (deep_wt_lv(j)%p(i,k) < 1.0) tr_la = tr(m)%t(i,j,k0a_lv(j)%p(i,k))
1293 if (deep_wt_rv(j)%p(i,k) < 1.0) tr_ra = tr(m)%t(i,j+1,k0a_rv(j)%p(i,k))
1294 tr_min_face = min(tr_min_face, tr_la, tr_lb, tr_ra, tr_rb)
1295 tr_max_face = max(tr_max_face, tr_la, tr_lb, tr_ra, tr_rb)
1300 klb = k0b_lv(j)%p(i,k) ; tr_lb = tr(m)%t(i,j,klb) ; tr_av_l = tr_lb
1301 if (deep_wt_lv(j)%p(i,k) < 1.0)
then 1302 kla = k0a_lv(j)%p(i,k) ; tr_la = tr(m)%t(i,j,kla)
1303 wt_b = deep_wt_lv(j)%p(i,k)
1304 tr_av_l = wt_b * tr_lb + (1.0-wt_b) * tr_la
1307 krb = k0b_rv(j)%p(i,k) ; tr_rb = tr(m)%t(i,j+1,krb) ; tr_av_r = tr_rb
1308 if (deep_wt_rv(j)%p(i,k) < 1.0)
then 1309 kra = k0a_rv(j)%p(i,k) ; tr_ra = tr(m)%t(i,j+1,kra)
1310 wt_b = deep_wt_rv(j)%p(i,k)
1311 tr_av_r = wt_b * tr_rb + (1.0-wt_b) * tr_ra
1314 h_l = hp_lv(j)%p(i,k) ; h_r = hp_rv(j)%p(i,k)
1315 tr_flux = i_maxitt * ((2.0 * h_l * h_r) / (h_l + h_r)) * &
1316 khdt_epi_y(i,j) * (tr_av_l - tr_av_r)
1317 tr_flux_3d(i,j,k) = tr_flux
1319 if (deep_wt_lv(j)%p(i,k) < 1.0)
then 1321 wt_b = deep_wt_lv(j)%p(i,k) ; wt_a = 1.0 - wt_b
1322 vol = hp_lv(j)%p(i,k) * g%areaT(i,j)
1326 if (tr_flux > 0.0)
then 1327 if (tr_la < tr_lb)
then ;
if (vol * (tr_la-tr_min_face) < tr_flux) &
1328 tr_adj_vert = -wt_a * min(tr_flux - vol * (tr_la-tr_min_face), &
1329 (vol*wt_b) * (tr_lb - tr_la))
1330 else ;
if (vol*(tr_lb-tr_min_face) < tr_flux) &
1331 tr_adj_vert = wt_b * min(tr_flux - vol * (tr_lb-tr_min_face), &
1332 (vol*wt_a) * (tr_la - tr_lb))
1334 elseif (tr_flux < 0.0)
then 1335 if (tr_la > tr_lb)
then ;
if (vol * (tr_max_face-tr_la) < -tr_flux) &
1336 tr_adj_vert = wt_a * min(-tr_flux - vol * (tr_max_face-tr_la), &
1337 (vol*wt_b) * (tr_la - tr_lb))
1338 else ;
if (vol*(tr_max_face-tr_lb) < -tr_flux) &
1339 tr_adj_vert = -wt_b * min(-tr_flux - vol * (tr_max_face-tr_lb), &
1340 (vol*wt_a)*(tr_lb - tr_la))
1343 tr_adj_vert_l(i,j,k) = tr_adj_vert
1346 if (deep_wt_rv(j)%p(i,k) < 1.0)
then 1348 wt_b = deep_wt_rv(j)%p(i,k) ; wt_a = 1.0 - wt_b
1349 vol = hp_rv(j)%p(i,k) * g%areaT(i,j+1)
1353 if (tr_flux < 0.0)
then 1354 if (tr_ra < tr_rb)
then ;
if (vol * (tr_ra-tr_min_face) < -tr_flux) &
1355 tr_adj_vert = -wt_a * min(-tr_flux - vol * (tr_ra-tr_min_face), &
1356 (vol*wt_b) * (tr_rb - tr_ra))
1357 else ;
if (vol*(tr_rb-tr_min_face) < (-tr_flux)) &
1358 tr_adj_vert = wt_b * min(-tr_flux - vol * (tr_rb-tr_min_face), &
1359 (vol*wt_a) * (tr_ra - tr_rb))
1361 elseif (tr_flux > 0.0)
then 1362 if (tr_ra > tr_rb)
then ;
if (vol * (tr_max_face-tr_ra) < tr_flux) &
1363 tr_adj_vert = wt_a * min(tr_flux - vol * (tr_max_face-tr_ra), &
1364 (vol*wt_b) * (tr_ra - tr_rb))
1365 else ;
if (vol*(tr_max_face-tr_rb) < tr_flux) &
1366 tr_adj_vert = -wt_b * min(tr_flux - vol * (tr_max_face-tr_rb), &
1367 (vol*wt_a)*(tr_rb - tr_ra))
1370 tr_adj_vert_r(i,j,k) = tr_adj_vert
1372 if (
associated(tr(m)%df2d_y)) &
1373 tr(m)%df2d_y(i,j) = tr(m)%df2d_y(i,j) + tr_flux * idt
1375 endif ;
enddo ;
enddo 1380 do i=is,ie ;
do j=js-1,je ;
if (g%mask2dCv(i,j) > 0.5)
then 1382 klb = k0b_lv(j)%p(i,k); krb = k0b_rv(j)%p(i,k)
1383 if (deep_wt_lv(j)%p(i,k) >= 1.0)
then 1384 tr_flux_conv(i,j,klb) = tr_flux_conv(i,j,klb) - tr_flux_3d(i,j,k)
1386 kla = k0a_lv(j)%p(i,k)
1387 wt_b = deep_wt_lv(j)%p(i,k) ; wt_a = 1.0 - wt_b
1388 tr_flux_conv(i,j,kla) = tr_flux_conv(i,j,kla) - (wt_a*tr_flux_3d(i,j,k) + tr_adj_vert_l(i,j,k))
1389 tr_flux_conv(i,j,klb) = tr_flux_conv(i,j,klb) - (wt_b*tr_flux_3d(i,j,k) - tr_adj_vert_l(i,j,k))
1391 if (deep_wt_rv(j)%p(i,k) >= 1.0)
then 1392 tr_flux_conv(i,j+1,krb) = tr_flux_conv(i,j+1,krb) + tr_flux_3d(i,j,k)
1394 kra = k0a_rv(j)%p(i,k)
1395 wt_b = deep_wt_rv(j)%p(i,k) ; wt_a = 1.0 - wt_b
1396 tr_flux_conv(i,j+1,kra) = tr_flux_conv(i,j+1,kra) + &
1397 (wt_a*tr_flux_3d(i,j,k) - tr_adj_vert_r(i,j,k))
1398 tr_flux_conv(i,j+1,krb) = tr_flux_conv(i,j+1,krb) + &
1399 (wt_b*tr_flux_3d(i,j,k) + tr_adj_vert_r(i,j,k))
1402 endif ;
enddo ;
enddo 1404 do k=1,pemax_krho ;
do j=js,je ;
do i=is,ie
1405 if ((g%mask2dT(i,j) > 0.5) .and. (h(i,j,k) > 0.0))
then 1406 tr(m)%t(i,j,k) = tr(m)%t(i,j,k) + tr_flux_conv(i,j,k) / &
1407 (h(i,j,k)*g%areaT(i,j))
1408 tr_flux_conv(i,j,k) = 0.0
1410 enddo ;
enddo ;
enddo 1416 deallocate(deep_wt_lu(j)%p) ;
deallocate(deep_wt_ru(j)%p)
1417 deallocate(hp_lu(j)%p) ;
deallocate(hp_ru(j)%p)
1418 deallocate(k0a_lu(j)%p) ;
deallocate(k0a_ru(j)%p)
1419 deallocate(k0b_lu(j)%p) ;
deallocate(k0b_ru(j)%p)
1423 deallocate(deep_wt_lv(j)%p) ;
deallocate(deep_wt_rv(j)%p)
1424 deallocate(hp_lv(j)%p) ;
deallocate(hp_rv(j)%p)
1425 deallocate(k0a_lv(j)%p) ;
deallocate(k0a_rv(j)%p)
1426 deallocate(k0b_lv(j)%p) ;
deallocate(k0b_rv(j)%p)
1429 end subroutine tracer_epipycnal_ml_diff
1433 subroutine tracer_hor_diff_init(Time, G, US, param_file, diag, EOS, diabatic_CSp, CS)
1434 type(time_type),
target,
intent(in) :: time
1437 type(
diag_ctrl),
target,
intent(inout) :: diag
1438 type(
eos_type),
target,
intent(in) :: eos
1439 type(
diabatic_cs),
pointer,
intent(in) :: diabatic_csp
1444 #include "version_variable.h" 1445 character(len=40) :: mdl =
"MOM_tracer_hor_diff" 1446 character(len=256) :: mesg
1448 if (
associated(cs))
then 1449 call mom_error(warning,
"tracer_hor_diff_init called with associated control structure.")
1455 cs%show_call_tree = calltree_showquery()
1459 call get_param(param_file, mdl,
"KHTR", cs%KhTr, &
1460 "The background along-isopycnal tracer diffusivity.", &
1461 units=
"m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1462 call get_param(param_file, mdl,
"KHTR_SLOPE_CFF", cs%KhTr_Slope_Cff, &
1463 "The scaling coefficient for along-isopycnal tracer "//&
1464 "diffusivity using a shear-based (Visbeck-like) "//&
1465 "parameterization. A non-zero value enables this param.", &
1466 units=
"nondim", default=0.0)
1467 call get_param(param_file, mdl,
"KHTR_MIN", cs%KhTr_Min, &
1468 "The minimum along-isopycnal tracer diffusivity.", &
1469 units=
"m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1470 call get_param(param_file, mdl,
"KHTR_MAX", cs%KhTr_Max, &
1471 "The maximum along-isopycnal tracer diffusivity.", &
1472 units=
"m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1473 call get_param(param_file, mdl,
"KHTR_PASSIVITY_COEFF", cs%KhTr_passivity_coeff, &
1474 "The coefficient that scales deformation radius over "//&
1475 "grid-spacing in passivity, where passivity is the ratio "//&
1476 "between along isopycnal mixing of tracers to thickness mixing. "//&
1477 "A non-zero value enables this parameterization.", &
1478 units=
"nondim", default=0.0)
1479 call get_param(param_file, mdl,
"KHTR_PASSIVITY_MIN", cs%KhTr_passivity_min, &
1480 "The minimum passivity which is the ratio between "//&
1481 "along isopycnal mixing of tracers to thickness mixing.", &
1482 units=
"nondim", default=0.5)
1483 call get_param(param_file, mdl,
"DIFFUSE_ML_TO_INTERIOR", cs%Diffuse_ML_interior, &
1484 "If true, enable epipycnal mixing between the surface "//&
1485 "boundary layer and the interior.", default=.false.)
1486 call get_param(param_file, mdl,
"CHECK_DIFFUSIVE_CFL", cs%check_diffusive_CFL, &
1487 "If true, use enough iterations the diffusion to ensure "//&
1488 "that the diffusive equivalent of the CFL limit is not "//&
1489 "violated. If false, always use the greater of 1 or "//&
1490 "MAX_TR_DIFFUSION_CFL iteration.", default=.false.)
1491 call get_param(param_file, mdl,
"MAX_TR_DIFFUSION_CFL", cs%max_diff_CFL, &
1492 "If positive, locally limit the along-isopycnal tracer "//&
1493 "diffusivity to keep the diffusive CFL locally at or "//&
1494 "below this value. The number of diffusive iterations "//&
1495 "is often this value or the next greater integer.", &
1496 units=
"nondim", default=-1.0)
1497 call get_param(param_file, mdl,
"RECALC_NEUTRAL_SURF", cs%recalc_neutral_surf, &
1498 "If true, then recalculate the neutral surfaces if the \n"//&
1499 "diffusive CFL is exceeded. If false, assume that the \n"//&
1500 "positions of the surfaces do not change \n", default = .false.)
1501 cs%ML_KhTR_scale = 1.0
1502 if (cs%Diffuse_ML_interior)
then 1503 call get_param(param_file, mdl,
"ML_KHTR_SCALE", cs%ML_KhTR_scale, &
1504 "With Diffuse_ML_interior, the ratio of the truly "//&
1505 "horizontal diffusivity in the mixed layer to the "//&
1506 "epipycnal diffusivity. The valid range is 0 to 1.", &
1507 units=
"nondim", default=1.0)
1510 cs%use_neutral_diffusion = neutral_diffusion_init(time, g, us, param_file, diag, eos, &
1511 diabatic_csp, cs%neutral_diffusion_CSp )
1512 if (cs%use_neutral_diffusion .and. cs%Diffuse_ML_interior)
call mom_error(fatal,
"MOM_tracer_hor_diff: "// &
1513 "USE_NEUTRAL_DIFFUSION and DIFFUSE_ML_TO_INTERIOR are mutually exclusive!")
1514 cs%use_lateral_boundary_diffusion = lateral_boundary_diffusion_init(time, g, param_file, diag, diabatic_csp, &
1515 cs%lateral_boundary_diffusion_CSp)
1516 if (cs%use_neutral_diffusion .and. cs%Diffuse_ML_interior)
call mom_error(fatal,
"MOM_tracer_hor_diff: "// &
1517 "USE_LATERAL_BOUNDARY_DIFFUSION and DIFFUSE_ML_TO_INTERIOR are mutually exclusive!")
1519 call get_param(param_file, mdl,
"DEBUG", cs%debug, default=.false.)
1521 id_clock_diffuse = cpu_clock_id(
'(Ocean diffuse tracer)', grain=clock_module)
1522 id_clock_epimix = cpu_clock_id(
'(Ocean epipycnal diffuse tracer)',grain=clock_module)
1523 id_clock_pass = cpu_clock_id(
'(Ocean tracer halo updates)', grain=clock_routine)
1524 id_clock_sync = cpu_clock_id(
'(Ocean tracer global synch)', grain=clock_routine)
1531 cs%id_KhTr_u = register_diag_field(
'ocean_model',
'KHTR_u', diag%axesCu1, time, &
1532 'Epipycnal tracer diffusivity at zonal faces of tracer cell',
'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
1533 cs%id_KhTr_v = register_diag_field(
'ocean_model',
'KHTR_v', diag%axesCv1, time, &
1534 'Epipycnal tracer diffusivity at meridional faces of tracer cell',
'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
1535 cs%id_KhTr_h = register_diag_field(
'ocean_model',
'KHTR_h', diag%axesT1, time, &
1536 'Epipycnal tracer diffusivity at tracer cell center',
'm2 s-1', conversion=us%L_to_m**2*us%s_to_T, &
1537 cmor_field_name=
'diftrelo', &
1538 cmor_standard_name=
'ocean_tracer_epineutral_laplacian_diffusivity', &
1539 cmor_long_name =
'Ocean Tracer Epineutral Laplacian Diffusivity')
1541 cs%id_khdt_x = register_diag_field(
'ocean_model',
'KHDT_x', diag%axesCu1, time, &
1542 'Epipycnal tracer diffusivity operator at zonal faces of tracer cell',
'm2', conversion=us%L_to_m**2)
1543 cs%id_khdt_y = register_diag_field(
'ocean_model',
'KHDT_y', diag%axesCv1, time, &
1544 'Epipycnal tracer diffusivity operator at meridional faces of tracer cell',
'm2', conversion=us%L_to_m**2)
1545 if (cs%check_diffusive_CFL)
then 1546 cs%id_CFL = register_diag_field(
'ocean_model',
'CFL_lateral_diff', diag%axesT1, time,&
1547 'Grid CFL number for lateral/neutral tracer diffusion',
'nondim')
1551 end subroutine tracer_hor_diff_init
1553 subroutine tracer_hor_diff_end(CS)
1556 call neutral_diffusion_end(cs%neutral_diffusion_CSp)
1557 if (
associated(cs))
deallocate(cs)
1559 end subroutine tracer_hor_diff_end
Calculates and applies diffusive fluxes as a parameterization of lateral mixing (non-neutral) by meso...
Ocean grid type. See mom_grid for details.
Calculates density of sea water from T, S and P.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Wraps the MPP cpu clock functions.
This routine drives the diabatic/dianeutral physics for MOM.
Main routine for lateral (along surface or neutral) diffusion of tracers.
Sets parameters for lateral boundary mixing module.
This type is used to exchange information related to the MEKE calculations.
The MOM6 facility to parse input files for runtime parameters.
Variable mixing coefficients.
Do a halo update on a pair of arrays representing the two components of a vector. ...
Variable mixing coefficients.
Make a diagnostic available for averaging or output.
Describes various unit conversion factors.
A type that can be used to create arrays of pointers to 2D arrays.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Type to carry basic tracer information.
A type that can be used to create arrays of pointers to 2D integer arrays.
Routines for error handling and I/O management.
Provides subroutines for quantities specific to the equation of state.
Control structure for this module.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Set up a group of halo updates.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
A column-wise toolbox for implementing neutral diffusion.
The control structure for the MOM_neutral_diffusion module.
Provides a transparent vertical ocean grid type and supporting routines.
The control structure for along-layer and epineutral tracer diffusion.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.
A control structure for the equation of state.