10 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
11 use mom_cpu_clock, only : clock_module_driver, clock_module, clock_routine
15 use mom_diabatic_aux, only : make_frazil, adjust_salt, differential_diffuse_t_s, tridiagts
16 use mom_diabatic_aux, only : find_uv_at_h, diagnosemldbydensitydifference, applyboundaryfluxesinout
23 use mom_diag_mediator, only : diag_copy_diag_to_storage, diag_copy_storage_to_diag
37 use mom_error_handler, only : mom_error, fatal, warning, calltree_showquery,mom_mesg
41 use mom_forcing_type, only : calculatebuoyancyflux2d, forcing_singlepointprint
42 use mom_geothermal, only : geothermal_entraining, geothermal_in_place
53 use mom_cvmix_kpp, only : kpp_nonlocaltransport_temp, kpp_nonlocaltransport_saln
63 use mom_time_manager, only : time_type, real_to_time,
operator(-),
operator(<=)
74 implicit none ;
private 76 #include <MOM_memory.h> 79 public diabatic_driver_init
80 public diabatic_driver_end
81 public extract_diabatic_member
83 public adiabatic_driver_init
94 logical :: use_legacy_diabatic
97 logical :: bulkmixedlayer
99 logical :: use_energetic_pbl
104 logical :: use_kappa_shear
106 logical :: use_cvmix_shear
108 logical :: use_cvmix_ddiff
109 logical :: use_cvmix_conv
111 logical :: use_sponge
115 logical :: use_geothermal
116 logical :: use_int_tides
118 logical :: epbl_is_additive
124 real :: uniform_test_cg
126 logical :: usealealgorithm
129 logical :: aggregate_fw_forcing
139 logical :: massless_match_targets
141 logical :: mix_boundary_tracers
152 real :: minimum_forcing_depth
154 real :: evap_cfl_limit = 0.8
156 integer :: halo_ts_diff = 0
158 logical :: usekpp = .false.
159 logical :: kppispassive
161 logical :: debugconservation
162 logical :: tracer_tridiag
164 logical :: debug_energy_req
166 real :: mlddensitydifference
169 real :: mld_en_vals(3)
172 integer :: id_cg1 = -1
173 integer,
allocatable,
dimension(:) :: id_cn
174 integer :: id_wd = -1, id_ea = -1, id_eb = -1
175 integer :: id_dudt_dia = -1, id_dvdt_dia = -1, id_ea_s = -1, id_eb_s = -1
176 integer :: id_ea_t = -1, id_eb_t = -1
177 integer :: id_kd_heat = -1, id_kd_salt = -1, id_kd_interface = -1, id_kd_epbl = -1
178 integer :: id_tdif = -1, id_tadv = -1, id_sdif = -1, id_sadv = -1
179 integer :: id_mld_003 = -1, id_mld_0125 = -1, id_mld_user = -1, id_mlotstsq = -1
180 integer :: id_mld_en1 = -1, id_mld_en2 = -1, id_mld_en3= -1
181 integer :: id_submln2 = -1
182 integer :: id_hf_dudt_dia_2d = -1, id_hf_dvdt_dia_2d = -1
185 integer :: id_u_predia = -1, id_v_predia = -1, id_h_predia = -1
186 integer :: id_t_predia = -1, id_s_predia = -1, id_e_predia = -1
188 integer :: id_diabatic_diff_temp_tend = -1
189 integer :: id_diabatic_diff_saln_tend = -1
190 integer :: id_diabatic_diff_heat_tend = -1
191 integer :: id_diabatic_diff_salt_tend = -1
192 integer :: id_diabatic_diff_heat_tend_2d = -1
193 integer :: id_diabatic_diff_salt_tend_2d = -1
194 integer :: id_diabatic_diff_h= -1
196 integer :: id_boundary_forcing_h = -1
197 integer :: id_boundary_forcing_h_tendency = -1
198 integer :: id_boundary_forcing_temp_tend = -1
199 integer :: id_boundary_forcing_saln_tend = -1
200 integer :: id_boundary_forcing_heat_tend = -1
201 integer :: id_boundary_forcing_salt_tend = -1
202 integer :: id_boundary_forcing_heat_tend_2d = -1
203 integer :: id_boundary_forcing_salt_tend_2d = -1
205 integer :: id_frazil_h = -1
206 integer :: id_frazil_temp_tend = -1
207 integer :: id_frazil_heat_tend = -1
208 integer :: id_frazil_heat_tend_2d = -1
211 logical :: diabatic_diff_tendency_diag = .false.
212 logical :: boundary_forcing_tendency_diag = .false.
213 logical :: frazil_tendency_diag = .false.
214 real,
allocatable,
dimension(:,:,:) :: frazil_heat_diag
215 real,
allocatable,
dimension(:,:,:) :: frazil_temp_diag
232 type(
kpp_cs),
pointer :: kpp_csp => null()
236 type(group_pass_type) :: pass_hold_eb_ea
237 type(group_pass_type) :: pass_kv
240 real,
allocatable,
dimension(:,:,:) :: kpp_nltheat
241 real,
allocatable,
dimension(:,:,:) :: kpp_nltscalar
242 real,
allocatable,
dimension(:,:,:) :: kpp_buoy_flux
243 real,
allocatable,
dimension(:,:) :: kpp_temp_flux
244 real,
allocatable,
dimension(:,:) :: kpp_salt_flux
246 type(time_type),
pointer :: time
250 integer :: id_clock_entrain, id_clock_mixedlayer, id_clock_set_diffusivity
251 integer :: id_clock_tracers, id_clock_tridiag, id_clock_pass, id_clock_sponge
252 integer :: id_clock_geothermal, id_clock_differential_diff, id_clock_remap
253 integer :: id_clock_kpp
260 subroutine diabatic(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
261 G, GV, US, CS, OBC, WAVES)
264 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
265 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
266 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
269 real,
dimension(:,:),
pointer :: Hml
270 type(
forcing),
intent(inout) :: fluxes
277 real,
intent(in) :: dt
278 type(time_type),
intent(in) :: Time_end
285 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
287 real,
dimension(SZI_(G),SZJ_(G),CS%nMode) :: &
289 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
290 integer :: i, j, k, m, is, ie, js, je, nz
291 logical :: showCallTree
293 if (g%ke == 1)
return 295 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
297 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_diabatic_driver: "// &
298 "Module must be initialized before it is used.")
299 if (dt == 0.0)
call mom_error(fatal,
"MOM_diabatic_driver: "// &
300 "diabatic was called with a zero length timestep.")
301 if (dt < 0.0)
call mom_error(fatal,
"MOM_diabatic_driver: "// &
302 "diabatic was called with a negative timestep.")
304 showcalltree = calltree_showquery()
308 if (cs%id_u_predia > 0)
call post_data(cs%id_u_predia, u, cs%diag)
309 if (cs%id_v_predia > 0)
call post_data(cs%id_v_predia, v, cs%diag)
310 if (cs%id_h_predia > 0)
call post_data(cs%id_h_predia, h, cs%diag)
311 if (cs%id_T_predia > 0)
call post_data(cs%id_T_predia, tv%T, cs%diag)
312 if (cs%id_S_predia > 0)
call post_data(cs%id_S_predia, tv%S, cs%diag)
313 if (cs%id_e_predia > 0)
then 314 call find_eta(h, tv, g, gv, us, eta, eta_to_m=1.0)
315 call post_data(cs%id_e_predia, eta, cs%diag)
319 call mom_state_chksum(
"Start of diabatic ", u, v, h, g, gv, us, haloshift=0)
320 call mom_forcing_chksum(
"Start of diabatic", fluxes, g, us, haloshift=0)
322 if (cs%debugConservation)
call mom_state_stats(
'Start of diabatic', u, v, h, tv%T, tv%S, g, gv, us)
324 if (cs%debug_energy_req) &
325 call diapyc_energy_req_test(h, dt, tv, g, gv, us, cs%diapyc_en_rec_CSp)
327 call cpu_clock_begin(id_clock_set_diffusivity)
328 call set_bbl_tke(u, v, h, fluxes, visc, g, gv, us, cs%set_diff_CSp, obc=obc)
329 call cpu_clock_end(id_clock_set_diffusivity)
334 if (
associated(tv%T) .AND.
associated(tv%frazil))
then 336 call enable_averages(0.5*dt, time_end - real_to_time(0.5*us%T_to_s*dt), cs%diag)
337 if (cs%frazil_tendency_diag)
then 338 do k=1,nz ;
do j=js,je ;
do i=is,ie
339 temp_diag(i,j,k) = tv%T(i,j,k)
340 enddo ;
enddo ;
enddo 343 if (
associated(fluxes%p_surf_full))
then 344 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp, fluxes%p_surf_full, halo=cs%halo_TS_diff)
346 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp, halo=cs%halo_TS_diff)
348 if (showcalltree)
call calltree_waypoint(
"done with 1st make_frazil (diabatic)")
350 if (cs%frazil_tendency_diag)
then 351 call diagnose_frazil_tendency(tv, h, temp_diag, 0.5*dt, g, gv, us, cs)
352 if (cs%id_frazil_h > 0)
call post_data(cs%id_frazil_h, h, cs%diag)
354 call disable_averaging(cs%diag)
356 if (cs%debugConservation)
call mom_state_stats(
'1st make_frazil', u, v, h, tv%T, tv%S, g, gv, us)
358 if (cs%use_int_tides)
then 360 call set_int_tide_input(u, v, h, tv, fluxes, cs%int_tide_input, dt, g, gv, us, &
361 cs%int_tide_input_CSp)
363 if (cs%uniform_test_cg > 0.0)
then 364 do m=1,cs%nMode ; cn_igw(:,:,m) = cs%uniform_test_cg ;
enddo 366 call wave_speeds(h, tv, g, gv, us, cs%nMode, cn_igw, full_halos=.true.)
369 call propagate_int_tide(h, tv, cn_igw, cs%int_tide_input%TKE_itidal_input, cs%int_tide_input%tideamp, &
370 cs%int_tide_input%Nb, dt, g, gv, us, cs%int_tide_CSp)
371 if (showcalltree)
call calltree_waypoint(
"done with propagate_int_tide (diabatic)")
374 if (cs%useALEalgorithm .and. cs%use_legacy_diabatic)
then 375 call diabatic_ale_legacy(u, v, h, tv, hml, fluxes, visc, adp, cdp, dt, time_end, &
376 g, gv, us, cs, waves)
377 elseif (cs%useALEalgorithm)
then 378 call diabatic_ale(u, v, h, tv, hml, fluxes, visc, adp, cdp, dt, time_end, &
379 g, gv, us, cs, waves)
381 call layered_diabatic(u, v, h, tv, hml, fluxes, visc, adp, cdp, dt, time_end, &
382 g, gv, us, cs, waves)
386 call cpu_clock_begin(id_clock_pass)
387 if (
associated(visc%Kv_shear)) &
388 call pass_var(visc%Kv_shear, g%Domain, to_all+omit_corners, halo=1)
389 call cpu_clock_end(id_clock_pass)
391 call disable_averaging(cs%diag)
395 if (
associated(tv%T) .AND.
associated(tv%frazil))
then 396 call enable_averages(0.5*dt, time_end, cs%diag)
397 if (cs%frazil_tendency_diag)
then 398 do k=1,nz ;
do j=js,je ;
do i=is,ie
399 temp_diag(i,j,k) = tv%T(i,j,k)
400 enddo ;
enddo ;
enddo 403 if (
associated(fluxes%p_surf_full))
then 404 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp, fluxes%p_surf_full)
406 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp)
409 if (cs%frazil_tendency_diag)
then 410 call diagnose_frazil_tendency(tv, h, temp_diag, 0.5*dt, g, gv, us, cs)
411 if (cs%id_frazil_h > 0 )
call post_data(cs%id_frazil_h, h, cs%diag)
414 if (showcalltree)
call calltree_waypoint(
"done with 2nd make_frazil (diabatic)")
415 if (cs%debugConservation)
call mom_state_stats(
'2nd make_frazil', u, v, h, tv%T, tv%S, g, gv, us)
416 call disable_averaging(cs%diag)
422 call enable_averages(dt, time_end, cs%diag)
423 if (cs%id_MLD_003 > 0 .or. cs%id_subMLN2 > 0 .or. cs%id_mlotstsq > 0)
then 424 call diagnosemldbydensitydifference(cs%id_MLD_003, h, tv, 0.03*us%kg_m3_to_R, g, gv, us, cs%diag, &
425 id_n2subml=cs%id_subMLN2, id_mldsq=cs%id_mlotstsq, dz_subml=cs%dz_subML_N2)
427 if (cs%id_MLD_0125 > 0)
then 428 call diagnosemldbydensitydifference(cs%id_MLD_0125, h, tv, 0.125*us%kg_m3_to_R, g, gv, us, cs%diag)
430 if (cs%id_MLD_user > 0)
then 431 call diagnosemldbydensitydifference(cs%id_MLD_user, h, tv, cs%MLDdensityDifference, g, gv, us, cs%diag)
433 if ((cs%id_MLD_EN1 > 0) .or. (cs%id_MLD_EN2 > 0) .or. (cs%id_MLD_EN3 > 0))
then 434 call diagnosemldbyenergy((/cs%id_MLD_EN1, cs%id_MLD_EN2, cs%id_MLD_EN3/),&
435 h, tv, g, gv, us, cs%MLD_EN_VALS, cs%diag)
437 if (cs%use_int_tides)
then 438 if (cs%id_cg1 > 0)
call post_data(cs%id_cg1, cn_igw(:,:,1),cs%diag)
439 do m=1,cs%nMode ;
if (cs%id_cn(m) > 0)
call post_data(cs%id_cn(m), cn_igw(:,:,m), cs%diag) ;
enddo 441 call disable_averaging(cs%diag)
443 if (cs%debugConservation)
call mom_state_stats(
'leaving diabatic', u, v, h, tv%T, tv%S, g, gv, us)
445 end subroutine diabatic
450 subroutine diabatic_ale_legacy(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
451 G, GV, US, CS, Waves)
455 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
456 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
457 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
460 real,
dimension(:,:),
pointer :: Hml
461 type(
forcing),
intent(inout) :: fluxes
468 real,
intent(in) :: dt
469 type(time_type),
intent(in) :: Time_end
474 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
493 real,
dimension(SZI_(G),SZJ_(G)) :: &
495 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: h_diag
496 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
497 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: saln_diag
498 real,
dimension(SZI_(G),SZJ_(G)) :: tendency_2d
502 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
507 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
target :: &
517 real,
allocatable,
dimension(:,:) :: &
518 hf_dudt_dia_2d, hf_dvdt_dia_2d
520 logical :: in_boundary(szi_(g))
540 real :: htot(szib_(g))
541 real :: b1(szib_(g)), d1(szib_(g))
542 real :: c1(szib_(g),szk_(g))
548 logical :: showCallTree
549 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, m
554 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
555 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
556 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect*h_neglect
557 kd_heat(:,:,:) = 0.0 ; kd_salt(:,:,:) = 0.0
559 showcalltree = calltree_showquery()
560 if (showcalltree)
call calltree_enter(
"diabatic_ALE_legacy(), MOM_diabatic_driver.F90")
563 call enable_averages(dt, time_end, cs%diag)
565 if (cs%use_geothermal)
then 566 call cpu_clock_begin(id_clock_geothermal)
567 call geothermal_in_place(h, tv, dt, g, gv, us, cs%geothermal_CSp, halo=cs%halo_TS_diff)
568 call cpu_clock_end(id_clock_geothermal)
569 if (showcalltree)
call calltree_waypoint(
"geothermal (diabatic)")
570 if (cs%debugConservation)
call mom_state_stats(
'geothermal', u, v, h, tv%T, tv%S, g, gv, us)
575 call diag_update_remap_grids(cs%diag)
580 if (
associated(cs%optics)) &
581 call set_pen_shortwave(cs%optics, fluxes, g, gv, us, cs%diabatic_aux_CSp, cs%opacity_CSp, cs%tracer_flow_CSp)
583 if (cs%debug)
call mom_state_chksum(
"before find_uv_at_h", u, v, h, g, gv, us, haloshift=0)
585 if (cs%use_kappa_shear .or. cs%use_CVMix_shear)
then 586 if (cs%use_geothermal)
then 589 eatr(:,:,:) = 0.0 ; ebtr(:,:,:) = 0.0
590 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us, eatr, ebtr)
592 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
594 if (showcalltree)
call calltree_waypoint(
"done with find_uv_at_h (diabatic)")
597 call cpu_clock_begin(id_clock_set_diffusivity)
601 call mom_state_chksum(
"before set_diffusivity", u, v, h, g, gv, us, haloshift=cs%halo_TS_diff)
602 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, &
603 visc, dt, g, gv, us, cs%set_diff_CSp, kd_int=kd_int)
604 call cpu_clock_end(id_clock_set_diffusivity)
605 if (showcalltree)
call calltree_waypoint(
"done with set_diffusivity (diabatic)")
611 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
612 kd_salt(i,j,k) = kd_int(i,j,k)
613 kd_heat(i,j,k) = kd_int(i,j,k)
614 enddo ;
enddo ;
enddo 616 if (
associated(visc%Kd_extra_S))
then 618 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
619 kd_salt(i,j,k) = kd_salt(i,j,k) + visc%Kd_extra_S(i,j,k)
620 enddo ;
enddo ;
enddo 622 if (
associated(visc%Kd_extra_T))
then 624 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
625 kd_heat(i,j,k) = kd_heat(i,j,k) + visc%Kd_extra_T(i,j,k)
626 enddo ;
enddo ;
enddo 631 call mom_state_chksum(
"after set_diffusivity ", u, v, h, g, gv, us, haloshift=0)
632 call mom_forcing_chksum(
"after set_diffusivity ", fluxes, g, us, haloshift=0)
633 call mom_thermovar_chksum(
"after set_diffusivity ", tv, g)
634 call hchksum(kd_int,
"after set_diffusivity Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
635 call hchksum(kd_heat,
"after set_diffusivity Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
636 call hchksum(kd_salt,
"after set_diffusivity Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
640 call cpu_clock_begin(id_clock_kpp)
648 call calculatebuoyancyflux2d(g, gv, us, fluxes, cs%optics, h, tv%T, tv%S, tv, &
649 cs%KPP_buoy_flux, cs%KPP_temp_flux, cs%KPP_salt_flux)
652 call kpp_compute_bld(cs%KPP_CSp, g, gv, us, h, tv%T, tv%S, u, v, tv, &
653 fluxes%ustar, cs%KPP_buoy_flux, waves=waves)
655 call kpp_calculate(cs%KPP_CSp, g, gv, us, h, fluxes%ustar, cs%KPP_buoy_flux, kd_heat, &
656 kd_salt, visc%Kv_shear, cs%KPP_NLTheat, cs%KPP_NLTscalar, waves=waves)
658 if (
associated(hml))
then 659 call kpp_get_bld(cs%KPP_CSp, hml(:,:), g, us)
660 call pass_var(hml, g%domain, halo=1)
662 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
665 if (.not.cs%KPPisPassive)
then 667 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
668 kd_int(i,j,k) = min( kd_salt(i,j,k), kd_heat(i,j,k) )
669 enddo ;
enddo ;
enddo 670 if (
associated(visc%Kd_extra_S))
then 672 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
673 visc%Kd_extra_S(i,j,k) = (kd_salt(i,j,k) - kd_int(i,j,k))
674 enddo ;
enddo ;
enddo 676 if (
associated(visc%Kd_extra_T))
then 678 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
679 visc%Kd_extra_T(i,j,k) = (kd_heat(i,j,k) - kd_int(i,j,k))
680 enddo ;
enddo ;
enddo 684 call cpu_clock_end(id_clock_kpp)
685 if (showcalltree)
call calltree_waypoint(
"done with KPP_calculate (diabatic)")
688 call mom_forcing_chksum(
"after KPP", fluxes, g, us, haloshift=0)
689 call mom_thermovar_chksum(
"after KPP", tv, g)
690 call hchksum(kd_heat,
"after KPP Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
691 call hchksum(kd_salt,
"after KPP Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
698 call cpu_clock_begin(id_clock_kpp)
700 call hchksum(cs%KPP_temp_flux,
"before KPP_applyNLT netHeat", g%HI, haloshift=0, scale=gv%H_to_m)
701 call hchksum(cs%KPP_salt_flux,
"before KPP_applyNLT netSalt", g%HI, haloshift=0, scale=gv%H_to_m)
702 call hchksum(cs%KPP_NLTheat,
"before KPP_applyNLT NLTheat", g%HI, haloshift=0)
703 call hchksum(cs%KPP_NLTscalar,
"before KPP_applyNLT NLTscalar", g%HI, haloshift=0)
707 call kpp_nonlocaltransport_temp(cs%KPP_CSp, g, gv, h, cs%KPP_NLTheat, cs%KPP_temp_flux, &
708 us%T_to_s*dt, tv%T, us%Q_to_J_kg*tv%C_p)
709 call kpp_nonlocaltransport_saln(cs%KPP_CSp, g, gv, h, cs%KPP_NLTscalar, cs%KPP_salt_flux, &
711 call cpu_clock_end(id_clock_kpp)
712 if (showcalltree)
call calltree_waypoint(
"done with KPP_applyNonLocalTransport (diabatic)")
713 if (cs%debugConservation)
call mom_state_stats(
'KPP_applyNonLocalTransport', u, v, h, tv%T, tv%S, g, gv, us)
716 call mom_state_chksum(
"after KPP_applyNLT ", u, v, h, g, gv, us, haloshift=0)
717 call mom_forcing_chksum(
"after KPP_applyNLT ", fluxes, g, us, haloshift=0)
718 call mom_thermovar_chksum(
"after KPP_applyNLT ", tv, g)
724 if (
associated(visc%Kd_extra_T) .and.
associated(visc%Kd_extra_S) .and.
associated(tv%T))
then 726 call cpu_clock_begin(id_clock_differential_diff)
727 call differential_diffuse_t_s(h, tv, visc, dt, g, gv)
728 call cpu_clock_end(id_clock_differential_diff)
730 if (showcalltree)
call calltree_waypoint(
"done with differential_diffuse_T_S (diabatic)")
731 if (cs%debugConservation)
call mom_state_stats(
'differential_diffuse_T_S', u, v, h, tv%T, tv%S, g, gv, us)
735 if (.not. cs%useKPP)
then 737 do k=2,nz ;
do j=js,je ;
do i=is,ie
738 kd_heat(i,j,k) = kd_heat(i,j,k) + visc%Kd_extra_T(i,j,k)
739 kd_salt(i,j,k) = kd_salt(i,j,k) + visc%Kd_extra_S(i,j,k)
740 enddo ;
enddo ;
enddo 746 if (cs%use_CVMix_conv)
then 747 call calculate_cvmix_conv(h, tv, g, gv, us, cs%CVMix_conv_csp, hml)
750 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
751 kd_int(i,j,k) = kd_int(i,j,k) + cs%CVMix_conv_csp%kd_conv(i,j,k)
752 visc%Kv_slow(i,j,k) = visc%Kv_slow(i,j,k) + cs%CVMix_conv_csp%kv_conv(i,j,k)
753 enddo ;
enddo ;
enddo 757 do j=js,je ;
do i=is,ie
761 do k=2,nz ;
do j=js,je ;
do i=is,ie
762 hval=1.0/(h_neglect + 0.5*(h(i,j,k-1) + h(i,j,k)))
763 ea_s(i,j,k) = (gv%Z_to_H**2) * dt * hval * kd_int(i,j,k)
764 eb_s(i,j,k-1) = ea_s(i,j,k)
765 ea_t(i,j,k-1) = ea_s(i,j,k-1) ; eb_t(i,j,k-1) = eb_s(i,j,k-1)
766 enddo ;
enddo ;
enddo 767 do j=js,je ;
do i=is,ie
769 ea_t(i,j,nz) = ea_s(i,j,nz) ; eb_t(i,j,nz) = eb_s(i,j,nz)
771 if (showcalltree)
call calltree_waypoint(
"done setting ea,eb from Kd_int (diabatic)")
774 call mom_forcing_chksum(
"after calc_entrain ", fluxes, g, us, haloshift=0)
775 call mom_thermovar_chksum(
"after calc_entrain ", tv, g)
776 call mom_state_chksum(
"after calc_entrain ", u, v, h, g, gv, us, haloshift=0)
777 call hchksum(ea_s,
"after calc_entrain ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
778 call hchksum(eb_s,
"after calc_entrain eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
782 if (cs%boundary_forcing_tendency_diag)
then 783 do k=1,nz ;
do j=js,je ;
do i=is,ie
784 h_diag(i,j,k) = h(i,j,k)
785 temp_diag(i,j,k) = tv%T(i,j,k)
786 saln_diag(i,j,k) = tv%S(i,j,k)
787 enddo ;
enddo ;
enddo 791 call cpu_clock_begin(id_clock_remap)
794 do k=1,nz ;
do j=js,je ;
do i=is,ie
795 h_prebound(i,j,k) = h(i,j,k)
796 enddo ;
enddo ;
enddo 797 if (cs%use_energetic_PBL)
then 799 skinbuoyflux(:,:) = 0.0
800 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
801 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, cs%evap_CFL_limit, &
802 cs%minimum_forcing_depth, ctke, dsv_dt, dsv_ds, skinbuoyflux=skinbuoyflux)
805 call hchksum(ea_t,
"after applyBoundaryFluxes ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
806 call hchksum(eb_t,
"after applyBoundaryFluxes eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
807 call hchksum(ea_s,
"after applyBoundaryFluxes ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
808 call hchksum(eb_s,
"after applyBoundaryFluxes eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
809 call hchksum(ctke,
"after applyBoundaryFluxes cTKE", g%HI, haloshift=0, &
810 scale=us%RZ3_T3_to_W_m2*us%T_to_s)
811 call hchksum(dsv_dt,
"after applyBoundaryFluxes dSV_dT", g%HI, haloshift=0, scale=us%kg_m3_to_R)
812 call hchksum(dsv_ds,
"after applyBoundaryFluxes dSV_dS", g%HI, haloshift=0, scale=us%kg_m3_to_R)
815 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
816 call energetic_pbl(h, u_h, v_h, tv, fluxes, dt, kd_epbl, g, gv, us, &
817 cs%energetic_PBL_CSp, dsv_dt, dsv_ds, ctke, skinbuoyflux, waves=waves)
819 if (
associated(hml))
then 820 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, hml(:,:), g, us)
821 call pass_var(hml, g%domain, halo=1)
823 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
824 elseif (
associated(visc%MLD))
then 825 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, visc%MLD, g, us)
826 call pass_var(visc%MLD, g%domain, halo=1)
830 do k=2,nz ;
do j=js,je ;
do i=is,ie
831 if (cs%ePBL_is_additive)
then 832 kd_add_here = kd_epbl(i,j,k)
833 visc%Kv_shear(i,j,k) = visc%Kv_shear(i,j,k) + cs%ePBL_Prandtl*kd_epbl(i,j,k)
835 kd_add_here = max(kd_epbl(i,j,k) - visc%Kd_shear(i,j,k), 0.0)
836 visc%Kv_shear(i,j,k) = max(visc%Kv_shear(i,j,k), cs%ePBL_Prandtl*kd_epbl(i,j,k))
839 ent_int = kd_add_here * (gv%Z_to_H**2 * dt) / &
840 (0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect)
841 eb_s(i,j,k-1) = eb_s(i,j,k-1) + ent_int
842 ea_s(i,j,k) = ea_s(i,j,k) + ent_int
843 kd_int(i,j,k) = kd_int(i,j,k) + kd_add_here
846 kd_heat(i,j,k) = kd_heat(i,j,k) + kd_int(i,j,k)
847 kd_salt(i,j,k) = kd_salt(i,j,k) + kd_int(i,j,k)
849 enddo ;
enddo ;
enddo 852 call hchksum(ea_t,
"after ePBL ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
853 call hchksum(eb_t,
"after ePBL eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
854 call hchksum(ea_s,
"after ePBL ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
855 call hchksum(eb_s,
"after ePBL eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
856 call hchksum(kd_epbl,
"after ePBL Kd_ePBL", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
860 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
861 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, &
862 cs%evap_CFL_limit, cs%minimum_forcing_depth)
869 if (cs%boundary_forcing_tendency_diag)
then 870 call diagnose_boundary_forcing_tendency(tv, h, temp_diag, saln_diag, h_diag, dt, g, gv, us, cs)
871 if (cs%id_boundary_forcing_h > 0)
call post_data(cs%id_boundary_forcing_h, h, cs%diag, alt_h = h_diag)
874 call diag_update_remap_grids(cs%diag)
875 call cpu_clock_end(id_clock_remap)
877 call mom_forcing_chksum(
"after applyBoundaryFluxes ", fluxes, g, us, haloshift=0)
878 call mom_thermovar_chksum(
"after applyBoundaryFluxes ", tv, g)
879 call mom_state_chksum(
"after applyBoundaryFluxes ", u, v, h, g, gv, us, haloshift=0)
881 if (showcalltree)
call calltree_waypoint(
"done with applyBoundaryFluxes (diabatic)")
882 if (cs%debugConservation)
call mom_state_stats(
'applyBoundaryFluxes', u, v, h, tv%T, tv%S, g, gv, us)
893 hold(i,j,1) = h(i,j,1)
895 hold(i,j,nz) = h(i,j,nz)
897 if (h(i,j,1) <= 0.0) h(i,j,1) = gv%Angstrom_H
898 if (h(i,j,nz) <= 0.0) h(i,j,nz) = gv%Angstrom_H
900 do k=2,nz-1 ;
do i=is,ie
901 hold(i,j,k) = h(i,j,k)
904 if (h(i,j,k) <= 0.0) h(i,j,k) = gv%Angstrom_H
908 call diag_update_remap_grids(cs%diag)
911 call mom_state_chksum(
"after negative check ", u, v, h, g, gv, us, haloshift=0)
912 call mom_forcing_chksum(
"after negative check ", fluxes, g, us, haloshift=0)
913 call mom_thermovar_chksum(
"after negative check ", tv, g)
915 if (showcalltree)
call calltree_waypoint(
"done with h=ea-eb (diabatic)")
916 if (cs%debugConservation)
call mom_state_stats(
'h=ea-eb', u, v, h, tv%T, tv%S, g, gv, us)
919 if (
associated(tv%T))
then 922 call hchksum(ea_t,
"before triDiagTS ea_t ", g%HI, haloshift=0, scale=gv%H_to_m)
923 call hchksum(eb_t,
"before triDiagTS eb_t ", g%HI, haloshift=0, scale=gv%H_to_m)
924 call hchksum(ea_s,
"before triDiagTS ea_s ", g%HI, haloshift=0, scale=gv%H_to_m)
925 call hchksum(eb_s,
"before triDiagTS eb_s ", g%HI, haloshift=0, scale=gv%H_to_m)
928 call cpu_clock_begin(id_clock_tridiag)
935 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
936 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
938 if (cs%diabatic_diff_tendency_diag)
then 939 do k=1,nz ;
do j=js,je ;
do i=is,ie
940 temp_diag(i,j,k) = tv%T(i,j,k)
941 saln_diag(i,j,k) = tv%S(i,j,k)
942 enddo ;
enddo ;
enddo 946 do k=1,nz ;
do j=js,je ;
do i=is,ie
947 ea_t(i,j,k) = ea_s(i,j,k) ; eb_t(i,j,k) = eb_s(i,j,k)
948 enddo ;
enddo ;
enddo 949 if (cs%tracer_tridiag)
then 950 call tracer_vertdiff(hold, ea_t, eb_t, dt, tv%T, g, gv)
951 call tracer_vertdiff(hold, ea_s, eb_s, dt, tv%S, g, gv)
953 call tridiagts(g, gv, is, ie, js, je, hold, ea_s, eb_s, tv%T, tv%S)
960 if (cs%diabatic_diff_tendency_diag)
then 961 call diagnose_diabatic_diff_tendency(tv, hold, temp_diag, saln_diag, dt, g, gv, us, cs)
962 if (cs%id_diabatic_diff_h > 0)
call post_data(cs%id_diabatic_diff_h, hold, cs%diag, alt_h = hold)
965 call cpu_clock_end(id_clock_tridiag)
967 if (showcalltree)
call calltree_waypoint(
"done with triDiagTS (diabatic)")
971 if (cs%debugConservation)
call mom_state_stats(
'triDiagTS', u, v, h, tv%T, tv%S, g, gv, us)
974 call mom_state_chksum(
"after mixed layer ", u, v, h, g, gv, us, haloshift=0)
975 call mom_thermovar_chksum(
"after mixed layer ", tv, g)
980 if (
associated(adp%du_dt_dia) .or.
associated(adp%dv_dt_dia)) &
982 call pass_var(h, g%domain, to_west+to_south+omit_corners, halo=1)
983 call diag_update_remap_grids(cs%diag)
987 if ((cs%id_Tdif > 0) .or. (cs%id_Tadv > 0))
then 988 do j=js,je ;
do i=is,ie
989 tdif_flx(i,j,1) = 0.0 ; tdif_flx(i,j,nz+1) = 0.0
990 tadv_flx(i,j,1) = 0.0 ; tadv_flx(i,j,nz+1) = 0.0
993 do k=2,nz ;
do j=js,je ;
do i=is,ie
994 tdif_flx(i,j,k) = (idt * 0.5*(ea_t(i,j,k) + eb_t(i,j,k-1))) * &
995 (tv%T(i,j,k-1) - tv%T(i,j,k))
996 tadv_flx(i,j,k) = (idt * (ea_t(i,j,k) - eb_t(i,j,k-1))) * &
997 0.5*(tv%T(i,j,k-1) + tv%T(i,j,k))
998 enddo ;
enddo ;
enddo 1000 if ((cs%id_Sdif > 0) .or. (cs%id_Sadv > 0))
then 1001 do j=js,je ;
do i=is,ie
1002 sdif_flx(i,j,1) = 0.0 ; sdif_flx(i,j,nz+1) = 0.0
1003 sadv_flx(i,j,1) = 0.0 ; sadv_flx(i,j,nz+1) = 0.0
1006 do k=2,nz ;
do j=js,je ;
do i=is,ie
1007 sdif_flx(i,j,k) = (idt * 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))) * &
1008 (tv%S(i,j,k-1) - tv%S(i,j,k))
1009 sadv_flx(i,j,k) = (idt * (ea_s(i,j,k) - eb_s(i,j,k-1))) * &
1010 0.5*(tv%S(i,j,k-1) + tv%S(i,j,k))
1011 enddo ;
enddo ;
enddo 1015 call cpu_clock_begin(id_clock_tracers)
1017 if (cs%mix_boundary_tracers)
then 1018 tr_ea_bbl = gv%Z_to_H * sqrt(dt*cs%Kd_BBL_tr)
1022 ebtr(i,j,nz) = eb_s(i,j,nz)
1024 in_boundary(i) = (g%mask2dT(i,j) > 0.0)
1026 do k=nz,2,-1 ;
do i=is,ie
1027 if (in_boundary(i))
then 1028 htot(i) = htot(i) + h(i,j,k)
1037 add_ent = ((dt * cs%Kd_min_tr) * gv%Z_to_H**2) * &
1038 ((h(i,j,k-1)+h(i,j,k)+h_neglect) / &
1039 (h(i,j,k-1)*h(i,j,k)+h_neglect2)) - &
1040 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))
1041 if (htot(i) < tr_ea_bbl)
then 1042 add_ent = max(0.0, add_ent, &
1043 (tr_ea_bbl - htot(i)) - min(ea_s(i,j,k),eb_s(i,j,k-1)))
1044 elseif (add_ent < 0.0)
then 1045 add_ent = 0.0 ; in_boundary(i) = .false.
1048 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1049 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1051 ebtr(i,j,k-1) = eb_s(i,j,k-1) ; eatr(i,j,k) = ea_s(i,j,k)
1054 if (
associated(visc%Kd_extra_S))
then ;
if (visc%Kd_extra_S(i,j,k) > 0.0)
then 1055 add_ent = ((dt * visc%Kd_extra_S(i,j,k)) * gv%Z_to_H**2) / &
1056 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + h_neglect)
1057 ebtr(i,j,k-1) = ebtr(i,j,k-1) + add_ent
1058 eatr(i,j,k) = eatr(i,j,k) + add_ent
1061 do i=is,ie ; eatr(i,j,1) = ea_s(i,j,1) ;
enddo 1067 call call_tracer_column_fns(h_prebound, h, ea_s, eb_s, fluxes, hml, dt, g, gv, us, tv, &
1068 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1069 evap_cfl_limit = cs%evap_CFL_limit, &
1070 minimum_forcing_depth=cs%minimum_forcing_depth)
1072 elseif (
associated(visc%Kd_extra_S))
then 1074 do j=js,je ;
do i=is,ie
1075 ebtr(i,j,nz) = eb_s(i,j,nz) ; eatr(i,j,1) = ea_s(i,j,1)
1078 do k=nz,2,-1 ;
do j=js,je ;
do i=is,ie
1079 if (visc%Kd_extra_S(i,j,k) > 0.0)
then 1080 add_ent = ((dt * visc%Kd_extra_S(i,j,k)) * gv%Z_to_H**2) / &
1081 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + h_neglect)
1085 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1086 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1087 enddo ;
enddo ;
enddo 1090 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1091 cs%optics, cs%tracer_flow_CSp, cs%debug,&
1092 evap_cfl_limit = cs%evap_CFL_limit, &
1093 minimum_forcing_depth=cs%minimum_forcing_depth)
1096 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1097 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1098 evap_cfl_limit = cs%evap_CFL_limit, &
1099 minimum_forcing_depth=cs%minimum_forcing_depth)
1102 call cpu_clock_end(id_clock_tracers)
1105 if (cs%use_sponge)
then 1106 call cpu_clock_begin(id_clock_sponge)
1107 if (
associated(cs%ALE_sponge_CSp))
then 1108 call apply_ale_sponge(h, dt, g, gv, us, cs%ALE_sponge_CSp, cs%Time)
1111 call cpu_clock_end(id_clock_sponge)
1114 call mom_thermovar_chksum(
"apply_sponge ", tv, g)
1118 call disable_averaging(cs%diag)
1120 call enable_averages(dt, time_end, cs%diag)
1122 if (cs%id_Kd_interface > 0)
call post_data(cs%id_Kd_interface, kd_int, cs%diag)
1123 if (cs%id_Kd_heat > 0)
call post_data(cs%id_Kd_heat, kd_heat, cs%diag)
1124 if (cs%id_Kd_salt > 0)
call post_data(cs%id_Kd_salt, kd_salt, cs%diag)
1125 if (cs%id_Kd_ePBL > 0)
call post_data(cs%id_Kd_ePBL, kd_epbl, cs%diag)
1127 if (cs%id_ea > 0)
call post_data(cs%id_ea, ea_s, cs%diag)
1128 if (cs%id_eb > 0)
call post_data(cs%id_eb, eb_s, cs%diag)
1129 if (cs%id_ea_t > 0)
call post_data(cs%id_ea_t, ea_t, cs%diag)
1130 if (cs%id_eb_t > 0)
call post_data(cs%id_eb_t, eb_t, cs%diag)
1131 if (cs%id_ea_s > 0)
call post_data(cs%id_ea_s, ea_s, cs%diag)
1132 if (cs%id_eb_s > 0)
call post_data(cs%id_eb_s, eb_s, cs%diag)
1134 if (cs%id_dudt_dia > 0)
call post_data(cs%id_dudt_dia, adp%du_dt_dia, cs%diag)
1135 if (cs%id_dvdt_dia > 0)
call post_data(cs%id_dvdt_dia, adp%dv_dt_dia, cs%diag)
1136 if (cs%id_wd > 0)
call post_data(cs%id_wd, cdp%diapyc_vel, cs%diag)
1138 if (cs%id_Tdif > 0)
call post_data(cs%id_Tdif, tdif_flx, cs%diag)
1139 if (cs%id_Tadv > 0)
call post_data(cs%id_Tadv, tadv_flx, cs%diag)
1140 if (cs%id_Sdif > 0)
call post_data(cs%id_Sdif, sdif_flx, cs%diag)
1141 if (cs%id_Sadv > 0)
call post_data(cs%id_Sadv, sadv_flx, cs%diag)
1144 if (cs%id_hf_dudt_dia_2d > 0)
then 1145 allocate(hf_dudt_dia_2d(g%IsdB:g%IedB,g%jsd:g%jed))
1146 hf_dudt_dia_2d(:,:) = 0.0
1147 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1148 hf_dudt_dia_2d(i,j) = hf_dudt_dia_2d(i,j) + adp%du_dt_dia(i,j,k) * adp%diag_hfrac_u(i,j,k)
1149 enddo ;
enddo ;
enddo 1150 call post_data(cs%id_hf_dudt_dia_2d, hf_dudt_dia_2d, cs%diag)
1151 deallocate(hf_dudt_dia_2d)
1154 if (cs%id_hf_dvdt_dia_2d > 0)
then 1155 allocate(hf_dvdt_dia_2d(g%isd:g%ied,g%JsdB:g%JedB))
1156 hf_dvdt_dia_2d(:,:) = 0.0
1157 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1158 hf_dvdt_dia_2d(i,j) = hf_dvdt_dia_2d(i,j) + adp%dv_dt_dia(i,j,k) * adp%diag_hfrac_v(i,j,k)
1159 enddo ;
enddo ;
enddo 1160 call post_data(cs%id_hf_dvdt_dia_2d, hf_dvdt_dia_2d, cs%diag)
1161 deallocate(hf_dvdt_dia_2d)
1164 call disable_averaging(cs%diag)
1166 if (showcalltree)
call calltree_leave(
"diabatic_ALE_legacy()")
1168 end subroutine diabatic_ale_legacy
1173 subroutine diabatic_ale(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
1174 G, GV, US, CS, Waves)
1178 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
1179 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
1180 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
1183 real,
dimension(:,:),
pointer :: Hml
1184 type(
forcing),
intent(inout) :: fluxes
1191 real,
intent(in) :: dt
1192 type(time_type),
intent(in) :: Time_end
1197 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
1212 real,
dimension(SZI_(G),SZJ_(G)) :: &
1214 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: h_diag
1215 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
1216 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: saln_diag
1217 real,
dimension(SZI_(G),SZJ_(G)) :: tendency_2d
1221 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
1226 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
target :: &
1236 real,
allocatable,
dimension(:,:) :: &
1237 hf_dudt_dia_2d, hf_dvdt_dia_2d
1239 logical :: in_boundary(szi_(g))
1259 real :: htot(szib_(g))
1260 real :: b1(szib_(g)), d1(szib_(g))
1261 real :: c1(szib_(g),szk_(g))
1267 logical :: showCallTree
1268 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, m
1270 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1271 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1272 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect*h_neglect
1273 kd_heat(:,:,:) = 0.0 ; kd_salt(:,:,:) = 0.0
1274 ea_s(:,:,:) = 0.0; eb_s(:,:,:) = 0.0; ea_t(:,:,:) = 0.0; eb_t(:,:,:) = 0.0
1276 showcalltree = calltree_showquery()
1277 if (showcalltree)
call calltree_enter(
"diabatic_ALE(), MOM_diabatic_driver.F90")
1279 if (.not. (cs%useALEalgorithm))
call mom_error(fatal,
"MOM_diabatic_driver: "// &
1280 "The ALE algorithm must be enabled when using MOM_diabatic_driver.")
1283 call enable_averages(dt, time_end, cs%diag)
1285 if (cs%use_geothermal)
then 1286 call cpu_clock_begin(id_clock_geothermal)
1287 call geothermal_in_place(h, tv, dt, g, gv, us, cs%geothermal_CSp, halo=cs%halo_TS_diff)
1288 call cpu_clock_end(id_clock_geothermal)
1289 if (showcalltree)
call calltree_waypoint(
"geothermal (diabatic)")
1290 if (cs%debugConservation)
call mom_state_stats(
'geothermal', u, v, h, tv%T, tv%S, g, gv, us)
1295 call diag_update_remap_grids(cs%diag)
1300 if (
associated(cs%optics)) &
1301 call set_pen_shortwave(cs%optics, fluxes, g, gv, us, cs%diabatic_aux_CSp, cs%opacity_CSp, cs%tracer_flow_CSp)
1303 if (cs%debug)
call mom_state_chksum(
"before find_uv_at_h", u, v, h, g, gv, us, haloshift=0)
1305 if (cs%use_kappa_shear .or. cs%use_CVMix_shear)
then 1306 if (cs%use_geothermal)
then 1309 eatr(:,:,:) = 0.0 ; ebtr(:,:,:) = 0.0
1310 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us, eatr, ebtr)
1312 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
1314 if (showcalltree)
call calltree_waypoint(
"done with find_uv_at_h (diabatic)")
1317 call cpu_clock_begin(id_clock_set_diffusivity)
1321 call mom_state_chksum(
"before set_diffusivity", u, v, h, g, gv, us, haloshift=cs%halo_TS_diff)
1322 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, visc, dt, g, gv, us, &
1323 cs%set_diff_CSp, kd_int=kd_int)
1324 call cpu_clock_end(id_clock_set_diffusivity)
1325 if (showcalltree)
call calltree_waypoint(
"done with set_diffusivity (diabatic)")
1328 call mom_state_chksum(
"after set_diffusivity ", u, v, h, g, gv, us, haloshift=0)
1329 call mom_forcing_chksum(
"after set_diffusivity ", fluxes, g, us, haloshift=0)
1330 call mom_thermovar_chksum(
"after set_diffusivity ", tv, g)
1331 call hchksum(kd_int,
"after set_diffusivity Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1337 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1338 kd_salt(i,j,k) = kd_int(i,j,k)
1339 kd_heat(i,j,k) = kd_int(i,j,k)
1340 enddo ;
enddo ;
enddo 1342 if (
associated(visc%Kd_extra_S))
then 1344 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1345 kd_salt(i,j,k) = kd_salt(i,j,k) + visc%Kd_extra_S(i,j,k)
1346 enddo ;
enddo ;
enddo 1348 if (
associated(visc%Kd_extra_T))
then 1350 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1351 kd_heat(i,j,k) = kd_heat(i,j,k) + visc%Kd_extra_T(i,j,k)
1352 enddo ;
enddo ;
enddo 1356 call hchksum(kd_heat,
"after set_diffusivity Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1357 call hchksum(kd_salt,
"after set_diffusivity Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1361 call cpu_clock_begin(id_clock_kpp)
1363 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1364 visc%Kv_shear(i,j,k) = visc%Kv_shear(i,j,k) + visc%Kv_slow(i,j,k)
1365 enddo ;
enddo ;
enddo 1372 call calculatebuoyancyflux2d(g, gv, us, fluxes, cs%optics, h, tv%T, tv%S, tv, &
1373 cs%KPP_buoy_flux, cs%KPP_temp_flux, cs%KPP_salt_flux)
1376 call kpp_compute_bld(cs%KPP_CSp, g, gv, us, h, tv%T, tv%S, u, v, tv, &
1377 fluxes%ustar, cs%KPP_buoy_flux, waves=waves)
1379 call kpp_calculate(cs%KPP_CSp, g, gv, us, h, fluxes%ustar, cs%KPP_buoy_flux, kd_heat, &
1380 kd_salt, visc%Kv_shear, cs%KPP_NLTheat, cs%KPP_NLTscalar, waves=waves)
1382 if (
associated(hml))
then 1383 call kpp_get_bld(cs%KPP_CSp, hml(:,:), g, us)
1384 call pass_var(hml, g%domain, halo=1)
1386 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
1389 call cpu_clock_end(id_clock_kpp)
1390 if (showcalltree)
call calltree_waypoint(
"done with KPP_calculate (diabatic)")
1393 call mom_forcing_chksum(
"after KPP", fluxes, g, us, haloshift=0)
1394 call mom_thermovar_chksum(
"after KPP", tv, g)
1395 call hchksum(kd_heat,
"after KPP Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1396 call hchksum(kd_salt,
"after KPP Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1403 call cpu_clock_begin(id_clock_kpp)
1405 call hchksum(cs%KPP_temp_flux,
"before KPP_applyNLT netHeat", g%HI, haloshift=0, scale=gv%H_to_m)
1406 call hchksum(cs%KPP_salt_flux,
"before KPP_applyNLT netSalt", g%HI, haloshift=0, scale=gv%H_to_m)
1407 call hchksum(cs%KPP_NLTheat,
"before KPP_applyNLT NLTheat", g%HI, haloshift=0)
1408 call hchksum(cs%KPP_NLTscalar,
"before KPP_applyNLT NLTscalar", g%HI, haloshift=0)
1412 call kpp_nonlocaltransport_temp(cs%KPP_CSp, g, gv, h, cs%KPP_NLTheat, cs%KPP_temp_flux, &
1413 us%T_to_s*dt, tv%T, us%Q_to_J_kg*tv%C_p)
1414 call kpp_nonlocaltransport_saln(cs%KPP_CSp, g, gv, h, cs%KPP_NLTscalar, cs%KPP_salt_flux, &
1416 call cpu_clock_end(id_clock_kpp)
1417 if (showcalltree)
call calltree_waypoint(
"done with KPP_applyNonLocalTransport (diabatic)")
1418 if (cs%debugConservation)
call mom_state_stats(
'KPP_applyNonLocalTransport', u, v, h, tv%T, tv%S, g, gv, us)
1421 call mom_state_chksum(
"after KPP_applyNLT ", u, v, h, g, gv, us, haloshift=0)
1422 call mom_forcing_chksum(
"after KPP_applyNLT ", fluxes, g, us, haloshift=0)
1423 call mom_thermovar_chksum(
"after KPP_applyNLT ", tv, g)
1429 if (
associated(visc%Kd_extra_T) .and.
associated(visc%Kd_extra_S) .and.
associated(tv%T) .and. &
1430 (.not.cs%use_CVMix_ddiff))
then 1432 call cpu_clock_begin(id_clock_differential_diff)
1433 call differential_diffuse_t_s(h, tv, visc, dt, g, gv)
1434 call cpu_clock_end(id_clock_differential_diff)
1436 if (showcalltree)
call calltree_waypoint(
"done with differential_diffuse_T_S (diabatic)")
1437 if (cs%debugConservation)
call mom_state_stats(
'differential_diffuse_T_S', u, v, h, tv%T, tv%S, g, gv, us)
1441 if (.not. cs%useKPP)
then 1443 do k=2,nz ;
do j=js,je ;
do i=is,ie
1444 kd_heat(i,j,k) = kd_heat(i,j,k) + visc%Kd_extra_T(i,j,k)
1445 kd_salt(i,j,k) = kd_salt(i,j,k) + visc%Kd_extra_S(i,j,k)
1446 enddo ;
enddo ;
enddo 1452 if (cs%use_CVMix_conv)
then 1453 call calculate_cvmix_conv(h, tv, g, gv, us, cs%CVMix_conv_csp, hml)
1456 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1457 kd_heat(i,j,k) = kd_heat(i,j,k) + cs%CVMix_conv_csp%kd_conv(i,j,k)
1458 kd_salt(i,j,k) = kd_salt(i,j,k) + cs%CVMix_conv_csp%kd_conv(i,j,k)
1460 visc%Kv_shear(i,j,k) = visc%Kv_shear(i,j,k) + cs%CVMix_conv_csp%kv_conv(i,j,k)
1462 visc%Kv_slow(i,j,k) = visc%Kv_slow(i,j,k) + cs%CVMix_conv_csp%kv_conv(i,j,k)
1464 enddo ;
enddo ;
enddo 1468 if (cs%boundary_forcing_tendency_diag)
then 1469 do k=1,nz ;
do j=js,je ;
do i=is,ie
1470 h_diag(i,j,k) = h(i,j,k)
1471 temp_diag(i,j,k) = tv%T(i,j,k)
1472 saln_diag(i,j,k) = tv%S(i,j,k)
1473 enddo ;
enddo ;
enddo 1477 call cpu_clock_begin(id_clock_remap)
1480 do k=1,nz ;
do j=js,je ;
do i=is,ie
1481 h_prebound(i,j,k) = h(i,j,k)
1482 enddo ;
enddo ;
enddo 1483 if (cs%use_energetic_PBL)
then 1485 skinbuoyflux(:,:) = 0.0
1486 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
1487 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, cs%evap_CFL_limit, &
1488 cs%minimum_forcing_depth, ctke, dsv_dt, dsv_ds, skinbuoyflux=skinbuoyflux)
1491 call hchksum(ea_t,
"after applyBoundaryFluxes ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
1492 call hchksum(eb_t,
"after applyBoundaryFluxes eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
1493 call hchksum(ea_s,
"after applyBoundaryFluxes ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
1494 call hchksum(eb_s,
"after applyBoundaryFluxes eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
1495 call hchksum(ctke,
"after applyBoundaryFluxes cTKE", g%HI, haloshift=0, &
1496 scale=us%RZ3_T3_to_W_m2*us%T_to_s)
1497 call hchksum(dsv_dt,
"after applyBoundaryFluxes dSV_dT", g%HI, haloshift=0, scale=us%kg_m3_to_R)
1498 call hchksum(dsv_ds,
"after applyBoundaryFluxes dSV_dS", g%HI, haloshift=0, scale=us%kg_m3_to_R)
1501 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
1502 call energetic_pbl(h, u_h, v_h, tv, fluxes, dt, kd_epbl, g, gv, us, &
1503 cs%energetic_PBL_CSp, dsv_dt, dsv_ds, ctke, skinbuoyflux, waves=waves)
1505 if (
associated(hml))
then 1506 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, hml(:,:), g, us)
1507 call pass_var(hml, g%domain, halo=1)
1509 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
1510 elseif (
associated(visc%MLD))
then 1511 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, visc%MLD, g, us)
1512 call pass_var(visc%MLD, g%domain, halo=1)
1516 do k=2,nz ;
do j=js,je ;
do i=is,ie
1517 if (cs%ePBL_is_additive)
then 1518 kd_add_here = kd_epbl(i,j,k)
1519 visc%Kv_shear(i,j,k) = visc%Kv_shear(i,j,k) + cs%ePBL_Prandtl*kd_epbl(i,j,k)
1521 kd_add_here = max(kd_epbl(i,j,k) - visc%Kd_shear(i,j,k), 0.0)
1522 visc%Kv_shear(i,j,k) = max(visc%Kv_shear(i,j,k), cs%ePBL_Prandtl*kd_epbl(i,j,k))
1525 kd_heat(i,j,k) = kd_heat(i,j,k) + kd_add_here
1526 kd_salt(i,j,k) = kd_salt(i,j,k) + kd_add_here
1528 enddo ;
enddo ;
enddo 1531 call hchksum(ea_t,
"after ePBL ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
1532 call hchksum(eb_t,
"after ePBL eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
1533 call hchksum(ea_s,
"after ePBL ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
1534 call hchksum(eb_s,
"after ePBL eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
1535 call hchksum(kd_epbl,
"after ePBL Kd_ePBL", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1539 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
1540 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, &
1541 cs%evap_CFL_limit, cs%minimum_forcing_depth)
1548 if (cs%boundary_forcing_tendency_diag)
then 1549 call diagnose_boundary_forcing_tendency(tv, h, temp_diag, saln_diag, h_diag, dt, g, gv, us, cs)
1550 if (cs%id_boundary_forcing_h > 0)
call post_data(cs%id_boundary_forcing_h, h, cs%diag, alt_h = h_diag)
1553 call diag_update_remap_grids(cs%diag)
1554 call cpu_clock_end(id_clock_remap)
1556 call mom_forcing_chksum(
"after applyBoundaryFluxes ", fluxes, g, us, haloshift=0)
1557 call mom_thermovar_chksum(
"after applyBoundaryFluxes ", tv, g)
1558 call mom_state_chksum(
"after applyBoundaryFluxes ", u, v, h, g, gv, us, haloshift=0)
1560 if (showcalltree)
call calltree_waypoint(
"done with applyBoundaryFluxes (diabatic)")
1561 if (cs%debugConservation)
call mom_state_stats(
'applyBoundaryFluxes', u, v, h, tv%T, tv%S, g, gv, us)
1563 if (showcalltree)
call calltree_waypoint(
"done with h=ea-eb (diabatic)")
1564 if (cs%debugConservation)
call mom_state_stats(
'h=ea-eb', u, v, h, tv%T, tv%S, g, gv, us)
1567 if (
associated(tv%T))
then 1570 call hchksum(ea_t,
"before triDiagTS ea_t ", g%HI, haloshift=0, scale=gv%H_to_m)
1571 call hchksum(eb_t,
"before triDiagTS eb_t ", g%HI, haloshift=0, scale=gv%H_to_m)
1572 call hchksum(ea_s,
"before triDiagTS ea_s ", g%HI, haloshift=0, scale=gv%H_to_m)
1573 call hchksum(eb_s,
"before triDiagTS eb_s ", g%HI, haloshift=0, scale=gv%H_to_m)
1576 call cpu_clock_begin(id_clock_tridiag)
1583 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
1584 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
1586 if (cs%diabatic_diff_tendency_diag)
then 1587 do k=1,nz ;
do j=js,je ;
do i=is,ie
1588 temp_diag(i,j,k) = tv%T(i,j,k)
1589 saln_diag(i,j,k) = tv%S(i,j,k)
1590 enddo ;
enddo ;
enddo 1594 do j=js,je ;
do i=is,ie
1595 ea_t(i,j,1) = 0.; ea_s(i,j,1) = 0.
1599 do k=2,nz ;
do j=js,je ;
do i=is,ie
1600 hval = 1.0 / (h_neglect + 0.5*(h(i,j,k-1) + h(i,j,k)))
1601 ea_t(i,j,k) = (gv%Z_to_H**2) * dt * hval * kd_heat(i,j,k)
1602 eb_t(i,j,k-1) = ea_t(i,j,k)
1603 ea_s(i,j,k) = (gv%Z_to_H**2) * dt * hval * kd_salt(i,j,k)
1604 eb_s(i,j,k-1) = ea_s(i,j,k)
1605 enddo ;
enddo ;
enddo 1606 do j=js,je ;
do i=is,ie
1607 eb_t(i,j,nz) = 0. ; eb_s(i,j,nz) = 0.
1609 if (showcalltree)
call calltree_waypoint(
"done setting ea_t,ea_s,eb_t,eb_s from Kd_heat" //&
1610 "and Kd_salt (diabatic)")
1616 ea_t(i,js-1,k) = 0.0 ; eb_t(i,js-1,k) = 0.0
1617 ea_s(i,js-1,k) = 0.0 ; eb_s(i,js-1,k) = 0.0
1618 ea_t(i,je+1,k) = 0.0 ; eb_t(i,je+1,k) = 0.0
1619 ea_s(i,je+1,k) = 0.0 ; eb_s(i,je+1,k) = 0.0
1622 ea_t(is-1,j,k) = 0.0 ; eb_t(is-1,j,k) = 0.0
1623 ea_s(is-1,j,k) = 0.0 ; eb_s(is-1,j,k) = 0.0
1624 ea_t(ie+1,j,k) = 0.0 ; eb_t(ie+1,j,k) = 0.0
1625 ea_s(ie+1,j,k) = 0.0 ; eb_s(ie+1,j,k) = 0.0
1630 call tracer_vertdiff(h, ea_t, eb_t, dt, tv%T, g, gv)
1631 call tracer_vertdiff(h, ea_s, eb_s, dt, tv%S, g, gv)
1635 if (cs%diabatic_diff_tendency_diag)
then 1636 call diagnose_diabatic_diff_tendency(tv, h, temp_diag, saln_diag, dt, g, gv, us, cs)
1638 call cpu_clock_end(id_clock_tridiag)
1640 if (showcalltree)
call calltree_waypoint(
"done with triDiagTS (diabatic)")
1644 if (cs%debugConservation)
call mom_state_stats(
'triDiagTS', u, v, h, tv%T, tv%S, g, gv, us)
1647 call mom_state_chksum(
"after mixed layer ", u, v, h, g, gv, us, haloshift=0)
1648 call mom_thermovar_chksum(
"after mixed layer ", tv, g)
1653 call diag_update_remap_grids(cs%diag)
1657 if ((cs%id_Tdif > 0) .or. (cs%id_Tadv > 0))
then 1658 do j=js,je ;
do i=is,ie
1659 tdif_flx(i,j,1) = 0.0 ; tdif_flx(i,j,nz+1) = 0.0
1660 tadv_flx(i,j,1) = 0.0 ; tadv_flx(i,j,nz+1) = 0.0
1663 do k=2,nz ;
do j=js,je ;
do i=is,ie
1664 tdif_flx(i,j,k) = (idt * 0.5*(ea_t(i,j,k) + eb_t(i,j,k-1))) * &
1665 (tv%T(i,j,k-1) - tv%T(i,j,k))
1666 tadv_flx(i,j,k) = (idt * (ea_t(i,j,k) - eb_t(i,j,k-1))) * &
1667 0.5*(tv%T(i,j,k-1) + tv%T(i,j,k))
1668 enddo ;
enddo ;
enddo 1670 if ((cs%id_Sdif > 0) .or. (cs%id_Sadv > 0))
then 1671 do j=js,je ;
do i=is,ie
1672 sdif_flx(i,j,1) = 0.0 ; sdif_flx(i,j,nz+1) = 0.0
1673 sadv_flx(i,j,1) = 0.0 ; sadv_flx(i,j,nz+1) = 0.0
1676 do k=2,nz ;
do j=js,je ;
do i=is,ie
1677 sdif_flx(i,j,k) = (idt * 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))) * &
1678 (tv%S(i,j,k-1) - tv%S(i,j,k))
1679 sadv_flx(i,j,k) = (idt * (ea_s(i,j,k) - eb_s(i,j,k-1))) * &
1680 0.5*(tv%S(i,j,k-1) + tv%S(i,j,k))
1681 enddo ;
enddo ;
enddo 1685 call cpu_clock_begin(id_clock_tracers)
1687 if (cs%mix_boundary_tracers)
then 1688 tr_ea_bbl = gv%Z_to_H * sqrt(dt*cs%Kd_BBL_tr)
1692 ebtr(i,j,nz) = eb_s(i,j,nz)
1694 in_boundary(i) = (g%mask2dT(i,j) > 0.0)
1696 do k=nz,2,-1 ;
do i=is,ie
1697 if (in_boundary(i))
then 1698 htot(i) = htot(i) + h(i,j,k)
1707 add_ent = ((dt * cs%Kd_min_tr) * gv%Z_to_H**2) * &
1708 ((h(i,j,k-1)+h(i,j,k)+h_neglect) / &
1709 (h(i,j,k-1)*h(i,j,k)+h_neglect2)) - &
1710 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))
1711 if (htot(i) < tr_ea_bbl)
then 1712 add_ent = max(0.0, add_ent, &
1713 (tr_ea_bbl - htot(i)) - min(ea_s(i,j,k),eb_s(i,j,k-1)))
1714 elseif (add_ent < 0.0)
then 1715 add_ent = 0.0 ; in_boundary(i) = .false.
1718 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1719 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1721 ebtr(i,j,k-1) = eb_s(i,j,k-1) ; eatr(i,j,k) = ea_s(i,j,k)
1724 if (
associated(visc%Kd_extra_S))
then ;
if (visc%Kd_extra_S(i,j,k) > 0.0)
then 1725 add_ent = ((dt * visc%Kd_extra_S(i,j,k)) * gv%Z_to_H**2) / &
1726 (0.5 * (h(i,j,k-1) + h(i,j,k)) + &
1728 ebtr(i,j,k-1) = ebtr(i,j,k-1) + add_ent
1729 eatr(i,j,k) = eatr(i,j,k) + add_ent
1732 do i=is,ie ; eatr(i,j,1) = ea_s(i,j,1) ;
enddo 1738 call call_tracer_column_fns(h_prebound, h, ea_s, eb_s, fluxes, hml, dt, g, gv, us, tv, &
1739 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1740 evap_cfl_limit = cs%evap_CFL_limit, &
1741 minimum_forcing_depth=cs%minimum_forcing_depth)
1743 elseif (
associated(visc%Kd_extra_S))
then 1745 do j=js,je ;
do i=is,ie
1746 ebtr(i,j,nz) = eb_s(i,j,nz) ; eatr(i,j,1) = ea_s(i,j,1)
1749 do k=nz,2,-1 ;
do j=js,je ;
do i=is,ie
1750 if (visc%Kd_extra_S(i,j,k) > 0.0)
then 1751 add_ent = ((dt * visc%Kd_extra_S(i,j,k)) * gv%Z_to_H**2) / &
1752 (0.5 * (h(i,j,k-1) + h(i,j,k)) + &
1757 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1758 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1759 enddo ;
enddo ;
enddo 1762 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1763 cs%optics, cs%tracer_flow_CSp, cs%debug,&
1764 evap_cfl_limit = cs%evap_CFL_limit, &
1765 minimum_forcing_depth=cs%minimum_forcing_depth)
1768 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1769 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1770 evap_cfl_limit = cs%evap_CFL_limit, &
1771 minimum_forcing_depth=cs%minimum_forcing_depth)
1774 call cpu_clock_end(id_clock_tracers)
1777 if (cs%use_sponge)
then 1778 call cpu_clock_begin(id_clock_sponge)
1779 if (
associated(cs%ALE_sponge_CSp))
then 1780 call apply_ale_sponge(h, dt, g, gv, us, cs%ALE_sponge_CSp, cs%Time)
1783 call cpu_clock_end(id_clock_sponge)
1786 call mom_thermovar_chksum(
"apply_sponge ", tv, g)
1790 call cpu_clock_begin(id_clock_pass)
1791 if (g%symmetric)
then ; dir_flag = to_all+omit_corners
1792 else ; dir_flag = to_west+to_south+omit_corners ;
endif 1797 call do_group_pass(cs%pass_hold_eb_ea, g%Domain)
1799 if (
associated(visc%Kv_slow)) &
1800 call pass_var(visc%Kv_slow, g%Domain, to_all+omit_corners, halo=1)
1801 call cpu_clock_end(id_clock_pass)
1803 call disable_averaging(cs%diag)
1806 call enable_averages(dt, time_end, cs%diag)
1808 if (cs%id_Kd_interface > 0)
call post_data(cs%id_Kd_interface, kd_int, cs%diag)
1809 if (cs%id_Kd_heat > 0)
call post_data(cs%id_Kd_heat, kd_heat, cs%diag)
1810 if (cs%id_Kd_salt > 0)
call post_data(cs%id_Kd_salt, kd_salt, cs%diag)
1811 if (cs%id_Kd_ePBL > 0)
call post_data(cs%id_Kd_ePBL, kd_epbl, cs%diag)
1813 if (cs%id_ea_t > 0)
call post_data(cs%id_ea_t, ea_t, cs%diag)
1814 if (cs%id_eb_t > 0)
call post_data(cs%id_eb_t, eb_t, cs%diag)
1815 if (cs%id_ea_s > 0)
call post_data(cs%id_ea_s, ea_s, cs%diag)
1816 if (cs%id_eb_s > 0)
call post_data(cs%id_eb_s, eb_s, cs%diag)
1818 if (cs%id_dudt_dia > 0)
call post_data(cs%id_dudt_dia, adp%du_dt_dia, cs%diag)
1819 if (cs%id_dvdt_dia > 0)
call post_data(cs%id_dvdt_dia, adp%dv_dt_dia, cs%diag)
1821 if (cs%id_Tdif > 0)
call post_data(cs%id_Tdif, tdif_flx, cs%diag)
1822 if (cs%id_Tadv > 0)
call post_data(cs%id_Tadv, tadv_flx, cs%diag)
1823 if (cs%id_Sdif > 0)
call post_data(cs%id_Sdif, sdif_flx, cs%diag)
1824 if (cs%id_Sadv > 0)
call post_data(cs%id_Sadv, sadv_flx, cs%diag)
1827 if (cs%id_hf_dudt_dia_2d > 0)
then 1828 allocate(hf_dudt_dia_2d(g%IsdB:g%IedB,g%jsd:g%jed))
1829 hf_dudt_dia_2d(:,:) = 0.0
1830 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1831 hf_dudt_dia_2d(i,j) = hf_dudt_dia_2d(i,j) + adp%du_dt_dia(i,j,k) * adp%diag_hfrac_u(i,j,k)
1832 enddo ;
enddo ;
enddo 1833 call post_data(cs%id_hf_dudt_dia_2d, hf_dudt_dia_2d, cs%diag)
1834 deallocate(hf_dudt_dia_2d)
1837 if (cs%id_hf_dvdt_dia_2d > 0)
then 1838 allocate(hf_dvdt_dia_2d(g%isd:g%ied,g%JsdB:g%JedB))
1839 hf_dvdt_dia_2d(:,:) = 0.0
1840 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1841 hf_dvdt_dia_2d(i,j) = hf_dvdt_dia_2d(i,j) + adp%dv_dt_dia(i,j,k) * adp%diag_hfrac_v(i,j,k)
1842 enddo ;
enddo ;
enddo 1843 call post_data(cs%id_hf_dvdt_dia_2d, hf_dvdt_dia_2d, cs%diag)
1844 deallocate(hf_dvdt_dia_2d)
1847 call disable_averaging(cs%diag)
1849 if (showcalltree)
call calltree_leave(
"diabatic_ALE()")
1851 end subroutine diabatic_ale
1855 subroutine layered_diabatic(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
1856 G, GV, US, CS, WAVES)
1860 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
1861 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
1862 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
1865 real,
dimension(:,:),
pointer :: Hml
1866 type(
forcing),
intent(inout) :: fluxes
1873 real,
intent(in) :: dt
1874 type(time_type),
intent(in) :: Time_end
1878 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
1893 real,
dimension(SZI_(G),SZJ_(G)) :: &
1896 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: h_diag
1897 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
1898 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: saln_diag
1899 real,
dimension(SZI_(G),SZJ_(G)) :: tendency_2d
1903 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
target :: &
1909 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
target :: &
1919 real,
allocatable,
dimension(:,:) :: &
1920 hf_dudt_dia_2d, hf_dvdt_dia_2d
1923 real,
pointer,
dimension(:,:,:) :: &
1929 integer :: kb(szi_(g),szj_(g))
1932 real :: p_ref_cv(szi_(g))
1935 logical :: in_boundary(szi_(g))
1955 real :: htot(szib_(g))
1956 real :: b1(szib_(g)), d1(szib_(g))
1957 real :: c1(szib_(g),szk_(g))
1964 logical :: showCallTree
1965 integer,
dimension(2) :: EOSdom
1966 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, nkmb, m, halo
1968 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1969 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1970 nkmb = gv%nk_rho_varies
1971 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect*h_neglect
1972 kd_heat(:,:,:) = 0.0 ; kd_salt(:,:,:) = 0.0
1975 showcalltree = calltree_showquery()
1976 if (showcalltree)
call calltree_enter(
"layered_diabatic(), MOM_diabatic_driver.F90")
1979 eaml => eatr ; ebml => ebtr
1982 call enable_averages(dt, time_end, cs%diag)
1984 if ((cs%ML_mix_first > 0.0) .or. cs%use_geothermal)
then 1985 halo = cs%halo_TS_diff
1987 do k=1,nz ;
do j=js-halo,je+halo ;
do i=is-halo,ie+halo
1988 h_orig(i,j,k) = h(i,j,k) ; eaml(i,j,k) = 0.0 ; ebml(i,j,k) = 0.0
1989 enddo ;
enddo ;
enddo 1992 if (cs%use_geothermal)
then 1993 call cpu_clock_begin(id_clock_geothermal)
1994 call geothermal_entraining(h, tv, dt, eaml, ebml, g, gv, us, cs%geothermal_CSp, halo=cs%halo_TS_diff)
1995 call cpu_clock_end(id_clock_geothermal)
1996 if (showcalltree)
call calltree_waypoint(
"geothermal (diabatic)")
1997 if (cs%debugConservation)
call mom_state_stats(
'geothermal', u, v, h, tv%T, tv%S, g, gv, us)
2002 call diag_update_remap_grids(cs%diag)
2007 if (
associated(cs%optics)) &
2008 call set_pen_shortwave(cs%optics, fluxes, g, gv, us, cs%diabatic_aux_CSp, cs%opacity_CSp, cs%tracer_flow_CSp)
2010 if (cs%bulkmixedlayer)
then 2011 if (cs%debug)
call mom_forcing_chksum(
"Before mixedlayer", fluxes, g, us, haloshift=0)
2013 if (cs%ML_mix_first > 0.0)
then 2021 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
2023 call cpu_clock_begin(id_clock_mixedlayer)
2024 if (cs%ML_mix_first < 1.0)
then 2026 call bulkmixedlayer(h, u_h, v_h, tv, fluxes, dt*cs%ML_mix_first, &
2027 eaml,ebml, g, gv, us, cs%bulkmixedlayer_CSp, cs%optics, &
2028 hml, cs%aggregate_FW_forcing, dt, last_call=.false.)
2030 call bulkmixedlayer(h, u_h, v_h, tv, fluxes, dt, eaml, ebml, &
2031 g, gv, us, cs%bulkmixedlayer_CSp, cs%optics, &
2032 hml, cs%aggregate_FW_forcing, dt, last_call=.true.)
2041 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
2042 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
2043 call cpu_clock_end(id_clock_mixedlayer)
2045 call mom_state_chksum(
"After mixedlayer ", u, v, h, g, gv, us, haloshift=0)
2046 call mom_forcing_chksum(
"After mixedlayer", fluxes, g, us, haloshift=0)
2048 if (showcalltree)
call calltree_waypoint(
"done with 1st bulkmixedlayer (diabatic)")
2049 if (cs%debugConservation)
call mom_state_stats(
'1st bulkmixedlayer', u, v, h, tv%T, tv%S, g, gv, us)
2054 call mom_state_chksum(
"before find_uv_at_h", u, v, h, g, gv, us, haloshift=0)
2055 if (cs%use_kappa_shear .or. cs%use_CVMix_shear)
then 2056 if ((cs%ML_mix_first > 0.0) .or. cs%use_geothermal)
then 2057 call find_uv_at_h(u, v, h_orig, u_h, v_h, g, gv, us, eaml, ebml)
2059 call hchksum(eaml,
"after find_uv_at_h eaml", g%HI, scale=gv%H_to_m)
2060 call hchksum(ebml,
"after find_uv_at_h ebml", g%HI, scale=gv%H_to_m)
2063 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
2065 if (showcalltree)
call calltree_waypoint(
"done with find_uv_at_h (diabatic)")
2068 call cpu_clock_begin(id_clock_set_diffusivity)
2071 if ((cs%halo_TS_diff > 0) .and. (cs%ML_mix_first > 0.0))
then 2072 if (
associated(tv%T))
call pass_var(tv%T, g%Domain, halo=cs%halo_TS_diff, complete=.false.)
2073 if (
associated(tv%T))
call pass_var(tv%S, g%Domain, halo=cs%halo_TS_diff, complete=.false.)
2074 call pass_var(h, g%domain, halo=cs%halo_TS_diff, complete=.true.)
2077 call mom_state_chksum(
"before set_diffusivity", u, v, h, g, gv, us, haloshift=cs%halo_TS_diff)
2078 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, &
2079 visc, dt, g, gv, us, cs%set_diff_CSp, kd_lay, kd_int)
2080 call cpu_clock_end(id_clock_set_diffusivity)
2081 if (showcalltree)
call calltree_waypoint(
"done with set_diffusivity (diabatic)")
2084 call mom_state_chksum(
"after set_diffusivity ", u, v, h, g, gv, us, haloshift=0)
2085 call mom_forcing_chksum(
"after set_diffusivity ", fluxes, g, us, haloshift=0)
2086 call mom_thermovar_chksum(
"after set_diffusivity ", tv, g)
2087 call hchksum(kd_lay,
"after set_diffusivity Kd_lay", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2088 call hchksum(kd_int,
"after set_diffusivity Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2093 call cpu_clock_begin(id_clock_kpp)
2099 call calculatebuoyancyflux2d(g, gv, us, fluxes, cs%optics, h, tv%T, tv%S, tv, &
2100 cs%KPP_buoy_flux, cs%KPP_temp_flux, cs%KPP_salt_flux)
2106 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2107 kd_salt(i,j,k) = kd_int(i,j,k)
2108 kd_heat(i,j,k) = kd_int(i,j,k)
2109 enddo ;
enddo ;
enddo 2111 if (
associated(visc%Kd_extra_S))
then 2113 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2114 kd_salt(i,j,k) = kd_salt(i,j,k) + visc%Kd_extra_S(i,j,k)
2115 enddo ;
enddo ;
enddo 2117 if (
associated(visc%Kd_extra_T))
then 2119 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2120 kd_heat(i,j,k) = kd_heat(i,j,k) + visc%Kd_extra_T(i,j,k)
2121 enddo ;
enddo ;
enddo 2124 call kpp_compute_bld(cs%KPP_CSp, g, gv, us, h, tv%T, tv%S, u, v, tv, &
2125 fluxes%ustar, cs%KPP_buoy_flux, waves=waves)
2127 call kpp_calculate(cs%KPP_CSp, g, gv, us, h, fluxes%ustar, cs%KPP_buoy_flux, kd_heat, &
2128 kd_salt, visc%Kv_shear, cs%KPP_NLTheat, cs%KPP_NLTscalar, waves=waves)
2130 if (
associated(hml))
then 2131 call kpp_get_bld(cs%KPP_CSp, hml(:,:), g, us)
2132 call pass_var(hml, g%domain, halo=1)
2134 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
2137 if (.not. cs%KPPisPassive)
then 2139 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2140 kd_int(i,j,k) = min( kd_salt(i,j,k), kd_heat(i,j,k) )
2141 enddo ;
enddo ;
enddo 2142 if (
associated(visc%Kd_extra_S))
then 2144 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2145 visc%Kd_extra_S(i,j,k) = (kd_salt(i,j,k) - kd_int(i,j,k))
2146 enddo ;
enddo ;
enddo 2148 if (
associated(visc%Kd_extra_T))
then 2150 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2151 visc%Kd_extra_T(i,j,k) = (kd_heat(i,j,k) - kd_int(i,j,k))
2152 enddo ;
enddo ;
enddo 2156 call cpu_clock_end(id_clock_kpp)
2157 if (showcalltree)
call calltree_waypoint(
"done with KPP_calculate (diabatic)")
2160 call mom_forcing_chksum(
"after KPP", fluxes, g, us, haloshift=0)
2161 call mom_thermovar_chksum(
"after KPP", tv, g)
2162 call hchksum(kd_lay,
"after KPP Kd_lay", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2163 call hchksum(kd_int,
"after KPP Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2169 if (cs%use_CVMix_conv)
then 2170 call calculate_cvmix_conv(h, tv, g, gv, us, cs%CVMix_conv_csp, hml)
2172 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2173 kd_int(i,j,k) = kd_int(i,j,k) + cs%CVMix_conv_csp%kd_conv(i,j,k)
2174 visc%Kv_slow(i,j,k) = visc%Kv_slow(i,j,k) + cs%CVMix_conv_csp%kv_conv(i,j,k)
2175 enddo ;
enddo ;
enddo 2179 call cpu_clock_begin(id_clock_kpp)
2181 call hchksum(cs%KPP_temp_flux,
"before KPP_applyNLT netHeat", g%HI, haloshift=0, scale=gv%H_to_m)
2182 call hchksum(cs%KPP_salt_flux,
"before KPP_applyNLT netSalt", g%HI, haloshift=0, scale=gv%H_to_m)
2183 call hchksum(cs%KPP_NLTheat,
"before KPP_applyNLT NLTheat", g%HI, haloshift=0)
2184 call hchksum(cs%KPP_NLTscalar,
"before KPP_applyNLT NLTscalar", g%HI, haloshift=0)
2188 call kpp_nonlocaltransport_temp(cs%KPP_CSp, g, gv, h, cs%KPP_NLTheat, cs%KPP_temp_flux, &
2189 us%T_to_s*dt, tv%T, us%Q_to_J_kg*tv%C_p)
2190 call kpp_nonlocaltransport_saln(cs%KPP_CSp, g, gv, h, cs%KPP_NLTscalar, cs%KPP_salt_flux, &
2192 call cpu_clock_end(id_clock_kpp)
2193 if (showcalltree)
call calltree_waypoint(
"done with KPP_applyNonLocalTransport (diabatic)")
2194 if (cs%debugConservation)
call mom_state_stats(
'KPP_applyNonLocalTransport', u, v, h, tv%T, tv%S, g, gv, us)
2197 call mom_state_chksum(
"after KPP_applyNLT ", u, v, h, g, gv, us, haloshift=0)
2198 call mom_forcing_chksum(
"after KPP_applyNLT ", fluxes, g, us, haloshift=0)
2199 call mom_thermovar_chksum(
"after KPP_applyNLT ", tv, g)
2205 if (
associated(visc%Kd_extra_T) .and.
associated(visc%Kd_extra_S) .and.
associated(tv%T))
then 2207 call cpu_clock_begin(id_clock_differential_diff)
2208 call differential_diffuse_t_s(h, tv, visc, dt, g, gv)
2209 call cpu_clock_end(id_clock_differential_diff)
2210 if (showcalltree)
call calltree_waypoint(
"done with differential_diffuse_T_S (diabatic)")
2211 if (cs%debugConservation)
call mom_state_stats(
'differential_diffuse_T_S', u, v, h, tv%T, tv%S, g, gv, us)
2215 if (.not. cs%useKPP)
then 2217 do k=2,nz ;
do j=js,je ;
do i=is,ie
2218 kd_heat(i,j,k) = kd_heat(i,j,k) + visc%Kd_extra_T(i,j,k)
2219 kd_salt(i,j,k) = kd_salt(i,j,k) + visc%Kd_extra_S(i,j,k)
2220 enddo ;
enddo ;
enddo 2227 call cpu_clock_begin(id_clock_entrain)
2230 call entrainment_diffusive(h, tv, fluxes, dt, g, gv, us, cs%entrain_diffusive_CSp, &
2231 ea, eb, kb, kd_lay=kd_lay, kd_int=kd_int)
2232 call cpu_clock_end(id_clock_entrain)
2233 if (showcalltree)
call calltree_waypoint(
"done with Entrainment_diffusive (diabatic)")
2236 call mom_forcing_chksum(
"after calc_entrain ", fluxes, g, us, haloshift=0)
2237 call mom_thermovar_chksum(
"after calc_entrain ", tv, g)
2238 call mom_state_chksum(
"after calc_entrain ", u, v, h, g, gv, us, haloshift=0)
2239 call hchksum(ea,
"after calc_entrain ea", g%HI, haloshift=0, scale=gv%H_to_m)
2240 call hchksum(eb,
"after calc_entrain eb", g%HI, haloshift=0, scale=gv%H_to_m)
2244 if (cs%boundary_forcing_tendency_diag)
then 2245 do k=1,nz ;
do j=js,je ;
do i=is,ie
2246 h_diag(i,j,k) = h(i,j,k)
2247 temp_diag(i,j,k) = tv%T(i,j,k)
2248 saln_diag(i,j,k) = tv%S(i,j,k)
2249 enddo ;
enddo ;
enddo 2263 hold(i,j,1) = h(i,j,1)
2264 h(i,j,1) = h(i,j,1) + (eb(i,j,1) - ea(i,j,2))
2265 hold(i,j,nz) = h(i,j,nz)
2266 h(i,j,nz) = h(i,j,nz) + (ea(i,j,nz) - eb(i,j,nz-1))
2267 if (h(i,j,1) <= 0.0)
then 2268 h(i,j,1) = gv%Angstrom_H
2270 if (h(i,j,nz) <= 0.0)
then 2271 h(i,j,nz) = gv%Angstrom_H
2274 do k=2,nz-1 ;
do i=is,ie
2275 hold(i,j,k) = h(i,j,k)
2276 h(i,j,k) = h(i,j,k) + ((ea(i,j,k) - eb(i,j,k-1)) + &
2277 (eb(i,j,k) - ea(i,j,k+1)))
2278 if (h(i,j,k) <= 0.0)
then 2279 h(i,j,k) = gv%Angstrom_H
2284 call diag_update_remap_grids(cs%diag)
2287 call mom_state_chksum(
"after negative check ", u, v, h, g, gv, us, haloshift=0)
2288 call mom_forcing_chksum(
"after negative check ", fluxes, g, us, haloshift=0)
2289 call mom_thermovar_chksum(
"after negative check ", tv, g)
2291 if (showcalltree)
call calltree_waypoint(
"done with h=ea-eb (diabatic)")
2292 if (cs%debugConservation)
call mom_state_stats(
'h=ea-eb', u, v, h, tv%T, tv%S, g, gv, us)
2298 if (cs%bulkmixedlayer)
then 2300 if (
associated(tv%T))
then 2301 call cpu_clock_begin(id_clock_tridiag)
2306 if (cs%massless_match_targets)
then 2310 h_tr = hold(i,j,1) + h_neglect
2311 b1(i) = 1.0 / (h_tr + eb(i,j,1))
2312 d1(i) = h_tr * b1(i)
2313 tv%T(i,j,1) = b1(i) * (h_tr*tv%T(i,j,1))
2314 tv%S(i,j,1) = b1(i) * (h_tr*tv%S(i,j,1))
2316 do k=2,nkmb ;
do i=is,ie
2317 c1(i,k) = eb(i,j,k-1) * b1(i)
2318 h_tr = hold(i,j,k) + h_neglect
2319 b_denom_1 = h_tr + d1(i)*ea(i,j,k)
2320 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
2321 if (k<nkmb) d1(i) = b_denom_1 * b1(i)
2322 tv%T(i,j,k) = b1(i) * (h_tr*tv%T(i,j,k) + ea(i,j,k)*tv%T(i,j,k-1))
2323 tv%S(i,j,k) = b1(i) * (h_tr*tv%S(i,j,k) + ea(i,j,k)*tv%S(i,j,k-1))
2326 do k=nkmb+1,nz ;
do i=is,ie
2327 if (k == kb(i,j))
then 2328 c1(i,k) = eb(i,j,k-1) * b1(i)
2329 d1(i) = (((eb(i,j,nkmb)-eb(i,j,k-1)) + hold(i,j,nkmb) + h_neglect) + &
2330 d1(i)*ea(i,j,nkmb)) * b1(i)
2331 h_tr = hold(i,j,k) + h_neglect
2332 b_denom_1 = h_tr + d1(i)*ea(i,j,k)
2333 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
2334 d1(i) = b_denom_1 * b1(i)
2335 tv%T(i,j,k) = b1(i) * (h_tr*tv%T(i,j,k) + ea(i,j,k)*tv%T(i,j,nkmb))
2336 tv%S(i,j,k) = b1(i) * (h_tr*tv%S(i,j,k) + ea(i,j,k)*tv%S(i,j,nkmb))
2337 elseif (k > kb(i,j))
then 2338 c1(i,k) = eb(i,j,k-1) * b1(i)
2339 h_tr = hold(i,j,k) + h_neglect
2340 b_denom_1 = h_tr + d1(i)*ea(i,j,k)
2341 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
2342 d1(i) = b_denom_1 * b1(i)
2343 tv%T(i,j,k) = b1(i) * (h_tr*tv%T(i,j,k) + ea(i,j,k)*tv%T(i,j,k-1))
2344 tv%S(i,j,k) = b1(i) * (h_tr*tv%S(i,j,k) + ea(i,j,k)*tv%S(i,j,k-1))
2345 elseif (eb(i,j,k) < eb(i,j,k-1))
then 2350 tv%T(i,j,nkmb) = tv%T(i,j,nkmb) + b1(i) * (eb(i,j,k-1) - eb(i,j,k)) * tv%T(i,j,k)
2351 tv%S(i,j,nkmb) = tv%S(i,j,nkmb) + b1(i) * (eb(i,j,k-1) - eb(i,j,k)) * tv%S(i,j,k)
2355 do k=nz-1,nkmb,-1 ;
do i=is,ie
2356 if (k >= kb(i,j))
then 2357 tv%T(i,j,k) = tv%T(i,j,k) + c1(i,k+1)*tv%T(i,j,k+1)
2358 tv%S(i,j,k) = tv%S(i,j,k) + c1(i,k+1)*tv%S(i,j,k+1)
2361 do i=is,ie ;
if (kb(i,j) <= nz)
then 2362 tv%T(i,j,nkmb) = tv%T(i,j,nkmb) + c1(i,kb(i,j))*tv%T(i,j,kb(i,j))
2363 tv%S(i,j,nkmb) = tv%S(i,j,nkmb) + c1(i,kb(i,j))*tv%S(i,j,kb(i,j))
2365 do k=nkmb-1,1,-1 ;
do i=is,ie
2366 tv%T(i,j,k) = tv%T(i,j,k) + c1(i,k+1)*tv%T(i,j,k+1)
2367 tv%S(i,j,k) = tv%S(i,j,k) + c1(i,k+1)*tv%S(i,j,k+1)
2374 if (cs%tracer_tridiag)
then 2375 call tracer_vertdiff(hold, ea, eb, dt, tv%T, g, gv)
2376 call tracer_vertdiff(hold, ea, eb, dt, tv%S, g, gv)
2378 call tridiagts(g, gv, is, ie, js, je, hold, ea, eb, tv%T, tv%S)
2381 call cpu_clock_end(id_clock_tridiag)
2384 if (cs%debugConservation)
call mom_state_stats(
'BML tridiag', u, v, h, tv%T, tv%S, g, gv, us)
2386 if ((cs%ML_mix_first > 0.0) .or. cs%use_geothermal)
then 2390 do k=1,nz ;
do j=js,je ;
do i=is,ie
2391 hold(i,j,k) = h_orig(i,j,k)
2392 ea(i,j,k) = ea(i,j,k) + eaml(i,j,k)
2393 eb(i,j,k) = eb(i,j,k) + ebml(i,j,k)
2394 enddo ;
enddo ;
enddo 2396 call hchksum(ea,
"after ea = ea + eaml", g%HI, haloshift=0, scale=gv%H_to_m)
2397 call hchksum(eb,
"after eb = eb + ebml", g%HI, haloshift=0, scale=gv%H_to_m)
2401 if (cs%ML_mix_first < 1.0)
then 2410 call find_uv_at_h(u, v, hold, u_h, v_h, g, gv, us, ea, eb)
2411 if (cs%debug)
call mom_state_chksum(
"find_uv_at_h1 ", u, v, h, g, gv, us, haloshift=0)
2413 dt_mix = min(dt, dt*(1.0 - cs%ML_mix_first))
2414 call cpu_clock_begin(id_clock_mixedlayer)
2416 call bulkmixedlayer(h, u_h, v_h, tv, fluxes, dt_mix, ea, eb, &
2417 g, gv, us, cs%bulkmixedlayer_CSp, cs%optics, &
2418 hml, cs%aggregate_FW_forcing, dt, last_call=.true.)
2426 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
2427 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
2429 call cpu_clock_end(id_clock_mixedlayer)
2430 if (showcalltree)
call calltree_waypoint(
"done with 2nd bulkmixedlayer (diabatic)")
2431 if (cs%debugConservation)
call mom_state_stats(
'2nd bulkmixedlayer', u, v, h, tv%T, tv%S, g, gv, us)
2437 if (
associated(tv%T))
then 2440 call hchksum(ea,
"before triDiagTS ea ", g%HI, haloshift=0, scale=gv%H_to_m)
2441 call hchksum(eb,
"before triDiagTS eb ", g%HI, haloshift=0, scale=gv%H_to_m)
2443 call cpu_clock_begin(id_clock_tridiag)
2451 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
2452 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
2454 if (cs%diabatic_diff_tendency_diag)
then 2455 do k=1,nz ;
do j=js,je ;
do i=is,ie
2456 temp_diag(i,j,k) = tv%T(i,j,k)
2457 saln_diag(i,j,k) = tv%S(i,j,k)
2458 enddo ;
enddo ;
enddo 2462 if (cs%tracer_tridiag)
then 2463 call tracer_vertdiff(hold, ea, eb, dt, tv%T, g, gv)
2464 call tracer_vertdiff(hold, ea, eb, dt, tv%S, g, gv)
2466 call tridiagts(g, gv, is, ie, js, je, hold, ea, eb, tv%T, tv%S)
2472 if (cs%diabatic_diff_tendency_diag)
then 2473 call diagnose_diabatic_diff_tendency(tv, hold, temp_diag, saln_diag, dt, g, gv, us, cs)
2474 if (cs%id_diabatic_diff_h > 0)
call post_data(cs%id_diabatic_diff_h, hold, cs%diag, alt_h=hold)
2477 call cpu_clock_end(id_clock_tridiag)
2478 if (showcalltree)
call calltree_waypoint(
"done with triDiagTS (diabatic)")
2481 if (cs%debugConservation)
call mom_state_stats(
'triDiagTS', u, v, h, tv%T, tv%S, g, gv, us)
2486 call mom_state_chksum(
"after mixed layer ", u, v, h, g, gv, us, haloshift=0)
2487 call mom_thermovar_chksum(
"after mixed layer ", tv, g)
2488 call hchksum(ea,
"after mixed layer ea", g%HI, scale=gv%H_to_m)
2489 call hchksum(eb,
"after mixed layer eb", g%HI, scale=gv%H_to_m)
2492 call cpu_clock_begin(id_clock_remap)
2493 call regularize_layers(h, tv, dt, ea, eb, g, gv, us, cs%regularize_layers_CSp)
2494 call cpu_clock_end(id_clock_remap)
2495 if (showcalltree)
call calltree_waypoint(
"done with regularize_layers (diabatic)")
2496 if (cs%debugConservation)
call mom_state_stats(
'regularize_layers', u, v, h, tv%T, tv%S, g, gv, us)
2500 if (
associated(adp%du_dt_dia) .or.
associated(adp%dv_dt_dia)) &
2502 call pass_var(h, g%domain, to_west+to_south+omit_corners, halo=1)
2503 call diag_update_remap_grids(cs%diag)
2507 if ((cs%id_Tdif > 0) .or. (cs%id_Tadv > 0))
then 2508 do j=js,je ;
do i=is,ie
2509 tdif_flx(i,j,1) = 0.0 ; tdif_flx(i,j,nz+1) = 0.0
2510 tadv_flx(i,j,1) = 0.0 ; tadv_flx(i,j,nz+1) = 0.0
2513 do k=2,nz ;
do j=js,je ;
do i=is,ie
2514 tdif_flx(i,j,k) = (idt * 0.5*(ea(i,j,k) + eb(i,j,k-1))) * &
2515 (tv%T(i,j,k-1) - tv%T(i,j,k))
2516 tadv_flx(i,j,k) = (idt * (ea(i,j,k) - eb(i,j,k-1))) * &
2517 0.5*(tv%T(i,j,k-1) + tv%T(i,j,k))
2518 enddo ;
enddo ;
enddo 2520 if ((cs%id_Sdif > 0) .or. (cs%id_Sadv > 0))
then 2521 do j=js,je ;
do i=is,ie
2522 sdif_flx(i,j,1) = 0.0 ; sdif_flx(i,j,nz+1) = 0.0
2523 sadv_flx(i,j,1) = 0.0 ; sadv_flx(i,j,nz+1) = 0.0
2526 do k=2,nz ;
do j=js,je ;
do i=is,ie
2527 sdif_flx(i,j,k) = (idt * 0.5*(ea(i,j,k) + eb(i,j,k-1))) * &
2528 (tv%S(i,j,k-1) - tv%S(i,j,k))
2529 sadv_flx(i,j,k) = (idt * (ea(i,j,k) - eb(i,j,k-1))) * &
2530 0.5*(tv%S(i,j,k-1) + tv%S(i,j,k))
2531 enddo ;
enddo ;
enddo 2535 call cpu_clock_begin(id_clock_tracers)
2536 if (cs%mix_boundary_tracers)
then 2537 tr_ea_bbl = gv%Z_to_H * sqrt(dt*cs%Kd_BBL_tr)
2541 ebtr(i,j,nz) = eb(i,j,nz)
2543 in_boundary(i) = (g%mask2dT(i,j) > 0.0)
2545 do k=nz,2,-1 ;
do i=is,ie
2546 if (in_boundary(i))
then 2547 htot(i) = htot(i) + h(i,j,k)
2556 add_ent = ((dt * cs%Kd_min_tr) * gv%Z_to_H**2) * &
2557 ((h(i,j,k-1)+h(i,j,k)+h_neglect) / &
2558 (h(i,j,k-1)*h(i,j,k)+h_neglect2)) - &
2559 0.5*(ea(i,j,k) + eb(i,j,k-1))
2560 if (htot(i) < tr_ea_bbl)
then 2561 add_ent = max(0.0, add_ent, &
2562 (tr_ea_bbl - htot(i)) - min(ea(i,j,k),eb(i,j,k-1)))
2563 elseif (add_ent < 0.0)
then 2564 add_ent = 0.0 ; in_boundary(i) = .false.
2567 ebtr(i,j,k-1) = eb(i,j,k-1) + add_ent
2568 eatr(i,j,k) = ea(i,j,k) + add_ent
2570 ebtr(i,j,k-1) = eb(i,j,k-1) ; eatr(i,j,k) = ea(i,j,k)
2572 if (
associated(visc%Kd_extra_S))
then ;
if (visc%Kd_extra_S(i,j,k) > 0.0)
then 2573 add_ent = ((dt * visc%Kd_extra_S(i,j,k)) * gv%Z_to_H**2) / &
2574 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + &
2576 ebtr(i,j,k-1) = ebtr(i,j,k-1) + add_ent
2577 eatr(i,j,k) = eatr(i,j,k) + add_ent
2580 do i=is,ie ; eatr(i,j,1) = ea(i,j,1) ;
enddo 2584 call call_tracer_column_fns(hold, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
2585 cs%optics, cs%tracer_flow_CSp, cs%debug)
2587 elseif (
associated(visc%Kd_extra_S))
then 2589 do j=js,je ;
do i=is,ie
2590 ebtr(i,j,nz) = eb(i,j,nz) ; eatr(i,j,1) = ea(i,j,1)
2593 do k=nz,2,-1 ;
do j=js,je ;
do i=is,ie
2594 if (visc%Kd_extra_S(i,j,k) > 0.0)
then 2595 add_ent = ((dt * visc%Kd_extra_S(i,j,k)) * gv%Z_to_H**2) / &
2596 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + &
2601 ebtr(i,j,k-1) = eb(i,j,k-1) + add_ent
2602 eatr(i,j,k) = ea(i,j,k) + add_ent
2603 enddo ;
enddo ;
enddo 2605 call call_tracer_column_fns(hold, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
2606 cs%optics, cs%tracer_flow_CSp, cs%debug)
2609 call call_tracer_column_fns(hold, h, ea, eb, fluxes, hml, dt, g, gv, us, tv, &
2610 cs%optics, cs%tracer_flow_CSp, cs%debug)
2614 call cpu_clock_end(id_clock_tracers)
2617 if (cs%use_sponge)
then 2618 call cpu_clock_begin(id_clock_sponge)
2620 if (cs%bulkmixedlayer .and.
associated(tv%eqn_of_state))
then 2621 do i=is,ie ; p_ref_cv(i) = tv%P_Ref ;
enddo 2622 eosdom(:) = eos_domain(g%HI)
2626 tv%eqn_of_state, eosdom)
2628 call apply_sponge(h, dt, g, gv, us, ea, eb, cs%sponge_CSp, rcv_ml)
2630 call apply_sponge(h, dt, g, gv, us, ea, eb, cs%sponge_CSp)
2632 call cpu_clock_end(id_clock_sponge)
2635 call mom_thermovar_chksum(
"apply_sponge ", tv, g)
2640 if (
associated(cdp%diapyc_vel))
then 2643 do k=2,nz ;
do i=is,ie
2644 cdp%diapyc_vel(i,j,k) = us%s_to_T*idt * (ea(i,j,k) - eb(i,j,k-1))
2647 cdp%diapyc_vel(i,j,1) = 0.0
2648 cdp%diapyc_vel(i,j,nz+1) = 0.0
2656 if (cs%bulkmixedlayer)
then 2658 call hchksum(ea,
"before net flux rearrangement ea", g%HI, scale=gv%H_to_m)
2659 call hchksum(eb,
"before net flux rearrangement eb", g%HI, scale=gv%H_to_m)
2663 do k=2,gv%nkml ;
do i=is,ie
2664 net_ent = ea(i,j,k) - eb(i,j,k-1)
2665 ea(i,j,k) = max(net_ent, 0.0)
2666 eb(i,j,k-1) = max(-net_ent, 0.0)
2670 call hchksum(ea,
"after net flux rearrangement ea", g%HI, scale=gv%H_to_m)
2671 call hchksum(eb,
"after net flux rearrangement eb", g%HI, scale=gv%H_to_m)
2679 hold(i,js-1,k) = gv%Angstrom_H ; ea(i,js-1,k) = 0.0 ; eb(i,js-1,k) = 0.0
2680 hold(i,je+1,k) = gv%Angstrom_H ; ea(i,je+1,k) = 0.0 ; eb(i,je+1,k) = 0.0
2683 hold(is-1,j,k) = gv%Angstrom_H ; ea(is-1,j,k) = 0.0 ; eb(is-1,j,k) = 0.0
2684 hold(ie+1,j,k) = gv%Angstrom_H ; ea(ie+1,j,k) = 0.0 ; eb(ie+1,j,k) = 0.0
2688 call cpu_clock_begin(id_clock_pass)
2689 if (g%symmetric)
then ; dir_flag = to_all+omit_corners
2690 else ; dir_flag = to_west+to_south+omit_corners ;
endif 2694 call do_group_pass(cs%pass_hold_eb_ea, g%Domain)
2695 call cpu_clock_end(id_clock_pass)
2701 call mom_state_chksum(
"before u/v tridiag ", u, v, h, g, gv, us, haloshift=0)
2702 call hchksum(ea,
"before u/v tridiag ea", g%HI, scale=gv%H_to_m)
2703 call hchksum(eb,
"before u/v tridiag eb", g%HI, scale=gv%H_to_m)
2704 call hchksum(hold,
"before u/v tridiag hold", g%HI, scale=gv%H_to_m)
2706 call cpu_clock_begin(id_clock_tridiag)
2707 idt_accel = 1.0 / dt
2711 if (
associated(adp%du_dt_dia)) adp%du_dt_dia(i,j,1) = u(i,j,1)
2712 hval = (hold(i,j,1) + hold(i+1,j,1)) + (ea(i,j,1) + ea(i+1,j,1)) + h_neglect
2713 b1(i) = 1.0 / (hval + (eb(i,j,1) + eb(i+1,j,1)))
2714 d1(i) = hval * b1(i)
2715 u(i,j,1) = b1(i) * (hval * u(i,j,1))
2717 do k=2,nz ;
do i=isq,ieq
2718 if (
associated(adp%du_dt_dia)) adp%du_dt_dia(i,j,k) = u(i,j,k)
2719 c1(i,k) = (eb(i,j,k-1)+eb(i+1,j,k-1)) * b1(i)
2720 eaval = ea(i,j,k) + ea(i+1,j,k)
2721 hval = hold(i,j,k) + hold(i+1,j,k) + h_neglect
2722 b1(i) = 1.0 / ((eb(i,j,k) + eb(i+1,j,k)) + (hval + d1(i)*eaval))
2723 d1(i) = (hval + d1(i)*eaval) * b1(i)
2724 u(i,j,k) = (hval*u(i,j,k) + eaval*u(i,j,k-1))*b1(i)
2726 do k=nz-1,1,-1 ;
do i=isq,ieq
2727 u(i,j,k) = u(i,j,k) + c1(i,k+1)*u(i,j,k+1)
2728 if (
associated(adp%du_dt_dia)) &
2729 adp%du_dt_dia(i,j,k) = (u(i,j,k) - adp%du_dt_dia(i,j,k)) * idt_accel
2731 if (
associated(adp%du_dt_dia))
then 2733 adp%du_dt_dia(i,j,nz) = (u(i,j,nz)-adp%du_dt_dia(i,j,nz)) * idt_accel
2738 call mom_state_chksum(
"aft 1st loop tridiag ", u, v, h, g, gv, us, haloshift=0)
2743 if (
associated(adp%dv_dt_dia)) adp%dv_dt_dia(i,j,1) = v(i,j,1)
2744 hval = (hold(i,j,1) + hold(i,j+1,1)) + (ea(i,j,1) + ea(i,j+1,1)) + h_neglect
2745 b1(i) = 1.0 / (hval + (eb(i,j,1) + eb(i,j+1,1)))
2746 d1(i) = hval * b1(i)
2747 v(i,j,1) = b1(i) * (hval * v(i,j,1))
2749 do k=2,nz ;
do i=is,ie
2750 if (
associated(adp%dv_dt_dia)) adp%dv_dt_dia(i,j,k) = v(i,j,k)
2751 c1(i,k) = (eb(i,j,k-1)+eb(i,j+1,k-1)) * b1(i)
2752 eaval = ea(i,j,k) + ea(i,j+1,k)
2753 hval = hold(i,j,k) + hold(i,j+1,k) + h_neglect
2754 b1(i) = 1.0 / ((eb(i,j,k) + eb(i,j+1,k)) + (hval + d1(i)*eaval))
2755 d1(i) = (hval + d1(i)*eaval) * b1(i)
2756 v(i,j,k) = (hval*v(i,j,k) + eaval*v(i,j,k-1))*b1(i)
2758 do k=nz-1,1,-1 ;
do i=is,ie
2759 v(i,j,k) = v(i,j,k) + c1(i,k+1)*v(i,j,k+1)
2760 if (
associated(adp%dv_dt_dia)) &
2761 adp%dv_dt_dia(i,j,k) = (v(i,j,k) - adp%dv_dt_dia(i,j,k)) * idt_accel
2763 if (
associated(adp%dv_dt_dia))
then 2765 adp%dv_dt_dia(i,j,nz) = (v(i,j,nz)-adp%dv_dt_dia(i,j,nz)) * idt_accel
2769 call cpu_clock_end(id_clock_tridiag)
2771 call mom_state_chksum(
"after u/v tridiag ", u, v, h, g, gv, us, haloshift=0)
2774 call disable_averaging(cs%diag)
2776 call enable_averages(dt, time_end, cs%diag)
2778 if (cs%id_Kd_interface > 0)
call post_data(cs%id_Kd_interface, kd_int, cs%diag)
2779 if (cs%id_Kd_heat > 0)
call post_data(cs%id_Kd_heat, kd_heat, cs%diag)
2780 if (cs%id_Kd_salt > 0)
call post_data(cs%id_Kd_salt, kd_salt, cs%diag)
2781 if (cs%id_Kd_ePBL > 0)
call post_data(cs%id_Kd_ePBL, kd_epbl, cs%diag)
2783 if (cs%id_ea > 0)
call post_data(cs%id_ea, ea, cs%diag)
2784 if (cs%id_eb > 0)
call post_data(cs%id_eb, eb, cs%diag)
2786 if (cs%id_dudt_dia > 0)
call post_data(cs%id_dudt_dia, adp%du_dt_dia, cs%diag)
2787 if (cs%id_dvdt_dia > 0)
call post_data(cs%id_dvdt_dia, adp%dv_dt_dia, cs%diag)
2788 if (cs%id_wd > 0)
call post_data(cs%id_wd, cdp%diapyc_vel, cs%diag)
2790 if (cs%id_Tdif > 0)
call post_data(cs%id_Tdif, tdif_flx, cs%diag)
2791 if (cs%id_Tadv > 0)
call post_data(cs%id_Tadv, tadv_flx, cs%diag)
2792 if (cs%id_Sdif > 0)
call post_data(cs%id_Sdif, sdif_flx, cs%diag)
2793 if (cs%id_Sadv > 0)
call post_data(cs%id_Sadv, sadv_flx, cs%diag)
2796 if (cs%id_hf_dudt_dia_2d > 0)
then 2797 allocate(hf_dudt_dia_2d(g%IsdB:g%IedB,g%jsd:g%jed))
2798 hf_dudt_dia_2d(:,:) = 0.0
2799 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
2800 hf_dudt_dia_2d(i,j) = hf_dudt_dia_2d(i,j) + adp%du_dt_dia(i,j,k) * adp%diag_hfrac_u(i,j,k)
2801 enddo ;
enddo ;
enddo 2802 call post_data(cs%id_hf_dudt_dia_2d, hf_dudt_dia_2d, cs%diag)
2803 deallocate(hf_dudt_dia_2d)
2806 if (cs%id_hf_dvdt_dia_2d > 0)
then 2807 allocate(hf_dvdt_dia_2d(g%isd:g%ied,g%JsdB:g%JedB))
2808 hf_dvdt_dia_2d(:,:) = 0.0
2809 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
2810 hf_dvdt_dia_2d(i,j) = hf_dvdt_dia_2d(i,j) + adp%dv_dt_dia(i,j,k) * adp%diag_hfrac_v(i,j,k)
2811 enddo ;
enddo ;
enddo 2812 call post_data(cs%id_hf_dvdt_dia_2d, hf_dvdt_dia_2d, cs%diag)
2813 deallocate(hf_dvdt_dia_2d)
2816 call disable_averaging(cs%diag)
2818 if (showcalltree)
call calltree_leave(
"layered_diabatic()")
2820 end subroutine layered_diabatic
2824 subroutine extract_diabatic_member(CS, opacity_CSp, optics_CSp, evap_CFL_limit, minimum_forcing_depth, &
2825 KPP_CSp, energetic_PBL_CSp, diabatic_aux_CSp, diabatic_halo)
2828 type(
opacity_cs),
optional,
pointer :: opacity_CSp
2829 type(
optics_type),
optional,
pointer :: optics_CSp
2830 type(
kpp_cs),
optional,
pointer :: KPP_CSp
2832 real,
optional,
intent( out) :: evap_CFL_limit
2834 real,
optional,
intent( out) :: minimum_forcing_depth
2838 integer,
optional,
intent( out) :: diabatic_halo
2842 if (
present(opacity_csp)) opacity_csp => cs%opacity_CSp
2843 if (
present(optics_csp)) optics_csp => cs%optics
2844 if (
present(kpp_csp)) kpp_csp => cs%KPP_CSp
2845 if (
present(energetic_pbl_csp)) energetic_pbl_csp => cs%energetic_PBL_CSp
2848 if (
present(evap_cfl_limit)) evap_cfl_limit = cs%evap_CFL_limit
2849 if (
present(minimum_forcing_depth)) minimum_forcing_depth = cs%minimum_forcing_depth
2850 if (
present(diabatic_halo)) diabatic_halo = cs%halo_TS_diff
2852 end subroutine extract_diabatic_member
2855 subroutine adiabatic(h, tv, fluxes, dt, G, GV, US, CS)
2857 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2860 type(
forcing),
intent(inout) :: fluxes
2861 real,
intent(in) :: dt
2866 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: zeros
2870 call call_tracer_column_fns(h, h, zeros, zeros, fluxes, zeros(:,:,1), dt, g, gv, us, tv, &
2871 cs%optics, cs%tracer_flow_CSp, cs%debug)
2873 end subroutine adiabatic
2879 subroutine diagnose_diabatic_diff_tendency(tv, h, temp_old, saln_old, dt, G, GV, US, CS)
2883 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
2884 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: temp_old
2885 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: saln_old
2886 real,
intent(in) :: dt
2891 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: work_3d
2892 real,
dimension(SZI_(G),SZJ_(G)) :: work_2d
2894 real :: ppt2mks = 0.001
2895 integer :: i, j, k, is, ie, js, je, nz
2896 logical :: do_saln_tend
2898 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
2899 idt = 0.0 ;
if (dt > 0.0) idt = 1. / dt
2900 work_3d(:,:,:) = 0.0
2905 do k=1,nz ;
do j=js,je ;
do i=is,ie
2906 work_3d(i,j,k) = (tv%T(i,j,k)-temp_old(i,j,k))*idt
2907 enddo ;
enddo ;
enddo 2908 if (cs%id_diabatic_diff_temp_tend > 0)
then 2909 call post_data(cs%id_diabatic_diff_temp_tend, work_3d, cs%diag, alt_h=h)
2913 if (cs%id_diabatic_diff_heat_tend > 0 .or. cs%id_diabatic_diff_heat_tend_2d > 0)
then 2914 do k=1,nz ;
do j=js,je ;
do i=is,ie
2915 work_3d(i,j,k) = h(i,j,k)*gv%H_to_RZ * tv%C_p * work_3d(i,j,k)
2916 enddo ;
enddo ;
enddo 2917 if (cs%id_diabatic_diff_heat_tend > 0)
then 2918 call post_data(cs%id_diabatic_diff_heat_tend, work_3d, cs%diag, alt_h=h)
2920 if (cs%id_diabatic_diff_heat_tend_2d > 0)
then 2921 do j=js,je ;
do i=is,ie
2924 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
2927 call post_data(cs%id_diabatic_diff_heat_tend_2d, work_2d, cs%diag)
2932 do_saln_tend = cs%id_diabatic_diff_saln_tend > 0 &
2933 .or. cs%id_diabatic_diff_salt_tend > 0 &
2934 .or. cs%id_diabatic_diff_salt_tend_2d > 0
2936 if (do_saln_tend)
then 2937 do k=1,nz ;
do j=js,je ;
do i=is,ie
2938 work_3d(i,j,k) = (tv%S(i,j,k) - saln_old(i,j,k)) * idt
2939 enddo ;
enddo ;
enddo 2941 if (cs%id_diabatic_diff_saln_tend > 0) &
2942 call post_data(cs%id_diabatic_diff_saln_tend, work_3d, cs%diag, alt_h=h)
2945 if (cs%id_diabatic_diff_salt_tend > 0 .or. cs%id_diabatic_diff_salt_tend_2d > 0)
then 2946 do k=1,nz ;
do j=js,je ;
do i=is,ie
2947 work_3d(i,j,k) = h(i,j,k)*gv%H_to_RZ * ppt2mks * work_3d(i,j,k)
2948 enddo ;
enddo ;
enddo 2949 if (cs%id_diabatic_diff_salt_tend > 0)
then 2950 call post_data(cs%id_diabatic_diff_salt_tend, work_3d, cs%diag, alt_h=h)
2952 if (cs%id_diabatic_diff_salt_tend_2d > 0)
then 2953 do j=js,je ;
do i=is,ie
2956 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
2959 call post_data(cs%id_diabatic_diff_salt_tend_2d, work_2d, cs%diag)
2964 end subroutine diagnose_diabatic_diff_tendency
2971 subroutine diagnose_boundary_forcing_tendency(tv, h, temp_old, saln_old, h_old, &
2976 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2978 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2979 intent(in) :: temp_old
2980 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2981 intent(in) :: saln_old
2982 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2984 real,
intent(in) :: dt
2989 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: work_3d
2990 real,
dimension(SZI_(G),SZJ_(G)) :: work_2d
2992 real :: ppt2mks = 0.001
2993 integer :: i, j, k, is, ie, js, je, nz
2995 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
2996 idt = 0.0 ;
if (dt > 0.0) idt = 1. / dt
2997 work_3d(:,:,:) = 0.0
3001 if (cs%id_boundary_forcing_h_tendency > 0)
then 3002 do k=1,nz ;
do j=js,je ;
do i=is,ie
3003 work_3d(i,j,k) = (h(i,j,k) - h_old(i,j,k))*idt
3004 enddo ;
enddo ;
enddo 3005 call post_data(cs%id_boundary_forcing_h_tendency, work_3d, cs%diag, alt_h=h_old)
3009 if (cs%id_boundary_forcing_temp_tend > 0)
then 3010 do k=1,nz ;
do j=js,je ;
do i=is,ie
3011 work_3d(i,j,k) = (tv%T(i,j,k)-temp_old(i,j,k))*idt
3012 enddo ;
enddo ;
enddo 3013 call post_data(cs%id_boundary_forcing_temp_tend, work_3d, cs%diag, alt_h=h_old)
3017 if (cs%id_boundary_forcing_heat_tend > 0 .or. cs%id_boundary_forcing_heat_tend_2d > 0)
then 3018 do k=1,nz ;
do j=js,je ;
do i=is,ie
3019 work_3d(i,j,k) = gv%H_to_RZ * tv%C_p * idt * (h(i,j,k) * tv%T(i,j,k) - h_old(i,j,k) * temp_old(i,j,k))
3020 enddo ;
enddo ;
enddo 3021 if (cs%id_boundary_forcing_heat_tend > 0)
then 3022 call post_data(cs%id_boundary_forcing_heat_tend, work_3d, cs%diag, alt_h = h_old)
3024 if (cs%id_boundary_forcing_heat_tend_2d > 0)
then 3025 do j=js,je ;
do i=is,ie
3028 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
3031 call post_data(cs%id_boundary_forcing_heat_tend_2d, work_2d, cs%diag)
3036 if (cs%id_boundary_forcing_saln_tend > 0)
then 3037 do k=1,nz ;
do j=js,je ;
do i=is,ie
3038 work_3d(i,j,k) = (tv%S(i,j,k)-saln_old(i,j,k))*idt
3039 enddo ;
enddo ;
enddo 3040 call post_data(cs%id_boundary_forcing_saln_tend, work_3d, cs%diag, alt_h = h_old)
3044 if (cs%id_boundary_forcing_salt_tend > 0 .or. cs%id_boundary_forcing_salt_tend_2d > 0)
then 3045 do k=1,nz ;
do j=js,je ;
do i=is,ie
3046 work_3d(i,j,k) = gv%H_to_RZ * ppt2mks * idt * (h(i,j,k) * tv%S(i,j,k) - h_old(i,j,k) * saln_old(i,j,k))
3047 enddo ;
enddo ;
enddo 3048 if (cs%id_boundary_forcing_salt_tend > 0)
then 3049 call post_data(cs%id_boundary_forcing_salt_tend, work_3d, cs%diag, alt_h = h_old)
3051 if (cs%id_boundary_forcing_salt_tend_2d > 0)
then 3052 do j=js,je ;
do i=is,ie
3055 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
3058 call post_data(cs%id_boundary_forcing_salt_tend_2d, work_2d, cs%diag)
3062 end subroutine diagnose_boundary_forcing_tendency
3069 subroutine diagnose_frazil_tendency(tv, h, temp_old, dt, G, GV, US, CS)
3073 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
3074 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: temp_old
3075 real,
intent(in) :: dt
3079 real,
dimension(SZI_(G),SZJ_(G)) :: work_2d
3081 integer :: i, j, k, is, ie, js, je, nz
3083 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
3084 idt = 0.0 ;
if (dt > 0.0) idt = 1. / dt
3087 if (cs%id_frazil_temp_tend > 0)
then 3088 do k=1,nz ;
do j=js,je ;
do i=is,ie
3089 cs%frazil_temp_diag(i,j,k) = idt * (tv%T(i,j,k)-temp_old(i,j,k))
3090 enddo ;
enddo ;
enddo 3091 call post_data(cs%id_frazil_temp_tend, cs%frazil_temp_diag(:,:,:), cs%diag)
3095 if (cs%id_frazil_heat_tend > 0 .or. cs%id_frazil_heat_tend_2d > 0)
then 3096 do k=1,nz ;
do j=js,je ;
do i=is,ie
3097 cs%frazil_heat_diag(i,j,k) = gv%H_to_RZ * tv%C_p * h(i,j,k) * idt * (tv%T(i,j,k)-temp_old(i,j,k))
3098 enddo ;
enddo ;
enddo 3099 if (cs%id_frazil_heat_tend > 0)
call post_data(cs%id_frazil_heat_tend, cs%frazil_heat_diag(:,:,:), cs%diag)
3103 if (cs%id_frazil_heat_tend_2d > 0)
then 3104 do j=js,je ;
do i=is,ie
3107 work_2d(i,j) = work_2d(i,j) + cs%frazil_heat_diag(i,j,k)
3110 call post_data(cs%id_frazil_heat_tend_2d, work_2d, cs%diag)
3114 end subroutine diagnose_frazil_tendency
3120 subroutine adiabatic_driver_init(Time, G, param_file, diag, CS, &
3122 type(time_type),
intent(in) :: Time
3125 type(
diag_ctrl),
target,
intent(inout) :: diag
3131 #include "version_variable.h" 3132 character(len=40) :: mdl =
"MOM_diabatic_driver" 3134 if (
associated(cs))
then 3135 call mom_error(warning,
"adiabatic_driver_init called with an "// &
3136 "associated control structure.")
3138 else ;
allocate(cs) ;
endif 3141 if (
associated(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
3145 "The following parameters are used for diabatic processes.")
3147 end subroutine adiabatic_driver_init
3151 subroutine diabatic_driver_init(Time, G, GV, US, param_file, useALEalgorithm, diag, &
3152 ADp, CDp, CS, tracer_flow_CSp, sponge_CSp, &
3154 type(time_type),
target :: Time
3159 logical,
intent(in) :: useALEalgorithm
3160 type(
diag_ctrl),
target,
intent(inout) :: diag
3172 logical :: use_temperature, differentialDiffusion
3173 character(len=20) :: EN1, EN2, EN3
3176 #include "version_variable.h" 3177 character(len=40) :: mdl =
"MOM_diabatic_driver" 3178 character(len=48) :: thickness_units
3179 character(len=40) :: var_name
3180 character(len=160) :: var_descript
3181 integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, nz, nbands, m
3182 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
3183 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
3185 if (
associated(cs))
then 3186 call mom_error(warning,
"diabatic_driver_init called with an "// &
3187 "associated control structure.")
3196 if (
associated(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
3197 if (
associated(sponge_csp)) cs%sponge_CSp => sponge_csp
3198 if (
associated(ale_sponge_csp)) cs%ALE_sponge_CSp => ale_sponge_csp
3200 cs%useALEalgorithm = usealealgorithm
3201 cs%bulkmixedlayer = (gv%nkml > 0)
3205 "The following parameters are used for diabatic processes.", &
3206 log_to_all=.true., debugging=.true.)
3207 call get_param(param_file, mdl,
"USE_LEGACY_DIABATIC_DRIVER", cs%use_legacy_diabatic, &
3208 "If true, use a legacy version of the diabatic subroutine. "//&
3209 "This is temporary and is needed to avoid change in answers.", &
3211 call get_param(param_file, mdl,
"SPONGE", cs%use_sponge, &
3212 "If true, sponges may be applied anywhere in the domain. "//&
3213 "The exact location and properties of those sponges are "//&
3214 "specified via calls to initialize_sponge and possibly "//&
3215 "set_up_sponge_field.", default=.false.)
3216 call get_param(param_file, mdl,
"ENABLE_THERMODYNAMICS", use_temperature, &
3217 "If true, temperature and salinity are used as state "//&
3218 "variables.", default=.true.)
3219 call get_param(param_file, mdl,
"ENERGETICS_SFC_PBL", cs%use_energetic_PBL, &
3220 "If true, use an implied energetics planetary boundary "//&
3221 "layer scheme to determine the diffusivity and viscosity "//&
3222 "in the surface boundary layer.", default=.false.)
3223 call get_param(param_file, mdl,
"EPBL_IS_ADDITIVE", cs%ePBL_is_additive, &
3224 "If true, the diffusivity from ePBL is added to all "//&
3225 "other diffusivities. Otherwise, the larger of kappa-shear "//&
3226 "and ePBL diffusivities are used.", default=.true.)
3227 call get_param(param_file, mdl,
"PRANDTL_EPBL", cs%ePBL_Prandtl, &
3228 "The Prandtl number used by ePBL to convert vertical diffusivities into "//&
3229 "viscosities.", default=1.0, units=
"nondim", do_not_log=.not.cs%use_energetic_PBL)
3230 call get_param(param_file, mdl,
"DOUBLE_DIFFUSION", differentialdiffusion, &
3231 "If true, apply parameterization of double-diffusion.", &
3233 call get_param(param_file, mdl,
"USE_KPP", cs%use_KPP, &
3234 "If true, turns on the [CVMix] KPP scheme of Large et al., 1994, "//&
3235 "to calculate diffusivities and non-local transport in the OBL.", &
3236 default=.false., do_not_log=.true.)
3237 cs%use_CVMix_ddiff = cvmix_ddiff_is_used(param_file)
3239 if (cs%use_CVMix_ddiff .and. differentialdiffusion)
then 3240 call mom_error(fatal,
'diabatic_driver_init: '// &
3241 'Multiple double-diffusion options selected (DOUBLE_DIFFUSION and'//&
3242 'USE_CVMIX_DDIFF), please disable all but one option to proceed.')
3245 cs%use_kappa_shear = kappa_shear_is_used(param_file)
3246 cs%use_CVMix_shear = cvmix_shear_is_used(param_file)
3248 if (cs%bulkmixedlayer)
then 3249 call get_param(param_file, mdl,
"ML_MIX_FIRST", cs%ML_mix_first, &
3250 "The fraction of the mixed layer mixing that is applied "//&
3251 "before interior diapycnal mixing. 0 by default.", &
3252 units=
"nondim", default=0.0)
3253 call get_param(param_file, mdl,
"NKBL", cs%nkbl, default=2, do_not_log=.true.)
3255 cs%ML_mix_first = 0.0
3257 if (use_temperature)
then 3258 call get_param(param_file, mdl,
"DO_GEOTHERMAL", cs%use_geothermal, &
3259 "If true, apply geothermal heating.", default=.false.)
3261 cs%use_geothermal = .false.
3263 call get_param(param_file, mdl,
"INTERNAL_TIDES", cs%use_int_tides, &
3264 "If true, use the code that advances a separate set of "//&
3265 "equations for the internal tide energy density.", default=.false.)
3267 if (cs%use_int_tides)
then 3268 call get_param(param_file, mdl,
"INTERNAL_TIDE_MODES", cs%nMode, &
3269 "The number of distinct internal tide modes "//&
3270 "that will be calculated.", default=1, do_not_log=.true.)
3271 call get_param(param_file, mdl,
"UNIFORM_TEST_CG", cs%uniform_test_cg, &
3272 "If positive, a uniform group velocity of internal tide for test case", &
3273 default=-1., units=
"m s-1", scale=us%m_s_to_L_T)
3276 call get_param(param_file, mdl,
"MASSLESS_MATCH_TARGETS", &
3277 cs%massless_match_targets, &
3278 "If true, the temperature and salinity of massless layers "//&
3279 "are kept consistent with their target densities. "//&
3280 "Otherwise the properties of massless layers evolve "//&
3281 "diffusively to match massive neighboring layers.", &
3284 call get_param(param_file, mdl,
"AGGREGATE_FW_FORCING", cs%aggregate_FW_forcing, &
3285 "If true, the net incoming and outgoing fresh water fluxes are combined "//&
3286 "and applied as either incoming or outgoing depending on the sign of the net. "//&
3287 "If false, the net incoming fresh water flux is added to the model and "//&
3288 "thereafter the net outgoing is removed from the topmost non-vanished "//&
3289 "layers of the updated state.", default=.true.)
3291 call get_param(param_file, mdl,
"DEBUG", cs%debug, &
3292 "If true, write out verbose debugging data.", &
3293 default=.false., debuggingparam=.true.)
3294 call get_param(param_file, mdl,
"DEBUG_CONSERVATION", cs%debugConservation, &
3295 "If true, monitor conservation and extrema.", &
3296 default=.false., debuggingparam=.true.)
3298 call get_param(param_file, mdl,
"DEBUG_ENERGY_REQ", cs%debug_energy_req, &
3299 "If true, debug the energy requirements.", default=.false., do_not_log=.true.)
3300 call get_param(param_file, mdl,
"MIX_BOUNDARY_TRACERS", cs%mix_boundary_tracers, &
3301 "If true, mix the passive tracers in massless layers at "//&
3302 "the bottom into the interior as though a diffusivity of "//&
3303 "KD_MIN_TR were operating.", default=.true.)
3305 if (cs%mix_boundary_tracers)
then 3306 call get_param(param_file, mdl,
"KD", kd, default=0.0)
3307 call get_param(param_file, mdl,
"KD_MIN_TR", cs%Kd_min_tr, &
3308 "A minimal diffusivity that should always be applied to "//&
3309 "tracers, especially in massless layers near the bottom. "//&
3310 "The default is 0.1*KD.", units=
"m2 s-1", default=0.1*kd, scale=us%m2_s_to_Z2_T)
3311 call get_param(param_file, mdl,
"KD_BBL_TR", cs%Kd_BBL_tr, &
3312 "A bottom boundary layer tracer diffusivity that will "//&
3313 "allow for explicitly specified bottom fluxes. The "//&
3314 "entrainment at the bottom is at least sqrt(Kd_BBL_tr*dt) "//&
3315 "over the same distance.", units=
"m2 s-1", default=0., scale=us%m2_s_to_Z2_T)
3318 call get_param(param_file, mdl,
"TRACER_TRIDIAG", cs%tracer_tridiag, &
3319 "If true, use the passive tracer tridiagonal solver for T and S", &
3322 call get_param(param_file, mdl,
"MINIMUM_FORCING_DEPTH", cs%minimum_forcing_depth, &
3323 "The smallest depth over which forcing can be applied. This "//&
3324 "only takes effect when near-surface layers become thin "//&
3325 "relative to this scale, in which case the forcing tendencies "//&
3326 "scaled down by distributing the forcing over this depth scale.", &
3327 units=
"m", default=0.001, scale=gv%m_to_H)
3328 call get_param(param_file, mdl,
"EVAP_CFL_LIMIT", cs%evap_CFL_limit, &
3329 "The largest fraction of a layer than can be lost to forcing "//&
3330 "(e.g. evaporation, sea-ice formation) in one time-step. The unused "//&
3331 "mass loss is passed down through the column.", &
3332 units=
"nondim", default=0.8)
3336 thickness_units = get_thickness_units(gv)
3338 cs%id_ea_t = register_diag_field(
'ocean_model',
'ea_t', diag%axesTL, time, &
3339 'Layer (heat) entrainment from above per timestep',
'm', conversion=gv%H_to_m)
3340 cs%id_eb_t = register_diag_field(
'ocean_model',
'eb_t', diag%axesTL, time, &
3341 'Layer (heat) entrainment from below per timestep',
'm', conversion=gv%H_to_m)
3342 cs%id_ea_s = register_diag_field(
'ocean_model',
'ea_s', diag%axesTL, time, &
3343 'Layer (salt) entrainment from above per timestep',
'm', conversion=gv%H_to_m)
3344 cs%id_eb_s = register_diag_field(
'ocean_model',
'eb_s', diag%axesTL, time, &
3345 'Layer (salt) entrainment from below per timestep',
'm', conversion=gv%H_to_m)
3347 cs%id_ea = register_diag_field(
'ocean_model',
'ea', diag%axesTL, time, &
3348 'Layer entrainment from above per timestep',
'm', conversion=gv%H_to_m)
3349 cs%id_eb = register_diag_field(
'ocean_model',
'eb', diag%axesTL, time, &
3350 'Layer entrainment from below per timestep',
'm', conversion=gv%H_to_m)
3351 cs%id_wd = register_diag_field(
'ocean_model',
'wd', diag%axesTi, time, &
3352 'Diapycnal velocity',
'm s-1', conversion=gv%H_to_m)
3353 if (cs%id_wd > 0)
call safe_alloc_ptr(cdp%diapyc_vel,isd,ied,jsd,jed,nz+1)
3355 cs%id_dudt_dia = register_diag_field(
'ocean_model',
'dudt_dia', diag%axesCuL, time, &
3356 'Zonal Acceleration from Diapycnal Mixing',
'm s-2', conversion=us%L_T2_to_m_s2)
3357 cs%id_dvdt_dia = register_diag_field(
'ocean_model',
'dvdt_dia', diag%axesCvL, time, &
3358 'Meridional Acceleration from Diapycnal Mixing',
'm s-2', conversion=us%L_T2_to_m_s2)
3360 cs%id_hf_dudt_dia_2d = register_diag_field(
'ocean_model',
'hf_dudt_dia_2d', diag%axesCu1, time, &
3361 'Depth-sum Fractional Thickness-weighted Zonal Acceleration from Diapycnal Mixing',
'm s-2', &
3362 conversion=us%L_T2_to_m_s2)
3363 if (cs%id_hf_dudt_dia_2d > 0)
then 3364 call safe_alloc_ptr(adp%diag_hfrac_u,isdb,iedb,jsd,jed,nz)
3365 call safe_alloc_ptr(adp%du_dt_dia,isdb,iedb,jsd,jed,nz)
3368 cs%id_hf_dvdt_dia_2d = register_diag_field(
'ocean_model',
'hf_dvdt_dia_2d', diag%axesCv1, time, &
3369 'Depth-sum Fractional Thickness-weighted Meridional Acceleration from Diapycnal Mixing',
'm s-2', &
3370 conversion=us%L_T2_to_m_s2)
3371 if (cs%id_hf_dvdt_dia_2d > 0)
then 3372 call safe_alloc_ptr(adp%diag_hfrac_v,isd,ied,jsd,jedb,nz)
3373 call safe_alloc_ptr(adp%dv_dt_dia,isd,ied,jsdb,jedb,nz)
3376 if (cs%use_int_tides)
then 3377 cs%id_cg1 = register_diag_field(
'ocean_model',
'cn1', diag%axesT1, &
3378 time,
'First baroclinic mode (eigen) speed',
'm s-1', conversion=us%L_T_to_m_s)
3379 allocate(cs%id_cn(cs%nMode)) ; cs%id_cn(:) = -1
3381 write(var_name,
'("cn_mode",i1)') m
3382 write(var_descript,
'("Baroclinic (eigen) speed of mode ",i1)') m
3383 cs%id_cn(m) = register_diag_field(
'ocean_model',var_name, diag%axesT1, &
3384 time, var_descript,
'm s-1', conversion=us%L_T_to_m_s)
3385 call mom_mesg(
"Registering "//trim(var_name)//
", Described as: "//var_descript, 5)
3389 if (use_temperature)
then 3390 cs%id_Tdif = register_diag_field(
'ocean_model',
"Tflx_dia_diff", diag%axesTi, &
3391 time,
"Diffusive diapycnal temperature flux across interfaces", &
3392 "degC m s-1", conversion=gv%H_to_m*us%s_to_T)
3393 cs%id_Tadv = register_diag_field(
'ocean_model',
"Tflx_dia_adv", diag%axesTi, &
3394 time,
"Advective diapycnal temperature flux across interfaces", &
3395 "degC m s-1", conversion=gv%H_to_m*us%s_to_T)
3396 cs%id_Sdif = register_diag_field(
'ocean_model',
"Sflx_dia_diff", diag%axesTi, &
3397 time,
"Diffusive diapycnal salnity flux across interfaces", &
3398 "psu m s-1", conversion=gv%H_to_m*us%s_to_T)
3399 cs%id_Sadv = register_diag_field(
'ocean_model',
"Sflx_dia_adv", diag%axesTi, &
3400 time,
"Advective diapycnal salnity flux across interfaces", &
3401 "psu m s-1", conversion=gv%H_to_m*us%s_to_T)
3402 cs%id_MLD_003 = register_diag_field(
'ocean_model',
'MLD_003', diag%axesT1, time, &
3403 'Mixed layer depth (delta rho = 0.03)',
'm', conversion=us%Z_to_m, &
3404 cmor_field_name=
'mlotst', cmor_long_name=
'Ocean Mixed Layer Thickness Defined by Sigma T', &
3405 cmor_standard_name=
'ocean_mixed_layer_thickness_defined_by_sigma_t')
3406 cs%id_mlotstsq = register_diag_field(
'ocean_model',
'mlotstsq', diag%axesT1, time, &
3407 long_name=
'Square of Ocean Mixed Layer Thickness Defined by Sigma T', &
3408 standard_name=
'square_of_ocean_mixed_layer_thickness_defined_by_sigma_t', &
3409 units=
'm2', conversion=us%Z_to_m**2)
3410 cs%id_MLD_0125 = register_diag_field(
'ocean_model',
'MLD_0125', diag%axesT1, time, &
3411 'Mixed layer depth (delta rho = 0.125)',
'm', conversion=us%Z_to_m)
3412 call get_param(param_file, mdl,
"MLD_EN_VALS", cs%MLD_EN_VALS, &
3413 "The energy values used to compute MLDs. If not set (or all set to 0.), the "//&
3414 "default will overwrite to 25., 2500., 250000.",units=
'J/m2', default=0., &
3415 scale=us%kg_m3_to_R*us%m_to_Z**3*us%T_to_s**2)
3416 if ((cs%MLD_EN_VALS(1)==0.).and.(cs%MLD_EN_VALS(2)==0.).and.(cs%MLD_EN_VALS(3)==0.))
then 3417 cs%MLD_EN_VALS = (/25.*us%kg_m3_to_R*us%m_to_Z*us%m_to_L**2*us%T_to_s**2,&
3418 2500.*us%kg_m3_to_R*us%m_to_Z*us%m_to_L**2*us%T_to_s**2,&
3419 250000.*us%kg_m3_to_R*us%m_to_Z*us%m_to_L**2*us%T_to_s**2/)
3421 write(en1,
'(F10.2)') cs%MLD_EN_VALS(1)*us%R_to_kg_m3*us%Z_to_m*us%L_to_m**2*us%s_to_T**2
3422 write(en2,
'(F10.2)') cs%MLD_EN_VALS(2)*us%R_to_kg_m3*us%Z_to_m*us%L_to_m**2*us%s_to_T**2
3423 write(en3,
'(F10.2)') cs%MLD_EN_VALS(3)*us%R_to_kg_m3*us%Z_to_m*us%L_to_m**2*us%s_to_T**2
3424 cs%id_MLD_EN1 = register_diag_field(
'ocean_model',
'MLD_EN1', diag%axesT1, time, &
3425 'Mixed layer depth for energy value set to '//trim(en1)//
' J/m2 (Energy set by 1st MLD_EN_VALS)', &
3426 'm', conversion=us%Z_to_m)
3427 cs%id_MLD_EN2 = register_diag_field(
'ocean_model',
'MLD_EN2', diag%axesT1, time, &
3428 'Mixed layer depth for energy value set to '//trim(en2)//
' J/m2 (Energy set by 2nd MLD_EN_VALS)', &
3429 'm', conversion=us%Z_to_m)
3430 cs%id_MLD_EN3 = register_diag_field(
'ocean_model',
'MLD_EN3', diag%axesT1, time, &
3431 'Mixed layer depth for energy value set to '//trim(en3)//
' J/m2 (Energy set by 3rd MLD_EN_VALS)', &
3432 'm', conversion=us%Z_to_m)
3433 cs%id_subMLN2 = register_diag_field(
'ocean_model',
'subML_N2', diag%axesT1, time, &
3434 'Squared buoyancy frequency below mixed layer',
's-2', conversion=us%s_to_T**2)
3435 cs%id_MLD_user = register_diag_field(
'ocean_model',
'MLD_user', diag%axesT1, time, &
3436 'Mixed layer depth (used defined)',
'm', conversion=us%Z_to_m)
3438 call get_param(param_file, mdl,
"DIAG_MLD_DENSITY_DIFF", cs%MLDdensityDifference, &
3439 "The density difference used to determine a diagnostic mixed "//&
3440 "layer depth, MLD_user, following the definition of Levitus 1982. "//&
3441 "The MLD is the depth at which the density is larger than the "//&
3442 "surface density by the specified amount.", &
3443 units=
'kg/m3', default=0.1, scale=us%kg_m3_to_R)
3444 call get_param(param_file, mdl,
"DIAG_DEPTH_SUBML_N2", cs%dz_subML_N2, &
3445 "The distance over which to calculate a diagnostic of the "//&
3446 "stratification at the base of the mixed layer.", &
3447 units=
'm', default=50.0, scale=us%m_to_Z)
3449 if (cs%id_dudt_dia > 0)
call safe_alloc_ptr(adp%du_dt_dia,isdb,iedb,jsd,jed,nz)
3450 if (cs%id_dvdt_dia > 0)
call safe_alloc_ptr(adp%dv_dt_dia,isd,ied,jsdb,jedb,nz)
3453 cs%id_u_predia = register_diag_field(
'ocean_model',
'u_predia', diag%axesCuL, time, &
3454 'Zonal velocity before diabatic forcing',
'm s-1', conversion=us%L_T_to_m_s)
3455 cs%id_v_predia = register_diag_field(
'ocean_model',
'v_predia', diag%axesCvL, time, &
3456 'Meridional velocity before diabatic forcing',
'm s-1', conversion=us%L_T_to_m_s)
3457 cs%id_h_predia = register_diag_field(
'ocean_model',
'h_predia', diag%axesTL, time, &
3458 'Layer Thickness before diabatic forcing', &
3459 trim(thickness_units), conversion=gv%H_to_MKS, v_extensive=.true.)
3460 cs%id_e_predia = register_diag_field(
'ocean_model',
'e_predia', diag%axesTi, time, &
3461 'Interface Heights before diabatic forcing',
'm')
3462 if (use_temperature)
then 3463 cs%id_T_predia = register_diag_field(
'ocean_model',
'temp_predia', diag%axesTL, time, &
3464 'Potential Temperature',
'degC')
3465 cs%id_S_predia = register_diag_field(
'ocean_model',
'salt_predia', diag%axesTL, time, &
3471 cs%id_Kd_interface = register_diag_field(
'ocean_model',
'Kd_interface', diag%axesTi, time, &
3472 'Total diapycnal diffusivity at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
3473 if (cs%use_energetic_PBL)
then 3474 cs%id_Kd_ePBL = register_diag_field(
'ocean_model',
'Kd_ePBL', diag%axesTi, time, &
3475 'ePBL diapycnal diffusivity at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
3478 cs%id_Kd_heat = register_diag_field(
'ocean_model',
'Kd_heat', diag%axesTi, time, &
3479 'Total diapycnal diffusivity for heat at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s, &
3480 cmor_field_name=
'difvho', &
3481 cmor_standard_name=
'ocean_vertical_heat_diffusivity', &
3482 cmor_long_name=
'Ocean vertical heat diffusivity')
3483 cs%id_Kd_salt = register_diag_field(
'ocean_model',
'Kd_salt', diag%axesTi, time, &
3484 'Total diapycnal diffusivity for salt at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s, &
3485 cmor_field_name=
'difvso', &
3486 cmor_standard_name=
'ocean_vertical_salt_diffusivity', &
3487 cmor_long_name=
'Ocean vertical salt diffusivity')
3491 cs%useKPP = kpp_init(param_file, g, gv, us, diag, time, cs%KPP_CSp, passive=cs%KPPisPassive)
3493 allocate( cs%KPP_NLTheat(isd:ied,jsd:jed,nz+1) ) ; cs%KPP_NLTheat(:,:,:) = 0.
3494 allocate( cs%KPP_NLTscalar(isd:ied,jsd:jed,nz+1) ) ; cs%KPP_NLTscalar(:,:,:) = 0.
3497 allocate( cs%KPP_buoy_flux(isd:ied,jsd:jed,nz+1) ) ; cs%KPP_buoy_flux(:,:,:) = 0.
3498 allocate( cs%KPP_temp_flux(isd:ied,jsd:jed) ) ; cs%KPP_temp_flux(:,:) = 0.
3499 allocate( cs%KPP_salt_flux(isd:ied,jsd:jed) ) ; cs%KPP_salt_flux(:,:) = 0.
3502 if (cs%useKPP .and. differentialdiffusion)
then 3503 call mom_error(fatal,
'diabatic_driver_init: '// &
3504 'DOUBLE_DIFFUSION (old method) does not work with KPP. Please'//&
3505 'set DOUBLE_DIFFUSION=False and USE_CVMIX_DDIFF=True.')
3511 cs%id_diabatic_diff_h = register_diag_field(
'ocean_model',
'diabatic_diff_h', diag%axesTL, time, &
3512 long_name=
'Cell thickness used during diabatic diffusion', &
3513 units=
'm', conversion=gv%H_to_m, v_extensive=.true.)
3514 if (cs%useALEalgorithm)
then 3515 cs%id_diabatic_diff_temp_tend = register_diag_field(
'ocean_model', &
3516 'diabatic_diff_temp_tendency', diag%axesTL, time, &
3517 'Diabatic diffusion temperature tendency',
'degC s-1', conversion=us%s_to_T)
3518 if (cs%id_diabatic_diff_temp_tend > 0)
then 3519 cs%diabatic_diff_tendency_diag = .true.
3522 cs%id_diabatic_diff_saln_tend = register_diag_field(
'ocean_model',&
3523 'diabatic_diff_saln_tendency', diag%axesTL, time, &
3524 'Diabatic diffusion salinity tendency',
'psu s-1', conversion=us%s_to_T)
3525 if (cs%id_diabatic_diff_saln_tend > 0)
then 3526 cs%diabatic_diff_tendency_diag = .true.
3529 cs%id_diabatic_diff_heat_tend = register_diag_field(
'ocean_model', &
3530 'diabatic_heat_tendency', diag%axesTL, time, &
3531 'Diabatic diffusion heat tendency', &
3532 'W m-2', conversion=us%QRZ_T_to_W_m2, cmor_field_name=
'opottempdiff', &
3533 cmor_standard_name=
'tendency_of_sea_water_potential_temperature_expressed_as_heat_content_'// &
3534 'due_to_parameterized_dianeutral_mixing', &
3535 cmor_long_name=
'Tendency of sea water potential temperature expressed as heat content '// &
3536 'due to parameterized dianeutral mixing', &
3538 if (cs%id_diabatic_diff_heat_tend > 0)
then 3539 cs%diabatic_diff_tendency_diag = .true.
3542 cs%id_diabatic_diff_salt_tend = register_diag_field(
'ocean_model', &
3543 'diabatic_salt_tendency', diag%axesTL, time, &
3544 'Diabatic diffusion of salt tendency', &
3545 'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s, cmor_field_name=
'osaltdiff', &
3546 cmor_standard_name=
'tendency_of_sea_water_salinity_expressed_as_salt_content_'// &
3547 'due_to_parameterized_dianeutral_mixing', &
3548 cmor_long_name=
'Tendency of sea water salinity expressed as salt content '// &
3549 'due to parameterized dianeutral mixing', &
3551 if (cs%id_diabatic_diff_salt_tend > 0)
then 3552 cs%diabatic_diff_tendency_diag = .true.
3556 cs%id_diabatic_diff_heat_tend_2d = register_diag_field(
'ocean_model', &
3557 'diabatic_heat_tendency_2d', diag%axesT1, time, &
3558 'Depth integrated diabatic diffusion heat tendency', &
3559 'W m-2', conversion=us%QRZ_T_to_W_m2, cmor_field_name=
'opottempdiff_2d', &
3560 cmor_standard_name=
'tendency_of_sea_water_potential_temperature_expressed_as_heat_content_'//&
3561 'due_to_parameterized_dianeutral_mixing_depth_integrated', &
3562 cmor_long_name=
'Tendency of sea water potential temperature expressed as heat content '//&
3563 'due to parameterized dianeutral mixing depth integrated')
3564 if (cs%id_diabatic_diff_heat_tend_2d > 0)
then 3565 cs%diabatic_diff_tendency_diag = .true.
3569 cs%id_diabatic_diff_salt_tend_2d = register_diag_field(
'ocean_model', &
3570 'diabatic_salt_tendency_2d', diag%axesT1, time, &
3571 'Depth integrated diabatic diffusion salt tendency', &
3572 'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s, cmor_field_name=
'osaltdiff_2d', &
3573 cmor_standard_name=
'tendency_of_sea_water_salinity_expressed_as_salt_content_'// &
3574 'due_to_parameterized_dianeutral_mixing_depth_integrated', &
3575 cmor_long_name=
'Tendency of sea water salinity expressed as salt content '// &
3576 'due to parameterized dianeutral mixing depth integrated')
3577 if (cs%id_diabatic_diff_salt_tend_2d > 0)
then 3578 cs%diabatic_diff_tendency_diag = .true.
3584 cs%id_boundary_forcing_h = register_diag_field(
'ocean_model',
'boundary_forcing_h', diag%axesTL, time, &
3585 long_name=
'Cell thickness after applying boundary forcing', &
3586 units=
'm', conversion=gv%H_to_m, v_extensive=.true.)
3587 cs%id_boundary_forcing_h_tendency = register_diag_field(
'ocean_model', &
3588 'boundary_forcing_h_tendency', diag%axesTL, time, &
3589 'Cell thickness tendency due to boundary forcing',
'm s-1', &
3590 conversion=gv%H_to_m*us%s_to_T, v_extensive=.true.)
3591 if (cs%id_boundary_forcing_h_tendency > 0)
then 3592 cs%boundary_forcing_tendency_diag = .true.
3595 cs%id_boundary_forcing_temp_tend = register_diag_field(
'ocean_model',&
3596 'boundary_forcing_temp_tendency', diag%axesTL, time, &
3597 'Boundary forcing temperature tendency',
'degC s-1', conversion=us%s_to_T)
3598 if (cs%id_boundary_forcing_temp_tend > 0)
then 3599 cs%boundary_forcing_tendency_diag = .true.
3602 cs%id_boundary_forcing_saln_tend = register_diag_field(
'ocean_model',&
3603 'boundary_forcing_saln_tendency', diag%axesTL, time, &
3604 'Boundary forcing saln tendency',
'psu s-1', conversion=us%s_to_T)
3605 if (cs%id_boundary_forcing_saln_tend > 0)
then 3606 cs%boundary_forcing_tendency_diag = .true.
3609 cs%id_boundary_forcing_heat_tend = register_diag_field(
'ocean_model',&
3610 'boundary_forcing_heat_tendency', diag%axesTL, time, &
3611 'Boundary forcing heat tendency', &
3612 'W m-2', conversion=us%QRZ_T_to_W_m2, v_extensive = .true.)
3613 if (cs%id_boundary_forcing_heat_tend > 0)
then 3614 cs%boundary_forcing_tendency_diag = .true.
3617 cs%id_boundary_forcing_salt_tend = register_diag_field(
'ocean_model',&
3618 'boundary_forcing_salt_tendency', diag%axesTL, time, &
3619 'Boundary forcing salt tendency',
'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s, &
3620 v_extensive = .true.)
3621 if (cs%id_boundary_forcing_salt_tend > 0)
then 3622 cs%boundary_forcing_tendency_diag = .true.
3626 cs%id_boundary_forcing_heat_tend_2d = register_diag_field(
'ocean_model',&
3627 'boundary_forcing_heat_tendency_2d', diag%axesT1, time, &
3628 'Depth integrated boundary forcing of ocean heat', &
3629 'W m-2', conversion=us%QRZ_T_to_W_m2)
3630 if (cs%id_boundary_forcing_heat_tend_2d > 0)
then 3631 cs%boundary_forcing_tendency_diag = .true.
3635 cs%id_boundary_forcing_salt_tend_2d = register_diag_field(
'ocean_model',&
3636 'boundary_forcing_salt_tendency_2d', diag%axesT1, time, &
3637 'Depth integrated boundary forcing of ocean salt', &
3638 'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s)
3639 if (cs%id_boundary_forcing_salt_tend_2d > 0)
then 3640 cs%boundary_forcing_tendency_diag = .true.
3645 cs%id_frazil_h = register_diag_field(
'ocean_model',
'frazil_h', diag%axesTL, time, &
3646 long_name=
'Cell Thickness', standard_name=
'cell_thickness', &
3647 units=
'm', conversion=gv%H_to_m, v_extensive=.true.)
3650 cs%id_frazil_temp_tend = register_diag_field(
'ocean_model',&
3651 'frazil_temp_tendency', diag%axesTL, time, &
3652 'Temperature tendency due to frazil formation',
'degC s-1', conversion=us%s_to_T)
3653 if (cs%id_frazil_temp_tend > 0)
then 3654 cs%frazil_tendency_diag = .true.
3658 cs%id_frazil_heat_tend = register_diag_field(
'ocean_model',&
3659 'frazil_heat_tendency', diag%axesTL, time, &
3660 'Heat tendency due to frazil formation', &
3661 'W m-2', conversion=us%QRZ_T_to_W_m2, v_extensive=.true.)
3662 if (cs%id_frazil_heat_tend > 0)
then 3663 cs%frazil_tendency_diag = .true.
3667 cs%id_frazil_heat_tend_2d = register_diag_field(
'ocean_model',&
3668 'frazil_heat_tendency_2d', diag%axesT1, time, &
3669 'Depth integrated heat tendency due to frazil formation', &
3670 'W m-2', conversion=us%QRZ_T_to_W_m2)
3671 if (cs%id_frazil_heat_tend_2d > 0)
then 3672 cs%frazil_tendency_diag = .true.
3675 if (cs%frazil_tendency_diag)
then 3676 allocate(cs%frazil_temp_diag(isd:ied,jsd:jed,nz) ) ; cs%frazil_temp_diag(:,:,:) = 0.
3677 allocate(cs%frazil_heat_diag(isd:ied,jsd:jed,nz) ) ; cs%frazil_heat_diag(:,:,:) = 0.
3681 cs%use_CVMix_conv = cvmix_conv_init(time, g, gv, us, param_file, diag, cs%CVMix_conv_csp)
3683 call entrain_diffusive_init(time, g, gv, us, param_file, diag, cs%entrain_diffusive_CSp, &
3684 just_read_params=cs%useALEalgorithm)
3687 if (cs%use_geothermal) &
3688 call geothermal_init(time, g, gv, us, param_file, diag, cs%geothermal_CSp, usealealgorithm)
3691 if (cs%use_int_tides)
then 3692 call int_tide_input_init(time, g, gv, us, param_file, diag, cs%int_tide_input_CSp, &
3694 call internal_tides_init(time, g, gv, us, param_file, diag, cs%int_tide_CSp)
3698 call set_diffusivity_init(time, g, gv, us, param_file, diag, cs%set_diff_CSp, &
3699 cs%int_tide_CSp, cs%halo_TS_diff)
3702 id_clock_entrain = cpu_clock_id(
'(Ocean diabatic entrain)', grain=clock_module)
3703 if (cs%bulkmixedlayer) &
3704 id_clock_mixedlayer = cpu_clock_id(
'(Ocean mixed layer)', grain=clock_module)
3705 id_clock_remap = cpu_clock_id(
'(Ocean vert remap)', grain=clock_module)
3706 if (cs%use_geothermal) &
3707 id_clock_geothermal = cpu_clock_id(
'(Ocean geothermal)', grain=clock_routine)
3708 id_clock_set_diffusivity = cpu_clock_id(
'(Ocean set_diffusivity)', grain=clock_module)
3709 id_clock_kpp = cpu_clock_id(
'(Ocean KPP)', grain=clock_module)
3710 id_clock_tracers = cpu_clock_id(
'(Ocean tracer_columns)', grain=clock_module_driver+5)
3711 if (cs%use_sponge) &
3712 id_clock_sponge = cpu_clock_id(
'(Ocean sponges)', grain=clock_module)
3713 id_clock_tridiag = cpu_clock_id(
'(Ocean diabatic tridiag)', grain=clock_routine)
3714 id_clock_pass = cpu_clock_id(
'(Ocean diabatic message passing)', grain=clock_routine)
3715 id_clock_differential_diff = -1 ;
if (differentialdiffusion) &
3716 id_clock_differential_diff = cpu_clock_id(
'(Ocean differential diffusion)', grain=clock_routine)
3719 call diabatic_aux_init(time, g, gv, us, param_file, diag, cs%diabatic_aux_CSp, &
3720 cs%useALEalgorithm, cs%use_energetic_PBL)
3723 if (cs%bulkmixedlayer) &
3724 call bulkmixedlayer_init(time, g, gv, us, param_file, diag, cs%bulkmixedlayer_CSp)
3725 if (cs%use_energetic_PBL) &
3726 call energetic_pbl_init(time, g, gv, us, param_file, diag, cs%energetic_PBL_CSp)
3728 call regularize_layers_init(time, g, gv, param_file, diag, cs%regularize_layers_CSp)
3730 if (cs%debug_energy_req) &
3731 call diapyc_energy_req_init(time, g, gv, us, param_file, diag, cs%diapyc_en_rec_CSp)
3734 if (use_temperature)
then 3735 call get_param(param_file, mdl,
"PEN_SW_NBANDS", nbands, default=1)
3736 if (nbands > 0)
then 3738 call opacity_init(time, g, gv, us, param_file, diag, cs%opacity_CSp, cs%optics)
3743 call diag_grid_storage_init(cs%diag_grids_prev, g, diag)
3745 end subroutine diabatic_driver_init
3749 subroutine diabatic_driver_end(CS)
3752 if (.not.
associated(cs))
return 3754 call diabatic_aux_end(cs%diabatic_aux_CSp)
3756 call entrain_diffusive_end(cs%entrain_diffusive_CSp)
3757 call set_diffusivity_end(cs%set_diff_CSp)
3760 deallocate( cs%KPP_buoy_flux )
3761 deallocate( cs%KPP_temp_flux )
3762 deallocate( cs%KPP_salt_flux )
3763 deallocate( cs%KPP_NLTheat )
3764 deallocate( cs%KPP_NLTscalar )
3765 call kpp_end(cs%KPP_CSp)
3768 if (cs%use_CVMix_conv)
call cvmix_conv_end(cs%CVMix_conv_csp)
3770 if (cs%use_energetic_PBL) &
3771 call energetic_pbl_end(cs%energetic_PBL_CSp)
3772 if (cs%debug_energy_req) &
3773 call diapyc_energy_req_end(cs%diapyc_en_rec_CSp)
3775 if (
associated(cs%optics))
then 3776 call opacity_end(cs%opacity_CSp, cs%optics)
3777 deallocate(cs%optics)
3788 end subroutine diabatic_driver_end
Calculates the heights of sruface or all interfaces from layer thicknesses.
Implemented geothermal heating at the ocean bottom.
Wraps the FMS time manager functions.
This module implements boundary forcing for MOM6.
Interface to CVMix interior shear schemes.
This control structure has parameters for the MOM_internal_tides module.
Provides functions for some diabatic processes such as fraxil, brine rejection, tendency due to surfa...
Functions for calculating interface heights, including free surface height.
Ocean grid type. See mom_grid for details.
The control structure with parameters for the MOM_bulk_mixed_layer module.
Shear-dependent mixing following Jackson et al. 2008.
Calculates density of sea water from T, S and P.
Calculates the freezing point of sea water from T, S and P.
Provides regularization of layers in isopycnal mode.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Vertical viscosities, drag coefficients, and related fields.
Wraps the MPP cpu clock functions.
This control structure contains parameters for MOM_set_diffusivity.
This routine drives the diabatic/dianeutral physics for MOM.
Routines for calculating baroclinic wave speeds.
This control structure holds parameters for the MOM_energetic_PBL module.
The MOM6 facility to parse input files for runtime parameters.
Pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances.
This module contains the routines used to apply sponge layers when using the ALE mode.
ALE sponge control structure.
Orchestrates the registration and calling of tracer packages.
This control structure holds parameters used by the MOM_regularize_layers module. ...
Interface for surface waves.
Make a diagnostic available for averaging or output.
The control structure with paramters for the MOM_opacity module.
Describes various unit conversion factors.
Control structure for diabatic_aux.
Provides routines that do checksums of groups of MOM variables.
Routines used to calculate the opacity of the ocean.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
This control structure holds parameters for the MOM_diapyc_energy_req module.
Describes the decomposed MOM domain and has routines for communications across PEs.
Build mixed layer parameterization.
Calculates the energy requirements of mixing.
Pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances.
Diapycnal mixing and advection in isopycnal mode.
Subroutines that use the ray-tracing equations to propagate the internal tide energy density...
Routines for error handling and I/O management.
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
This control structure holds memory and parameters for the MOM_sponge module.
The control structure for orchestrating the calling of tracer packages.
The control structure holding parametes for the MOM_entrain_diffusive module.
Container for all surface wave related parameters.
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Implements sponge regions in isopycnal mode.
Control structure for containing KPP parameters/data.
Energetically consistent planetary boundary layer parameterization.
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.
Interface to CVMix double diffusion scheme.
Control structure including parameters for CVMix convection.
An overloaded interface to read various types of parameters.
Set up a group of halo updates.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
This type is used to store information about ocean optical properties.
Calculate vertical diffusivity from all mixing processes.
Controls where open boundary conditions are applied.
Provides a transparent vertical ocean grid type and supporting routines.
Interface to CVMix convection scheme.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Control structure for geothermal heating.
A structure for creating arrays of pointers to 3D arrays.
Provides the K-Profile Parameterization (KPP) of Large et al., 1994, via CVMix.
Do a halo update on an array.
Write out checksums of the MOM6 state variables.
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.