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 :: double_diffuse
112 logical :: use_sponge
116 logical :: use_geothermal
117 logical :: use_int_tides
119 logical :: epbl_is_additive
125 real :: uniform_test_cg
127 logical :: usealealgorithm
130 logical :: aggregate_fw_forcing
140 logical :: massless_match_targets
142 logical :: mix_boundary_tracers
153 real :: minimum_forcing_depth
155 real :: evap_cfl_limit = 0.8
157 integer :: halo_ts_diff = 0
159 logical :: usekpp = .false.
160 logical :: kppispassive
162 logical :: debugconservation
163 logical :: tracer_tridiag
165 logical :: debug_energy_req
167 real :: mlddensitydifference
170 real :: mld_en_vals(3)
173 integer :: id_cg1 = -1
174 integer,
allocatable,
dimension(:) :: id_cn
175 integer :: id_wd = -1, id_ea = -1, id_eb = -1
176 integer :: id_dudt_dia = -1, id_dvdt_dia = -1, id_ea_s = -1, id_eb_s = -1
177 integer :: id_ea_t = -1, id_eb_t = -1
178 integer :: id_kd_heat = -1, id_kd_salt = -1, id_kd_interface = -1, id_kd_epbl = -1
179 integer :: id_tdif = -1, id_tadv = -1, id_sdif = -1, id_sadv = -1
180 integer :: id_mld_003 = -1, id_mld_0125 = -1, id_mld_user = -1, id_mlotstsq = -1
181 integer :: id_mld_en1 = -1, id_mld_en2 = -1, id_mld_en3= -1
182 integer :: id_submln2 = -1
183 integer :: id_hf_dudt_dia_2d = -1, id_hf_dvdt_dia_2d = -1
186 integer :: id_u_predia = -1, id_v_predia = -1, id_h_predia = -1
187 integer :: id_t_predia = -1, id_s_predia = -1, id_e_predia = -1
189 integer :: id_diabatic_diff_temp_tend = -1
190 integer :: id_diabatic_diff_saln_tend = -1
191 integer :: id_diabatic_diff_heat_tend = -1
192 integer :: id_diabatic_diff_salt_tend = -1
193 integer :: id_diabatic_diff_heat_tend_2d = -1
194 integer :: id_diabatic_diff_salt_tend_2d = -1
195 integer :: id_diabatic_diff_h= -1
197 integer :: id_boundary_forcing_h = -1
198 integer :: id_boundary_forcing_h_tendency = -1
199 integer :: id_boundary_forcing_temp_tend = -1
200 integer :: id_boundary_forcing_saln_tend = -1
201 integer :: id_boundary_forcing_heat_tend = -1
202 integer :: id_boundary_forcing_salt_tend = -1
203 integer :: id_boundary_forcing_heat_tend_2d = -1
204 integer :: id_boundary_forcing_salt_tend_2d = -1
206 integer :: id_frazil_h = -1
207 integer :: id_frazil_temp_tend = -1
208 integer :: id_frazil_heat_tend = -1
209 integer :: id_frazil_heat_tend_2d = -1
212 logical :: diabatic_diff_tendency_diag = .false.
213 logical :: boundary_forcing_tendency_diag = .false.
214 logical :: frazil_tendency_diag = .false.
215 real,
allocatable,
dimension(:,:,:) :: frazil_heat_diag
216 real,
allocatable,
dimension(:,:,:) :: frazil_temp_diag
233 type(
kpp_cs),
pointer :: kpp_csp => null()
237 type(group_pass_type) :: pass_hold_eb_ea
238 type(group_pass_type) :: pass_kv
241 real,
allocatable,
dimension(:,:,:) :: kpp_nltheat
242 real,
allocatable,
dimension(:,:,:) :: kpp_nltscalar
243 real,
allocatable,
dimension(:,:,:) :: kpp_buoy_flux
244 real,
allocatable,
dimension(:,:) :: kpp_temp_flux
245 real,
allocatable,
dimension(:,:) :: kpp_salt_flux
247 type(time_type),
pointer :: time
251 integer :: id_clock_entrain, id_clock_mixedlayer, id_clock_set_diffusivity
252 integer :: id_clock_tracers, id_clock_tridiag, id_clock_pass, id_clock_sponge
253 integer :: id_clock_geothermal, id_clock_differential_diff, id_clock_remap
254 integer :: id_clock_kpp
261 subroutine diabatic(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
262 G, GV, US, CS, OBC, WAVES)
265 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
266 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
267 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
270 real,
dimension(:,:),
pointer :: hml
271 type(
forcing),
intent(inout) :: fluxes
278 real,
intent(in) :: dt
279 type(time_type),
intent(in) :: time_end
286 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
288 real,
dimension(SZI_(G),SZJ_(G),CS%nMode) :: &
290 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
291 integer :: i, j, k, m, is, ie, js, je, nz
292 logical :: showcalltree
294 if (g%ke == 1)
return
296 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
298 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_diabatic_driver: "// &
299 "Module must be initialized before it is used.")
300 if (dt == 0.0)
call mom_error(fatal,
"MOM_diabatic_driver: "// &
301 "diabatic was called with a zero length timestep.")
302 if (dt < 0.0)
call mom_error(fatal,
"MOM_diabatic_driver: "// &
303 "diabatic was called with a negative timestep.")
305 showcalltree = calltree_showquery()
309 if (cs%id_u_predia > 0)
call post_data(cs%id_u_predia, u, cs%diag)
310 if (cs%id_v_predia > 0)
call post_data(cs%id_v_predia, v, cs%diag)
311 if (cs%id_h_predia > 0)
call post_data(cs%id_h_predia, h, cs%diag)
312 if (cs%id_T_predia > 0)
call post_data(cs%id_T_predia, tv%T, cs%diag)
313 if (cs%id_S_predia > 0)
call post_data(cs%id_S_predia, tv%S, cs%diag)
314 if (cs%id_e_predia > 0)
then
315 call find_eta(h, tv, g, gv, us, eta, eta_to_m=1.0)
316 call post_data(cs%id_e_predia, eta, cs%diag)
320 call mom_state_chksum(
"Start of diabatic ", u, v, h, g, gv, us, haloshift=0)
321 call mom_forcing_chksum(
"Start of diabatic", fluxes, g, us, haloshift=0)
323 if (cs%debugConservation)
call mom_state_stats(
'Start of diabatic', u, v, h, tv%T, tv%S, g, gv, us)
325 if (cs%debug_energy_req) &
326 call diapyc_energy_req_test(h, dt, tv, g, gv, us, cs%diapyc_en_rec_CSp)
328 call cpu_clock_begin(id_clock_set_diffusivity)
329 call set_bbl_tke(u, v, h, fluxes, visc, g, gv, us, cs%set_diff_CSp, obc=obc)
330 call cpu_clock_end(id_clock_set_diffusivity)
335 if (
associated(tv%T) .AND.
associated(tv%frazil))
then
337 call enable_averages(0.5*dt, time_end - real_to_time(0.5*us%T_to_s*dt), cs%diag)
338 if (cs%frazil_tendency_diag)
then
339 do k=1,nz ;
do j=js,je ;
do i=is,ie
340 temp_diag(i,j,k) = tv%T(i,j,k)
341 enddo ;
enddo ;
enddo
344 if (
associated(fluxes%p_surf_full))
then
345 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp, fluxes%p_surf_full, halo=cs%halo_TS_diff)
347 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp, halo=cs%halo_TS_diff)
349 if (showcalltree)
call calltree_waypoint(
"done with 1st make_frazil (diabatic)")
351 if (cs%frazil_tendency_diag)
then
352 call diagnose_frazil_tendency(tv, h, temp_diag, 0.5*dt, g, gv, us, cs)
353 if (cs%id_frazil_h > 0)
call post_data(cs%id_frazil_h, h, cs%diag)
355 call disable_averaging(cs%diag)
357 if (cs%debugConservation)
call mom_state_stats(
'1st make_frazil', u, v, h, tv%T, tv%S, g, gv, us)
359 if (cs%use_int_tides)
then
361 call set_int_tide_input(u, v, h, tv, fluxes, cs%int_tide_input, dt, g, gv, us, &
362 cs%int_tide_input_CSp)
364 if (cs%uniform_test_cg > 0.0)
then
365 do m=1,cs%nMode ; cn_igw(:,:,m) = cs%uniform_test_cg ;
enddo
367 call wave_speeds(h, tv, g, gv, us, cs%nMode, cn_igw, full_halos=.true.)
370 call propagate_int_tide(h, tv, cn_igw, cs%int_tide_input%TKE_itidal_input, cs%int_tide_input%tideamp, &
371 cs%int_tide_input%Nb, dt, g, gv, us, cs%int_tide_CSp)
372 if (showcalltree)
call calltree_waypoint(
"done with propagate_int_tide (diabatic)")
375 if (cs%useALEalgorithm .and. cs%use_legacy_diabatic)
then
376 call diabatic_ale_legacy(u, v, h, tv, hml, fluxes, visc, adp, cdp, dt, time_end, &
377 g, gv, us, cs, waves)
378 elseif (cs%useALEalgorithm)
then
379 call diabatic_ale(u, v, h, tv, hml, fluxes, visc, adp, cdp, dt, time_end, &
380 g, gv, us, cs, waves)
382 call layered_diabatic(u, v, h, tv, hml, fluxes, visc, adp, cdp, dt, time_end, &
383 g, gv, us, cs, waves)
387 call cpu_clock_begin(id_clock_pass)
388 if (
associated(visc%Kv_shear)) &
389 call pass_var(visc%Kv_shear, g%Domain, to_all+omit_corners, halo=1)
390 call cpu_clock_end(id_clock_pass)
392 call disable_averaging(cs%diag)
396 if (
associated(tv%T) .AND.
associated(tv%frazil))
then
397 call enable_averages(0.5*dt, time_end, cs%diag)
398 if (cs%frazil_tendency_diag)
then
399 do k=1,nz ;
do j=js,je ;
do i=is,ie
400 temp_diag(i,j,k) = tv%T(i,j,k)
401 enddo ;
enddo ;
enddo
404 if (
associated(fluxes%p_surf_full))
then
405 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp, fluxes%p_surf_full)
407 call make_frazil(h, tv, g, gv, us, cs%diabatic_aux_CSp)
410 if (cs%frazil_tendency_diag)
then
411 call diagnose_frazil_tendency(tv, h, temp_diag, 0.5*dt, g, gv, us, cs)
412 if (cs%id_frazil_h > 0 )
call post_data(cs%id_frazil_h, h, cs%diag)
415 if (showcalltree)
call calltree_waypoint(
"done with 2nd make_frazil (diabatic)")
416 if (cs%debugConservation)
call mom_state_stats(
'2nd make_frazil', u, v, h, tv%T, tv%S, g, gv, us)
417 call disable_averaging(cs%diag)
423 call enable_averages(dt, time_end, cs%diag)
424 if (cs%id_MLD_003 > 0 .or. cs%id_subMLN2 > 0 .or. cs%id_mlotstsq > 0)
then
425 call diagnosemldbydensitydifference(cs%id_MLD_003, h, tv, 0.03*us%kg_m3_to_R, g, gv, us, cs%diag, &
426 id_n2subml=cs%id_subMLN2, id_mldsq=cs%id_mlotstsq, dz_subml=cs%dz_subML_N2)
428 if (cs%id_MLD_0125 > 0)
then
429 call diagnosemldbydensitydifference(cs%id_MLD_0125, h, tv, 0.125*us%kg_m3_to_R, g, gv, us, cs%diag)
431 if (cs%id_MLD_user > 0)
then
432 call diagnosemldbydensitydifference(cs%id_MLD_user, h, tv, cs%MLDdensityDifference, g, gv, us, cs%diag)
434 if ((cs%id_MLD_EN1 > 0) .or. (cs%id_MLD_EN2 > 0) .or. (cs%id_MLD_EN3 > 0))
then
435 call diagnosemldbyenergy((/cs%id_MLD_EN1, cs%id_MLD_EN2, cs%id_MLD_EN3/),&
436 h, tv, g, gv, us, cs%MLD_EN_VALS, cs%diag)
438 if (cs%use_int_tides)
then
439 if (cs%id_cg1 > 0)
call post_data(cs%id_cg1, cn_igw(:,:,1),cs%diag)
440 do m=1,cs%nMode ;
if (cs%id_cn(m) > 0)
call post_data(cs%id_cn(m), cn_igw(:,:,m), cs%diag) ;
enddo
442 call disable_averaging(cs%diag)
444 if (cs%debugConservation)
call mom_state_stats(
'leaving diabatic', u, v, h, tv%T, tv%S, g, gv, us)
446 end subroutine diabatic
451 subroutine diabatic_ale_legacy(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
452 G, GV, US, CS, Waves)
456 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
457 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
458 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
461 real,
dimension(:,:),
pointer :: Hml
462 type(
forcing),
intent(inout) :: fluxes
469 real,
intent(in) :: dt
470 type(time_type),
intent(in) :: Time_end
475 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
494 real,
dimension(SZI_(G),SZJ_(G)) :: &
496 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: h_diag
497 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
498 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: saln_diag
499 real,
dimension(SZI_(G),SZJ_(G)) :: tendency_2d
503 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
508 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
522 real,
allocatable,
dimension(:,:) :: &
523 hf_dudt_dia_2d, hf_dvdt_dia_2d
525 logical :: in_boundary(SZI_(G))
545 real :: htot(SZIB_(G))
546 real :: b1(SZIB_(G)), d1(SZIB_(G))
547 real :: c1(SZIB_(G),SZK_(G))
553 logical :: showCallTree
554 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, m
559 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
560 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
561 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect*h_neglect
562 kd_heat(:,:,:) = 0.0 ; kd_salt(:,:,:) = 0.0
564 showcalltree = calltree_showquery()
565 if (showcalltree)
call calltree_enter(
"diabatic_ALE_legacy(), MOM_diabatic_driver.F90")
568 call enable_averages(dt, time_end, cs%diag)
570 if (cs%use_geothermal)
then
571 call cpu_clock_begin(id_clock_geothermal)
572 call geothermal_in_place(h, tv, dt, g, gv, us, cs%geothermal_CSp, halo=cs%halo_TS_diff)
573 call cpu_clock_end(id_clock_geothermal)
574 if (showcalltree)
call calltree_waypoint(
"geothermal (diabatic)")
575 if (cs%debugConservation)
call mom_state_stats(
'geothermal', u, v, h, tv%T, tv%S, g, gv, us)
580 call diag_update_remap_grids(cs%diag)
585 if (
associated(cs%optics)) &
586 call set_pen_shortwave(cs%optics, fluxes, g, gv, us, cs%diabatic_aux_CSp, cs%opacity_CSp, cs%tracer_flow_CSp)
588 if (cs%debug)
call mom_state_chksum(
"before find_uv_at_h", u, v, h, g, gv, us, haloshift=0)
590 if (cs%use_kappa_shear .or. cs%use_CVMix_shear)
then
591 if (cs%use_geothermal)
then
594 eatr(:,:,:) = 0.0 ; ebtr(:,:,:) = 0.0
595 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us, eatr, ebtr)
597 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
599 if (showcalltree)
call calltree_waypoint(
"done with find_uv_at_h (diabatic)")
602 call cpu_clock_begin(id_clock_set_diffusivity)
606 call mom_state_chksum(
"before set_diffusivity", u, v, h, g, gv, us, haloshift=cs%halo_TS_diff)
607 if (cs%double_diffuse)
then
608 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, visc, dt, g, gv, us, cs%set_diff_CSp, &
609 kd_int=kd_int, kd_extra_t=kd_extra_t, kd_extra_s=kd_extra_s)
611 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, visc, dt, g, gv, us, &
612 cs%set_diff_CSp, kd_int=kd_int)
614 call cpu_clock_end(id_clock_set_diffusivity)
615 if (showcalltree)
call calltree_waypoint(
"done with set_diffusivity (diabatic)")
621 if (cs%double_diffuse)
then
623 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
624 kd_salt(i,j,k) = kd_int(i,j,k) + kd_extra_s(i,j,k)
625 kd_heat(i,j,k) = kd_int(i,j,k) + kd_extra_t(i,j,k)
626 enddo ;
enddo ;
enddo
629 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
630 kd_salt(i,j,k) = kd_int(i,j,k)
631 kd_heat(i,j,k) = kd_int(i,j,k)
632 enddo ;
enddo ;
enddo
637 call mom_state_chksum(
"after set_diffusivity ", u, v, h, g, gv, us, haloshift=0)
638 call mom_forcing_chksum(
"after set_diffusivity ", fluxes, g, us, haloshift=0)
639 call mom_thermovar_chksum(
"after set_diffusivity ", tv, g)
640 call hchksum(kd_int,
"after set_diffusivity Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
641 call hchksum(kd_heat,
"after set_diffusivity Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
642 call hchksum(kd_salt,
"after set_diffusivity Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
646 call cpu_clock_begin(id_clock_kpp)
654 call calculatebuoyancyflux2d(g, gv, us, fluxes, cs%optics, h, tv%T, tv%S, tv, &
655 cs%KPP_buoy_flux, cs%KPP_temp_flux, cs%KPP_salt_flux)
658 call kpp_compute_bld(cs%KPP_CSp, g, gv, us, h, tv%T, tv%S, u, v, tv, &
659 fluxes%ustar, cs%KPP_buoy_flux, waves=waves)
661 call kpp_calculate(cs%KPP_CSp, g, gv, us, h, fluxes%ustar, cs%KPP_buoy_flux, kd_heat, &
662 kd_salt, visc%Kv_shear, cs%KPP_NLTheat, cs%KPP_NLTscalar, waves=waves)
664 if (
associated(hml))
then
665 call kpp_get_bld(cs%KPP_CSp, hml(:,:), g, us)
666 call pass_var(hml, g%domain, halo=1)
668 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
671 if (.not.cs%KPPisPassive)
then
673 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
674 kd_int(i,j,k) = min( kd_salt(i,j,k), kd_heat(i,j,k) )
675 enddo ;
enddo ;
enddo
676 if (cs%double_diffuse)
then
678 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
679 kd_extra_s(i,j,k) = (kd_salt(i,j,k) - kd_int(i,j,k))
680 kd_extra_t(i,j,k) = (kd_heat(i,j,k) - kd_int(i,j,k))
681 enddo ;
enddo ;
enddo
685 call cpu_clock_end(id_clock_kpp)
686 if (showcalltree)
call calltree_waypoint(
"done with KPP_calculate (diabatic)")
689 call mom_forcing_chksum(
"after KPP", fluxes, g, us, haloshift=0)
690 call mom_thermovar_chksum(
"after KPP", tv, g)
691 call hchksum(kd_heat,
"after KPP Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
692 call hchksum(kd_salt,
"after KPP Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
699 call cpu_clock_begin(id_clock_kpp)
701 call hchksum(cs%KPP_temp_flux,
"before KPP_applyNLT netHeat", g%HI, haloshift=0, scale=gv%H_to_m)
702 call hchksum(cs%KPP_salt_flux,
"before KPP_applyNLT netSalt", g%HI, haloshift=0, scale=gv%H_to_m)
703 call hchksum(cs%KPP_NLTheat,
"before KPP_applyNLT NLTheat", g%HI, haloshift=0)
704 call hchksum(cs%KPP_NLTscalar,
"before KPP_applyNLT NLTscalar", g%HI, haloshift=0)
708 call kpp_nonlocaltransport_temp(cs%KPP_CSp, g, gv, h, cs%KPP_NLTheat, cs%KPP_temp_flux, &
709 us%T_to_s*dt, tv%T, us%Q_to_J_kg*tv%C_p)
710 call kpp_nonlocaltransport_saln(cs%KPP_CSp, g, gv, h, cs%KPP_NLTscalar, cs%KPP_salt_flux, &
712 call cpu_clock_end(id_clock_kpp)
713 if (showcalltree)
call calltree_waypoint(
"done with KPP_applyNonLocalTransport (diabatic)")
714 if (cs%debugConservation)
call mom_state_stats(
'KPP_applyNonLocalTransport', u, v, h, tv%T, tv%S, g, gv, us)
717 call mom_state_chksum(
"after KPP_applyNLT ", u, v, h, g, gv, us, haloshift=0)
718 call mom_forcing_chksum(
"after KPP_applyNLT ", fluxes, g, us, haloshift=0)
719 call mom_thermovar_chksum(
"after KPP_applyNLT ", tv, g)
725 if (cs%double_diffuse .and.
associated(tv%T))
then
727 call cpu_clock_begin(id_clock_differential_diff)
728 call differential_diffuse_t_s(h, tv%T, tv%S, kd_extra_t, kd_extra_s, dt, g, gv)
729 call cpu_clock_end(id_clock_differential_diff)
731 if (showcalltree)
call calltree_waypoint(
"done with differential_diffuse_T_S (diabatic)")
732 if (cs%debugConservation)
call mom_state_stats(
'differential_diffuse_T_S', u, v, h, tv%T, tv%S, g, gv, us)
736 if (.not. cs%useKPP)
then
738 do k=2,nz ;
do j=js,je ;
do i=is,ie
739 kd_heat(i,j,k) = kd_heat(i,j,k) + kd_extra_t(i,j,k)
740 kd_salt(i,j,k) = kd_salt(i,j,k) + kd_extra_s(i,j,k)
741 enddo ;
enddo ;
enddo
747 if (cs%use_CVMix_conv)
then
749 call calculate_cvmix_conv(h, tv, g, gv, us, cs%CVMix_conv_csp, hml, kd=kd_int, kv=visc%Kv_slow)
753 do j=js,je ;
do i=is,ie
757 do k=2,nz ;
do j=js,je ;
do i=is,ie
758 hval=1.0/(h_neglect + 0.5*(h(i,j,k-1) + h(i,j,k)))
759 ea_s(i,j,k) = (gv%Z_to_H**2) * dt * hval * kd_int(i,j,k)
760 eb_s(i,j,k-1) = ea_s(i,j,k)
761 ea_t(i,j,k-1) = ea_s(i,j,k-1) ; eb_t(i,j,k-1) = eb_s(i,j,k-1)
762 enddo ;
enddo ;
enddo
763 do j=js,je ;
do i=is,ie
765 ea_t(i,j,nz) = ea_s(i,j,nz) ; eb_t(i,j,nz) = eb_s(i,j,nz)
767 if (showcalltree)
call calltree_waypoint(
"done setting ea,eb from Kd_int (diabatic)")
770 call mom_forcing_chksum(
"after calc_entrain ", fluxes, g, us, haloshift=0)
771 call mom_thermovar_chksum(
"after calc_entrain ", tv, g)
772 call mom_state_chksum(
"after calc_entrain ", u, v, h, g, gv, us, haloshift=0)
773 call hchksum(ea_s,
"after calc_entrain ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
774 call hchksum(eb_s,
"after calc_entrain eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
778 if (cs%boundary_forcing_tendency_diag)
then
779 do k=1,nz ;
do j=js,je ;
do i=is,ie
780 h_diag(i,j,k) = h(i,j,k)
781 temp_diag(i,j,k) = tv%T(i,j,k)
782 saln_diag(i,j,k) = tv%S(i,j,k)
783 enddo ;
enddo ;
enddo
787 call cpu_clock_begin(id_clock_remap)
790 do k=1,nz ;
do j=js,je ;
do i=is,ie
791 h_prebound(i,j,k) = h(i,j,k)
792 enddo ;
enddo ;
enddo
793 if (cs%use_energetic_PBL)
then
795 skinbuoyflux(:,:) = 0.0
796 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
797 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, cs%evap_CFL_limit, &
798 cs%minimum_forcing_depth, ctke, dsv_dt, dsv_ds, skinbuoyflux=skinbuoyflux)
801 call hchksum(ea_t,
"after applyBoundaryFluxes ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
802 call hchksum(eb_t,
"after applyBoundaryFluxes eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
803 call hchksum(ea_s,
"after applyBoundaryFluxes ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
804 call hchksum(eb_s,
"after applyBoundaryFluxes eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
805 call hchksum(ctke,
"after applyBoundaryFluxes cTKE", g%HI, haloshift=0, &
806 scale=us%RZ3_T3_to_W_m2*us%T_to_s)
807 call hchksum(dsv_dt,
"after applyBoundaryFluxes dSV_dT", g%HI, haloshift=0, scale=us%kg_m3_to_R)
808 call hchksum(dsv_ds,
"after applyBoundaryFluxes dSV_dS", g%HI, haloshift=0, scale=us%kg_m3_to_R)
811 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
812 call energetic_pbl(h, u_h, v_h, tv, fluxes, dt, kd_epbl, g, gv, us, &
813 cs%energetic_PBL_CSp, dsv_dt, dsv_ds, ctke, skinbuoyflux, waves=waves)
815 if (
associated(hml))
then
816 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, hml(:,:), g, us)
817 call pass_var(hml, g%domain, halo=1)
819 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
820 elseif (
associated(visc%MLD))
then
821 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, visc%MLD, g, us)
822 call pass_var(visc%MLD, g%domain, halo=1)
826 do k=2,nz ;
do j=js,je ;
do i=is,ie
827 if (cs%ePBL_is_additive)
then
828 kd_add_here = kd_epbl(i,j,k)
829 visc%Kv_shear(i,j,k) = visc%Kv_shear(i,j,k) + cs%ePBL_Prandtl*kd_epbl(i,j,k)
831 kd_add_here = max(kd_epbl(i,j,k) - visc%Kd_shear(i,j,k), 0.0)
832 visc%Kv_shear(i,j,k) = max(visc%Kv_shear(i,j,k), cs%ePBL_Prandtl*kd_epbl(i,j,k))
835 ent_int = kd_add_here * (gv%Z_to_H**2 * dt) / &
836 (0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect)
837 eb_s(i,j,k-1) = eb_s(i,j,k-1) + ent_int
838 ea_s(i,j,k) = ea_s(i,j,k) + ent_int
839 kd_int(i,j,k) = kd_int(i,j,k) + kd_add_here
842 kd_heat(i,j,k) = kd_heat(i,j,k) + kd_int(i,j,k)
843 kd_salt(i,j,k) = kd_salt(i,j,k) + kd_int(i,j,k)
845 enddo ;
enddo ;
enddo
848 call hchksum(ea_t,
"after ePBL ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
849 call hchksum(eb_t,
"after ePBL eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
850 call hchksum(ea_s,
"after ePBL ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
851 call hchksum(eb_s,
"after ePBL eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
852 call hchksum(kd_epbl,
"after ePBL Kd_ePBL", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
856 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
857 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, &
858 cs%evap_CFL_limit, cs%minimum_forcing_depth)
865 if (cs%boundary_forcing_tendency_diag)
then
866 call diagnose_boundary_forcing_tendency(tv, h, temp_diag, saln_diag, h_diag, dt, g, gv, us, cs)
867 if (cs%id_boundary_forcing_h > 0)
call post_data(cs%id_boundary_forcing_h, h, cs%diag, alt_h = h_diag)
870 call diag_update_remap_grids(cs%diag)
871 call cpu_clock_end(id_clock_remap)
873 call mom_forcing_chksum(
"after applyBoundaryFluxes ", fluxes, g, us, haloshift=0)
874 call mom_thermovar_chksum(
"after applyBoundaryFluxes ", tv, g)
875 call mom_state_chksum(
"after applyBoundaryFluxes ", u, v, h, g, gv, us, haloshift=0)
877 if (showcalltree)
call calltree_waypoint(
"done with applyBoundaryFluxes (diabatic)")
878 if (cs%debugConservation)
call mom_state_stats(
'applyBoundaryFluxes', u, v, h, tv%T, tv%S, g, gv, us)
889 hold(i,j,1) = h(i,j,1)
891 hold(i,j,nz) = h(i,j,nz)
893 if (h(i,j,1) <= 0.0) h(i,j,1) = gv%Angstrom_H
894 if (h(i,j,nz) <= 0.0) h(i,j,nz) = gv%Angstrom_H
896 do k=2,nz-1 ;
do i=is,ie
897 hold(i,j,k) = h(i,j,k)
900 if (h(i,j,k) <= 0.0) h(i,j,k) = gv%Angstrom_H
904 call diag_update_remap_grids(cs%diag)
907 call mom_state_chksum(
"after negative check ", u, v, h, g, gv, us, haloshift=0)
908 call mom_forcing_chksum(
"after negative check ", fluxes, g, us, haloshift=0)
909 call mom_thermovar_chksum(
"after negative check ", tv, g)
911 if (showcalltree)
call calltree_waypoint(
"done with h=ea-eb (diabatic)")
912 if (cs%debugConservation)
call mom_state_stats(
'h=ea-eb', u, v, h, tv%T, tv%S, g, gv, us)
915 if (
associated(tv%T))
then
918 call hchksum(ea_t,
"before triDiagTS ea_t ", g%HI, haloshift=0, scale=gv%H_to_m)
919 call hchksum(eb_t,
"before triDiagTS eb_t ", g%HI, haloshift=0, scale=gv%H_to_m)
920 call hchksum(ea_s,
"before triDiagTS ea_s ", g%HI, haloshift=0, scale=gv%H_to_m)
921 call hchksum(eb_s,
"before triDiagTS eb_s ", g%HI, haloshift=0, scale=gv%H_to_m)
924 call cpu_clock_begin(id_clock_tridiag)
931 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
932 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
934 if (cs%diabatic_diff_tendency_diag)
then
935 do k=1,nz ;
do j=js,je ;
do i=is,ie
936 temp_diag(i,j,k) = tv%T(i,j,k)
937 saln_diag(i,j,k) = tv%S(i,j,k)
938 enddo ;
enddo ;
enddo
942 do k=1,nz ;
do j=js,je ;
do i=is,ie
943 ea_t(i,j,k) = ea_s(i,j,k) ; eb_t(i,j,k) = eb_s(i,j,k)
944 enddo ;
enddo ;
enddo
945 if (cs%tracer_tridiag)
then
946 call tracer_vertdiff(hold, ea_t, eb_t, dt, tv%T, g, gv)
947 call tracer_vertdiff(hold, ea_s, eb_s, dt, tv%S, g, gv)
949 call tridiagts(g, gv, is, ie, js, je, hold, ea_s, eb_s, tv%T, tv%S)
956 if (cs%diabatic_diff_tendency_diag)
then
957 call diagnose_diabatic_diff_tendency(tv, hold, temp_diag, saln_diag, dt, g, gv, us, cs)
958 if (cs%id_diabatic_diff_h > 0)
call post_data(cs%id_diabatic_diff_h, hold, cs%diag, alt_h = hold)
961 call cpu_clock_end(id_clock_tridiag)
963 if (showcalltree)
call calltree_waypoint(
"done with triDiagTS (diabatic)")
967 if (cs%debugConservation)
call mom_state_stats(
'triDiagTS', u, v, h, tv%T, tv%S, g, gv, us)
970 call mom_state_chksum(
"after mixed layer ", u, v, h, g, gv, us, haloshift=0)
971 call mom_thermovar_chksum(
"after mixed layer ", tv, g)
976 if (
associated(adp%du_dt_dia) .or.
associated(adp%dv_dt_dia)) &
978 call pass_var(h, g%domain, to_west+to_south+omit_corners, halo=1)
979 call diag_update_remap_grids(cs%diag)
983 if ((cs%id_Tdif > 0) .or. (cs%id_Tadv > 0))
then
984 do j=js,je ;
do i=is,ie
985 tdif_flx(i,j,1) = 0.0 ; tdif_flx(i,j,nz+1) = 0.0
986 tadv_flx(i,j,1) = 0.0 ; tadv_flx(i,j,nz+1) = 0.0
989 do k=2,nz ;
do j=js,je ;
do i=is,ie
990 tdif_flx(i,j,k) = (idt * 0.5*(ea_t(i,j,k) + eb_t(i,j,k-1))) * &
991 (tv%T(i,j,k-1) - tv%T(i,j,k))
992 tadv_flx(i,j,k) = (idt * (ea_t(i,j,k) - eb_t(i,j,k-1))) * &
993 0.5*(tv%T(i,j,k-1) + tv%T(i,j,k))
994 enddo ;
enddo ;
enddo
996 if ((cs%id_Sdif > 0) .or. (cs%id_Sadv > 0))
then
997 do j=js,je ;
do i=is,ie
998 sdif_flx(i,j,1) = 0.0 ; sdif_flx(i,j,nz+1) = 0.0
999 sadv_flx(i,j,1) = 0.0 ; sadv_flx(i,j,nz+1) = 0.0
1002 do k=2,nz ;
do j=js,je ;
do i=is,ie
1003 sdif_flx(i,j,k) = (idt * 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))) * &
1004 (tv%S(i,j,k-1) - tv%S(i,j,k))
1005 sadv_flx(i,j,k) = (idt * (ea_s(i,j,k) - eb_s(i,j,k-1))) * &
1006 0.5*(tv%S(i,j,k-1) + tv%S(i,j,k))
1007 enddo ;
enddo ;
enddo
1011 call cpu_clock_begin(id_clock_tracers)
1013 if (cs%mix_boundary_tracers)
then
1014 tr_ea_bbl = gv%Z_to_H * sqrt(dt*cs%Kd_BBL_tr)
1018 ebtr(i,j,nz) = eb_s(i,j,nz)
1020 in_boundary(i) = (g%mask2dT(i,j) > 0.0)
1022 do k=nz,2,-1 ;
do i=is,ie
1023 if (in_boundary(i))
then
1024 htot(i) = htot(i) + h(i,j,k)
1033 add_ent = ((dt * cs%Kd_min_tr) * gv%Z_to_H**2) * &
1034 ((h(i,j,k-1)+h(i,j,k)+h_neglect) / &
1035 (h(i,j,k-1)*h(i,j,k)+h_neglect2)) - &
1036 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))
1037 if (htot(i) < tr_ea_bbl)
then
1038 add_ent = max(0.0, add_ent, &
1039 (tr_ea_bbl - htot(i)) - min(ea_s(i,j,k),eb_s(i,j,k-1)))
1040 elseif (add_ent < 0.0)
then
1041 add_ent = 0.0 ; in_boundary(i) = .false.
1044 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1045 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1047 ebtr(i,j,k-1) = eb_s(i,j,k-1) ; eatr(i,j,k) = ea_s(i,j,k)
1050 if (cs%double_diffuse)
then ;
if (kd_extra_s(i,j,k) > 0.0)
then
1051 add_ent = ((dt * kd_extra_s(i,j,k)) * gv%Z_to_H**2) / &
1052 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + h_neglect)
1053 ebtr(i,j,k-1) = ebtr(i,j,k-1) + add_ent
1054 eatr(i,j,k) = eatr(i,j,k) + add_ent
1057 do i=is,ie ; eatr(i,j,1) = ea_s(i,j,1) ;
enddo
1064 call call_tracer_column_fns(h_prebound, h, ea_s, eb_s, fluxes, hml, dt, g, gv, us, tv, &
1065 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1066 evap_cfl_limit = cs%evap_CFL_limit, &
1067 minimum_forcing_depth=cs%minimum_forcing_depth)
1069 elseif (cs%double_diffuse)
then
1071 do j=js,je ;
do i=is,ie
1072 ebtr(i,j,nz) = eb_s(i,j,nz) ; eatr(i,j,1) = ea_s(i,j,1)
1075 do k=nz,2,-1 ;
do j=js,je ;
do i=is,ie
1076 if (kd_extra_s(i,j,k) > 0.0)
then
1077 add_ent = ((dt * kd_extra_s(i,j,k)) * gv%Z_to_H**2) / &
1078 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + h_neglect)
1082 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1083 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1084 enddo ;
enddo ;
enddo
1087 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1088 cs%optics, cs%tracer_flow_CSp, cs%debug,&
1089 evap_cfl_limit = cs%evap_CFL_limit, &
1090 minimum_forcing_depth=cs%minimum_forcing_depth)
1094 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1095 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1096 evap_cfl_limit = cs%evap_CFL_limit, &
1097 minimum_forcing_depth=cs%minimum_forcing_depth)
1100 call cpu_clock_end(id_clock_tracers)
1103 if (cs%use_sponge)
then
1104 call cpu_clock_begin(id_clock_sponge)
1105 if (
associated(cs%ALE_sponge_CSp))
then
1106 call apply_ale_sponge(h, dt, g, gv, us, cs%ALE_sponge_CSp, cs%Time)
1109 call cpu_clock_end(id_clock_sponge)
1112 call mom_thermovar_chksum(
"apply_sponge ", tv, g)
1116 call disable_averaging(cs%diag)
1118 call enable_averages(dt, time_end, cs%diag)
1120 if (cs%id_Kd_interface > 0)
call post_data(cs%id_Kd_interface, kd_int, cs%diag)
1121 if (cs%id_Kd_heat > 0)
call post_data(cs%id_Kd_heat, kd_heat, cs%diag)
1122 if (cs%id_Kd_salt > 0)
call post_data(cs%id_Kd_salt, kd_salt, cs%diag)
1123 if (cs%id_Kd_ePBL > 0)
call post_data(cs%id_Kd_ePBL, kd_epbl, cs%diag)
1125 if (cs%id_ea > 0)
call post_data(cs%id_ea, ea_s, cs%diag)
1126 if (cs%id_eb > 0)
call post_data(cs%id_eb, eb_s, cs%diag)
1127 if (cs%id_ea_t > 0)
call post_data(cs%id_ea_t, ea_t, cs%diag)
1128 if (cs%id_eb_t > 0)
call post_data(cs%id_eb_t, eb_t, cs%diag)
1129 if (cs%id_ea_s > 0)
call post_data(cs%id_ea_s, ea_s, cs%diag)
1130 if (cs%id_eb_s > 0)
call post_data(cs%id_eb_s, eb_s, cs%diag)
1132 if (cs%id_dudt_dia > 0)
call post_data(cs%id_dudt_dia, adp%du_dt_dia, cs%diag)
1133 if (cs%id_dvdt_dia > 0)
call post_data(cs%id_dvdt_dia, adp%dv_dt_dia, cs%diag)
1134 if (cs%id_wd > 0)
call post_data(cs%id_wd, cdp%diapyc_vel, cs%diag)
1136 if (cs%id_Tdif > 0)
call post_data(cs%id_Tdif, tdif_flx, cs%diag)
1137 if (cs%id_Tadv > 0)
call post_data(cs%id_Tadv, tadv_flx, cs%diag)
1138 if (cs%id_Sdif > 0)
call post_data(cs%id_Sdif, sdif_flx, cs%diag)
1139 if (cs%id_Sadv > 0)
call post_data(cs%id_Sadv, sadv_flx, cs%diag)
1142 if (cs%id_hf_dudt_dia_2d > 0)
then
1143 allocate(hf_dudt_dia_2d(g%IsdB:g%IedB,g%jsd:g%jed))
1144 hf_dudt_dia_2d(:,:) = 0.0
1145 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1146 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)
1147 enddo ;
enddo ;
enddo
1148 call post_data(cs%id_hf_dudt_dia_2d, hf_dudt_dia_2d, cs%diag)
1149 deallocate(hf_dudt_dia_2d)
1152 if (cs%id_hf_dvdt_dia_2d > 0)
then
1153 allocate(hf_dvdt_dia_2d(g%isd:g%ied,g%JsdB:g%JedB))
1154 hf_dvdt_dia_2d(:,:) = 0.0
1155 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1156 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)
1157 enddo ;
enddo ;
enddo
1158 call post_data(cs%id_hf_dvdt_dia_2d, hf_dvdt_dia_2d, cs%diag)
1159 deallocate(hf_dvdt_dia_2d)
1162 call disable_averaging(cs%diag)
1164 if (showcalltree)
call calltree_leave(
"diabatic_ALE_legacy()")
1166 end subroutine diabatic_ale_legacy
1171 subroutine diabatic_ale(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
1172 G, GV, US, CS, Waves)
1176 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
1177 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
1178 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
1181 real,
dimension(:,:),
pointer :: Hml
1182 type(
forcing),
intent(inout) :: fluxes
1189 real,
intent(in) :: dt
1190 type(time_type),
intent(in) :: Time_end
1195 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
1210 real,
dimension(SZI_(G),SZJ_(G)) :: &
1212 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: h_diag
1213 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
1214 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: saln_diag
1215 real,
dimension(SZI_(G),SZJ_(G)) :: tendency_2d
1219 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
1224 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
1238 real,
allocatable,
dimension(:,:) :: &
1239 hf_dudt_dia_2d, hf_dvdt_dia_2d
1241 logical :: in_boundary(SZI_(G))
1261 real :: htot(SZIB_(G))
1262 real :: b1(SZIB_(G)), d1(SZIB_(G))
1263 real :: c1(SZIB_(G),SZK_(G))
1269 logical :: showCallTree
1270 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, m
1272 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1273 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1274 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect*h_neglect
1275 kd_heat(:,:,:) = 0.0 ; kd_salt(:,:,:) = 0.0
1276 ea_s(:,:,:) = 0.0; eb_s(:,:,:) = 0.0; ea_t(:,:,:) = 0.0; eb_t(:,:,:) = 0.0
1278 showcalltree = calltree_showquery()
1279 if (showcalltree)
call calltree_enter(
"diabatic_ALE(), MOM_diabatic_driver.F90")
1281 if (.not. (cs%useALEalgorithm))
call mom_error(fatal,
"MOM_diabatic_driver: "// &
1282 "The ALE algorithm must be enabled when using MOM_diabatic_driver.")
1285 call enable_averages(dt, time_end, cs%diag)
1287 if (cs%use_geothermal)
then
1288 call cpu_clock_begin(id_clock_geothermal)
1289 call geothermal_in_place(h, tv, dt, g, gv, us, cs%geothermal_CSp, halo=cs%halo_TS_diff)
1290 call cpu_clock_end(id_clock_geothermal)
1291 if (showcalltree)
call calltree_waypoint(
"geothermal (diabatic)")
1292 if (cs%debugConservation)
call mom_state_stats(
'geothermal', u, v, h, tv%T, tv%S, g, gv, us)
1297 call diag_update_remap_grids(cs%diag)
1302 if (
associated(cs%optics)) &
1303 call set_pen_shortwave(cs%optics, fluxes, g, gv, us, cs%diabatic_aux_CSp, cs%opacity_CSp, cs%tracer_flow_CSp)
1305 if (cs%debug)
call mom_state_chksum(
"before find_uv_at_h", u, v, h, g, gv, us, haloshift=0)
1307 if (cs%use_kappa_shear .or. cs%use_CVMix_shear)
then
1308 if (cs%use_geothermal)
then
1311 eatr(:,:,:) = 0.0 ; ebtr(:,:,:) = 0.0
1312 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us, eatr, ebtr)
1314 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
1316 if (showcalltree)
call calltree_waypoint(
"done with find_uv_at_h (diabatic)")
1319 call cpu_clock_begin(id_clock_set_diffusivity)
1323 call mom_state_chksum(
"before set_diffusivity", u, v, h, g, gv, us, haloshift=cs%halo_TS_diff)
1324 if (cs%double_diffuse)
then
1325 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, visc, dt, g, gv, us, cs%set_diff_CSp, &
1326 kd_int=kd_int, kd_extra_t=kd_extra_t, kd_extra_s=kd_extra_s)
1328 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, visc, dt, g, gv, us, &
1329 cs%set_diff_CSp, kd_int=kd_int)
1331 call cpu_clock_end(id_clock_set_diffusivity)
1332 if (showcalltree)
call calltree_waypoint(
"done with set_diffusivity (diabatic)")
1335 call mom_state_chksum(
"after set_diffusivity ", u, v, h, g, gv, us, haloshift=0)
1336 call mom_forcing_chksum(
"after set_diffusivity ", fluxes, g, us, haloshift=0)
1337 call mom_thermovar_chksum(
"after set_diffusivity ", tv, g)
1338 call hchksum(kd_int,
"after set_diffusivity Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1344 if (cs%double_diffuse)
then
1346 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1347 kd_salt(i,j,k) = kd_int(i,j,k) + kd_extra_s(i,j,k)
1348 kd_heat(i,j,k) = kd_int(i,j,k) + kd_extra_t(i,j,k)
1349 enddo ;
enddo ;
enddo
1352 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1353 kd_salt(i,j,k) = kd_int(i,j,k)
1354 kd_heat(i,j,k) = kd_int(i,j,k)
1355 enddo ;
enddo ;
enddo
1359 call hchksum(kd_heat,
"after set_diffusivity Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1360 call hchksum(kd_salt,
"after set_diffusivity Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1364 call cpu_clock_begin(id_clock_kpp)
1366 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
1367 visc%Kv_shear(i,j,k) = visc%Kv_shear(i,j,k) + visc%Kv_slow(i,j,k)
1368 enddo ;
enddo ;
enddo
1375 call calculatebuoyancyflux2d(g, gv, us, fluxes, cs%optics, h, tv%T, tv%S, tv, &
1376 cs%KPP_buoy_flux, cs%KPP_temp_flux, cs%KPP_salt_flux)
1379 call kpp_compute_bld(cs%KPP_CSp, g, gv, us, h, tv%T, tv%S, u, v, tv, &
1380 fluxes%ustar, cs%KPP_buoy_flux, waves=waves)
1382 call kpp_calculate(cs%KPP_CSp, g, gv, us, h, fluxes%ustar, cs%KPP_buoy_flux, kd_heat, &
1383 kd_salt, visc%Kv_shear, cs%KPP_NLTheat, cs%KPP_NLTscalar, waves=waves)
1385 if (
associated(hml))
then
1386 call kpp_get_bld(cs%KPP_CSp, hml(:,:), g, us)
1387 call pass_var(hml, g%domain, halo=1)
1389 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
1392 call cpu_clock_end(id_clock_kpp)
1393 if (showcalltree)
call calltree_waypoint(
"done with KPP_calculate (diabatic)")
1396 call mom_forcing_chksum(
"after KPP", fluxes, g, us, haloshift=0)
1397 call mom_thermovar_chksum(
"after KPP", tv, g)
1398 call hchksum(kd_heat,
"after KPP Kd_heat", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1399 call hchksum(kd_salt,
"after KPP Kd_salt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1406 call cpu_clock_begin(id_clock_kpp)
1408 call hchksum(cs%KPP_temp_flux,
"before KPP_applyNLT netHeat", g%HI, haloshift=0, scale=gv%H_to_m)
1409 call hchksum(cs%KPP_salt_flux,
"before KPP_applyNLT netSalt", g%HI, haloshift=0, scale=gv%H_to_m)
1410 call hchksum(cs%KPP_NLTheat,
"before KPP_applyNLT NLTheat", g%HI, haloshift=0)
1411 call hchksum(cs%KPP_NLTscalar,
"before KPP_applyNLT NLTscalar", g%HI, haloshift=0)
1415 call kpp_nonlocaltransport_temp(cs%KPP_CSp, g, gv, h, cs%KPP_NLTheat, cs%KPP_temp_flux, &
1416 us%T_to_s*dt, tv%T, us%Q_to_J_kg*tv%C_p)
1417 call kpp_nonlocaltransport_saln(cs%KPP_CSp, g, gv, h, cs%KPP_NLTscalar, cs%KPP_salt_flux, &
1419 call cpu_clock_end(id_clock_kpp)
1420 if (showcalltree)
call calltree_waypoint(
"done with KPP_applyNonLocalTransport (diabatic)")
1421 if (cs%debugConservation)
call mom_state_stats(
'KPP_applyNonLocalTransport', u, v, h, tv%T, tv%S, g, gv, us)
1424 call mom_state_chksum(
"after KPP_applyNLT ", u, v, h, g, gv, us, haloshift=0)
1425 call mom_forcing_chksum(
"after KPP_applyNLT ", fluxes, g, us, haloshift=0)
1426 call mom_thermovar_chksum(
"after KPP_applyNLT ", tv, g)
1432 if (cs%double_diffuse .and.
associated(tv%T) .and. (.not.cs%use_CVMix_ddiff))
then
1434 call cpu_clock_begin(id_clock_differential_diff)
1435 call differential_diffuse_t_s(h, tv%T, tv%S, kd_extra_t, kd_extra_s, dt, g, gv)
1436 call cpu_clock_end(id_clock_differential_diff)
1438 if (showcalltree)
call calltree_waypoint(
"done with differential_diffuse_T_S (diabatic)")
1439 if (cs%debugConservation)
call mom_state_stats(
'differential_diffuse_T_S', u, v, h, tv%T, tv%S, g, gv, us)
1443 if (.not. cs%useKPP)
then
1445 do k=2,nz ;
do j=js,je ;
do i=is,ie
1446 kd_heat(i,j,k) = kd_heat(i,j,k) + kd_extra_t(i,j,k)
1447 kd_salt(i,j,k) = kd_salt(i,j,k) + kd_extra_s(i,j,k)
1448 enddo ;
enddo ;
enddo
1454 if (cs%use_CVMix_conv)
then
1457 call calculate_cvmix_conv(h, tv, g, gv, us, cs%CVMix_conv_csp, hml, kd=kd_heat, kv=visc%Kv_shear, kd_aux=kd_salt)
1459 call calculate_cvmix_conv(h, tv, g, gv, us, cs%CVMix_conv_csp, hml, kd=kd_heat, kv=visc%Kv_slow, kd_aux=kd_salt)
1464 if (cs%boundary_forcing_tendency_diag)
then
1465 do k=1,nz ;
do j=js,je ;
do i=is,ie
1466 h_diag(i,j,k) = h(i,j,k)
1467 temp_diag(i,j,k) = tv%T(i,j,k)
1468 saln_diag(i,j,k) = tv%S(i,j,k)
1469 enddo ;
enddo ;
enddo
1473 call cpu_clock_begin(id_clock_remap)
1476 do k=1,nz ;
do j=js,je ;
do i=is,ie
1477 h_prebound(i,j,k) = h(i,j,k)
1478 enddo ;
enddo ;
enddo
1479 if (cs%use_energetic_PBL)
then
1481 skinbuoyflux(:,:) = 0.0
1482 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
1483 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, cs%evap_CFL_limit, &
1484 cs%minimum_forcing_depth, ctke, dsv_dt, dsv_ds, skinbuoyflux=skinbuoyflux)
1487 call hchksum(ea_t,
"after applyBoundaryFluxes ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
1488 call hchksum(eb_t,
"after applyBoundaryFluxes eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
1489 call hchksum(ea_s,
"after applyBoundaryFluxes ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
1490 call hchksum(eb_s,
"after applyBoundaryFluxes eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
1491 call hchksum(ctke,
"after applyBoundaryFluxes cTKE", g%HI, haloshift=0, &
1492 scale=us%RZ3_T3_to_W_m2*us%T_to_s)
1493 call hchksum(dsv_dt,
"after applyBoundaryFluxes dSV_dT", g%HI, haloshift=0, scale=us%kg_m3_to_R)
1494 call hchksum(dsv_ds,
"after applyBoundaryFluxes dSV_dS", g%HI, haloshift=0, scale=us%kg_m3_to_R)
1497 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
1498 call energetic_pbl(h, u_h, v_h, tv, fluxes, dt, kd_epbl, g, gv, us, &
1499 cs%energetic_PBL_CSp, dsv_dt, dsv_ds, ctke, skinbuoyflux, waves=waves)
1501 if (
associated(hml))
then
1502 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, hml(:,:), g, us)
1503 call pass_var(hml, g%domain, halo=1)
1505 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
1506 elseif (
associated(visc%MLD))
then
1507 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, visc%MLD, g, us)
1508 call pass_var(visc%MLD, g%domain, halo=1)
1512 do k=2,nz ;
do j=js,je ;
do i=is,ie
1513 if (cs%ePBL_is_additive)
then
1514 kd_add_here = kd_epbl(i,j,k)
1515 visc%Kv_shear(i,j,k) = visc%Kv_shear(i,j,k) + cs%ePBL_Prandtl*kd_epbl(i,j,k)
1517 kd_add_here = max(kd_epbl(i,j,k) - visc%Kd_shear(i,j,k), 0.0)
1518 visc%Kv_shear(i,j,k) = max(visc%Kv_shear(i,j,k), cs%ePBL_Prandtl*kd_epbl(i,j,k))
1521 kd_heat(i,j,k) = kd_heat(i,j,k) + kd_add_here
1522 kd_salt(i,j,k) = kd_salt(i,j,k) + kd_add_here
1524 enddo ;
enddo ;
enddo
1527 call hchksum(ea_t,
"after ePBL ea_t", g%HI, haloshift=0, scale=gv%H_to_m)
1528 call hchksum(eb_t,
"after ePBL eb_t", g%HI, haloshift=0, scale=gv%H_to_m)
1529 call hchksum(ea_s,
"after ePBL ea_s", g%HI, haloshift=0, scale=gv%H_to_m)
1530 call hchksum(eb_s,
"after ePBL eb_s", g%HI, haloshift=0, scale=gv%H_to_m)
1531 call hchksum(kd_epbl,
"after ePBL Kd_ePBL", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
1535 call applyboundaryfluxesinout(cs%diabatic_aux_CSp, g, gv, us, dt, fluxes, cs%optics, &
1536 optics_nbands(cs%optics), h, tv, cs%aggregate_FW_forcing, &
1537 cs%evap_CFL_limit, cs%minimum_forcing_depth)
1544 if (cs%boundary_forcing_tendency_diag)
then
1545 call diagnose_boundary_forcing_tendency(tv, h, temp_diag, saln_diag, h_diag, dt, g, gv, us, cs)
1546 if (cs%id_boundary_forcing_h > 0)
call post_data(cs%id_boundary_forcing_h, h, cs%diag, alt_h = h_diag)
1549 call diag_update_remap_grids(cs%diag)
1550 call cpu_clock_end(id_clock_remap)
1552 call mom_forcing_chksum(
"after applyBoundaryFluxes ", fluxes, g, us, haloshift=0)
1553 call mom_thermovar_chksum(
"after applyBoundaryFluxes ", tv, g)
1554 call mom_state_chksum(
"after applyBoundaryFluxes ", u, v, h, g, gv, us, haloshift=0)
1556 if (showcalltree)
call calltree_waypoint(
"done with applyBoundaryFluxes (diabatic)")
1557 if (cs%debugConservation)
call mom_state_stats(
'applyBoundaryFluxes', u, v, h, tv%T, tv%S, g, gv, us)
1559 if (showcalltree)
call calltree_waypoint(
"done with h=ea-eb (diabatic)")
1560 if (cs%debugConservation)
call mom_state_stats(
'h=ea-eb', u, v, h, tv%T, tv%S, g, gv, us)
1563 if (
associated(tv%T))
then
1566 call hchksum(ea_t,
"before triDiagTS ea_t ", g%HI, haloshift=0, scale=gv%H_to_m)
1567 call hchksum(eb_t,
"before triDiagTS eb_t ", g%HI, haloshift=0, scale=gv%H_to_m)
1568 call hchksum(ea_s,
"before triDiagTS ea_s ", g%HI, haloshift=0, scale=gv%H_to_m)
1569 call hchksum(eb_s,
"before triDiagTS eb_s ", g%HI, haloshift=0, scale=gv%H_to_m)
1572 call cpu_clock_begin(id_clock_tridiag)
1579 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
1580 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
1582 if (cs%diabatic_diff_tendency_diag)
then
1583 do k=1,nz ;
do j=js,je ;
do i=is,ie
1584 temp_diag(i,j,k) = tv%T(i,j,k)
1585 saln_diag(i,j,k) = tv%S(i,j,k)
1586 enddo ;
enddo ;
enddo
1590 do j=js,je ;
do i=is,ie
1591 ea_t(i,j,1) = 0.; ea_s(i,j,1) = 0.
1595 do k=2,nz ;
do j=js,je ;
do i=is,ie
1596 hval = 1.0 / (h_neglect + 0.5*(h(i,j,k-1) + h(i,j,k)))
1597 ea_t(i,j,k) = (gv%Z_to_H**2) * dt * hval * kd_heat(i,j,k)
1598 eb_t(i,j,k-1) = ea_t(i,j,k)
1599 ea_s(i,j,k) = (gv%Z_to_H**2) * dt * hval * kd_salt(i,j,k)
1600 eb_s(i,j,k-1) = ea_s(i,j,k)
1601 enddo ;
enddo ;
enddo
1602 do j=js,je ;
do i=is,ie
1603 eb_t(i,j,nz) = 0. ; eb_s(i,j,nz) = 0.
1605 if (showcalltree)
call calltree_waypoint(
"done setting ea_t,ea_s,eb_t,eb_s from Kd_heat" //&
1606 "and Kd_salt (diabatic)")
1612 ea_t(i,js-1,k) = 0.0 ; eb_t(i,js-1,k) = 0.0
1613 ea_s(i,js-1,k) = 0.0 ; eb_s(i,js-1,k) = 0.0
1614 ea_t(i,je+1,k) = 0.0 ; eb_t(i,je+1,k) = 0.0
1615 ea_s(i,je+1,k) = 0.0 ; eb_s(i,je+1,k) = 0.0
1618 ea_t(is-1,j,k) = 0.0 ; eb_t(is-1,j,k) = 0.0
1619 ea_s(is-1,j,k) = 0.0 ; eb_s(is-1,j,k) = 0.0
1620 ea_t(ie+1,j,k) = 0.0 ; eb_t(ie+1,j,k) = 0.0
1621 ea_s(ie+1,j,k) = 0.0 ; eb_s(ie+1,j,k) = 0.0
1626 call tracer_vertdiff(h, ea_t, eb_t, dt, tv%T, g, gv)
1627 call tracer_vertdiff(h, ea_s, eb_s, dt, tv%S, g, gv)
1631 if (cs%diabatic_diff_tendency_diag)
then
1632 call diagnose_diabatic_diff_tendency(tv, h, temp_diag, saln_diag, dt, g, gv, us, cs)
1634 call cpu_clock_end(id_clock_tridiag)
1636 if (showcalltree)
call calltree_waypoint(
"done with triDiagTS (diabatic)")
1640 if (cs%debugConservation)
call mom_state_stats(
'triDiagTS', u, v, h, tv%T, tv%S, g, gv, us)
1643 call mom_state_chksum(
"after mixed layer ", u, v, h, g, gv, us, haloshift=0)
1644 call mom_thermovar_chksum(
"after mixed layer ", tv, g)
1649 call diag_update_remap_grids(cs%diag)
1653 if ((cs%id_Tdif > 0) .or. (cs%id_Tadv > 0))
then
1654 do j=js,je ;
do i=is,ie
1655 tdif_flx(i,j,1) = 0.0 ; tdif_flx(i,j,nz+1) = 0.0
1656 tadv_flx(i,j,1) = 0.0 ; tadv_flx(i,j,nz+1) = 0.0
1659 do k=2,nz ;
do j=js,je ;
do i=is,ie
1660 tdif_flx(i,j,k) = (idt * 0.5*(ea_t(i,j,k) + eb_t(i,j,k-1))) * &
1661 (tv%T(i,j,k-1) - tv%T(i,j,k))
1662 tadv_flx(i,j,k) = (idt * (ea_t(i,j,k) - eb_t(i,j,k-1))) * &
1663 0.5*(tv%T(i,j,k-1) + tv%T(i,j,k))
1664 enddo ;
enddo ;
enddo
1666 if ((cs%id_Sdif > 0) .or. (cs%id_Sadv > 0))
then
1667 do j=js,je ;
do i=is,ie
1668 sdif_flx(i,j,1) = 0.0 ; sdif_flx(i,j,nz+1) = 0.0
1669 sadv_flx(i,j,1) = 0.0 ; sadv_flx(i,j,nz+1) = 0.0
1672 do k=2,nz ;
do j=js,je ;
do i=is,ie
1673 sdif_flx(i,j,k) = (idt * 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))) * &
1674 (tv%S(i,j,k-1) - tv%S(i,j,k))
1675 sadv_flx(i,j,k) = (idt * (ea_s(i,j,k) - eb_s(i,j,k-1))) * &
1676 0.5*(tv%S(i,j,k-1) + tv%S(i,j,k))
1677 enddo ;
enddo ;
enddo
1681 call cpu_clock_begin(id_clock_tracers)
1683 if (cs%mix_boundary_tracers)
then
1684 tr_ea_bbl = gv%Z_to_H * sqrt(dt*cs%Kd_BBL_tr)
1688 ebtr(i,j,nz) = eb_s(i,j,nz)
1690 in_boundary(i) = (g%mask2dT(i,j) > 0.0)
1692 do k=nz,2,-1 ;
do i=is,ie
1693 if (in_boundary(i))
then
1694 htot(i) = htot(i) + h(i,j,k)
1703 add_ent = ((dt * cs%Kd_min_tr) * gv%Z_to_H**2) * &
1704 ((h(i,j,k-1)+h(i,j,k)+h_neglect) / &
1705 (h(i,j,k-1)*h(i,j,k)+h_neglect2)) - &
1706 0.5*(ea_s(i,j,k) + eb_s(i,j,k-1))
1707 if (htot(i) < tr_ea_bbl)
then
1708 add_ent = max(0.0, add_ent, &
1709 (tr_ea_bbl - htot(i)) - min(ea_s(i,j,k),eb_s(i,j,k-1)))
1710 elseif (add_ent < 0.0)
then
1711 add_ent = 0.0 ; in_boundary(i) = .false.
1714 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1715 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1717 ebtr(i,j,k-1) = eb_s(i,j,k-1) ; eatr(i,j,k) = ea_s(i,j,k)
1720 if (cs%double_diffuse)
then ;
if (kd_extra_s(i,j,k) > 0.0)
then
1721 add_ent = ((dt * kd_extra_s(i,j,k)) * gv%Z_to_H**2) / &
1722 (0.5 * (h(i,j,k-1) + h(i,j,k)) + &
1724 ebtr(i,j,k-1) = ebtr(i,j,k-1) + add_ent
1725 eatr(i,j,k) = eatr(i,j,k) + add_ent
1728 do i=is,ie ; eatr(i,j,1) = ea_s(i,j,1) ;
enddo
1735 call call_tracer_column_fns(h_prebound, h, ea_s, eb_s, fluxes, hml, dt, g, gv, us, tv, &
1736 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1737 evap_cfl_limit = cs%evap_CFL_limit, &
1738 minimum_forcing_depth=cs%minimum_forcing_depth)
1740 elseif (cs%double_diffuse)
then
1742 do j=js,je ;
do i=is,ie
1743 ebtr(i,j,nz) = eb_s(i,j,nz) ; eatr(i,j,1) = ea_s(i,j,1)
1746 do k=nz,2,-1 ;
do j=js,je ;
do i=is,ie
1747 if (kd_extra_s(i,j,k) > 0.0)
then
1748 add_ent = ((dt * kd_extra_s(i,j,k)) * gv%Z_to_H**2) / &
1749 (0.5 * (h(i,j,k-1) + h(i,j,k)) + &
1754 ebtr(i,j,k-1) = eb_s(i,j,k-1) + add_ent
1755 eatr(i,j,k) = ea_s(i,j,k) + add_ent
1756 enddo ;
enddo ;
enddo
1759 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1760 cs%optics, cs%tracer_flow_CSp, cs%debug,&
1761 evap_cfl_limit = cs%evap_CFL_limit, &
1762 minimum_forcing_depth=cs%minimum_forcing_depth)
1766 call call_tracer_column_fns(h_prebound, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
1767 cs%optics, cs%tracer_flow_CSp, cs%debug, &
1768 evap_cfl_limit = cs%evap_CFL_limit, &
1769 minimum_forcing_depth=cs%minimum_forcing_depth)
1772 call cpu_clock_end(id_clock_tracers)
1775 if (cs%use_sponge)
then
1776 call cpu_clock_begin(id_clock_sponge)
1777 if (
associated(cs%ALE_sponge_CSp))
then
1778 call apply_ale_sponge(h, dt, g, gv, us, cs%ALE_sponge_CSp, cs%Time)
1781 call cpu_clock_end(id_clock_sponge)
1784 call mom_thermovar_chksum(
"apply_sponge ", tv, g)
1788 call cpu_clock_begin(id_clock_pass)
1789 if (g%symmetric)
then ; dir_flag = to_all+omit_corners
1790 else ; dir_flag = to_west+to_south+omit_corners ;
endif
1795 call do_group_pass(cs%pass_hold_eb_ea, g%Domain)
1797 if (
associated(visc%Kv_slow)) &
1798 call pass_var(visc%Kv_slow, g%Domain, to_all+omit_corners, halo=1)
1799 call cpu_clock_end(id_clock_pass)
1801 call disable_averaging(cs%diag)
1804 call enable_averages(dt, time_end, cs%diag)
1806 if (cs%id_Kd_interface > 0)
call post_data(cs%id_Kd_interface, kd_int, cs%diag)
1807 if (cs%id_Kd_heat > 0)
call post_data(cs%id_Kd_heat, kd_heat, cs%diag)
1808 if (cs%id_Kd_salt > 0)
call post_data(cs%id_Kd_salt, kd_salt, cs%diag)
1809 if (cs%id_Kd_ePBL > 0)
call post_data(cs%id_Kd_ePBL, kd_epbl, cs%diag)
1811 if (cs%id_ea_t > 0)
call post_data(cs%id_ea_t, ea_t, cs%diag)
1812 if (cs%id_eb_t > 0)
call post_data(cs%id_eb_t, eb_t, cs%diag)
1813 if (cs%id_ea_s > 0)
call post_data(cs%id_ea_s, ea_s, cs%diag)
1814 if (cs%id_eb_s > 0)
call post_data(cs%id_eb_s, eb_s, cs%diag)
1816 if (cs%id_dudt_dia > 0)
call post_data(cs%id_dudt_dia, adp%du_dt_dia, cs%diag)
1817 if (cs%id_dvdt_dia > 0)
call post_data(cs%id_dvdt_dia, adp%dv_dt_dia, cs%diag)
1819 if (cs%id_Tdif > 0)
call post_data(cs%id_Tdif, tdif_flx, cs%diag)
1820 if (cs%id_Tadv > 0)
call post_data(cs%id_Tadv, tadv_flx, cs%diag)
1821 if (cs%id_Sdif > 0)
call post_data(cs%id_Sdif, sdif_flx, cs%diag)
1822 if (cs%id_Sadv > 0)
call post_data(cs%id_Sadv, sadv_flx, cs%diag)
1825 if (cs%id_hf_dudt_dia_2d > 0)
then
1826 allocate(hf_dudt_dia_2d(g%IsdB:g%IedB,g%jsd:g%jed))
1827 hf_dudt_dia_2d(:,:) = 0.0
1828 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1829 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)
1830 enddo ;
enddo ;
enddo
1831 call post_data(cs%id_hf_dudt_dia_2d, hf_dudt_dia_2d, cs%diag)
1832 deallocate(hf_dudt_dia_2d)
1835 if (cs%id_hf_dvdt_dia_2d > 0)
then
1836 allocate(hf_dvdt_dia_2d(g%isd:g%ied,g%JsdB:g%JedB))
1837 hf_dvdt_dia_2d(:,:) = 0.0
1838 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1839 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)
1840 enddo ;
enddo ;
enddo
1841 call post_data(cs%id_hf_dvdt_dia_2d, hf_dvdt_dia_2d, cs%diag)
1842 deallocate(hf_dvdt_dia_2d)
1845 call disable_averaging(cs%diag)
1847 if (showcalltree)
call calltree_leave(
"diabatic_ALE()")
1849 end subroutine diabatic_ale
1853 subroutine layered_diabatic(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end, &
1854 G, GV, US, CS, WAVES)
1858 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: u
1859 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: v
1860 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
1863 real,
dimension(:,:),
pointer :: Hml
1864 type(
forcing),
intent(inout) :: fluxes
1871 real,
intent(in) :: dt
1872 type(time_type),
intent(in) :: Time_end
1876 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
1891 real,
dimension(SZI_(G),SZJ_(G)) :: &
1894 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: h_diag
1895 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: temp_diag
1896 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: saln_diag
1897 real,
dimension(SZI_(G),SZJ_(G)) :: tendency_2d
1901 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
target :: &
1907 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
1921 real,
allocatable,
dimension(:,:) :: &
1922 hf_dudt_dia_2d, hf_dvdt_dia_2d
1925 real,
pointer,
dimension(:,:,:) :: &
1931 integer :: kb(SZI_(G),SZJ_(G))
1934 real :: p_ref_cv(SZI_(G))
1937 logical :: in_boundary(SZI_(G))
1957 real :: htot(SZIB_(G))
1958 real :: b1(SZIB_(G)), d1(SZIB_(G))
1959 real :: c1(SZIB_(G),SZK_(G))
1966 logical :: showCallTree
1967 integer,
dimension(2) :: EOSdom
1968 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, nkmb, m, halo
1970 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1971 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1972 nkmb = gv%nk_rho_varies
1973 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect*h_neglect
1974 kd_heat(:,:,:) = 0.0 ; kd_salt(:,:,:) = 0.0
1977 showcalltree = calltree_showquery()
1978 if (showcalltree)
call calltree_enter(
"layered_diabatic(), MOM_diabatic_driver.F90")
1981 eaml => eatr ; ebml => ebtr
1984 call enable_averages(dt, time_end, cs%diag)
1986 if ((cs%ML_mix_first > 0.0) .or. cs%use_geothermal)
then
1987 halo = cs%halo_TS_diff
1989 do k=1,nz ;
do j=js-halo,je+halo ;
do i=is-halo,ie+halo
1990 h_orig(i,j,k) = h(i,j,k) ; eaml(i,j,k) = 0.0 ; ebml(i,j,k) = 0.0
1991 enddo ;
enddo ;
enddo
1994 if (cs%use_geothermal)
then
1995 call cpu_clock_begin(id_clock_geothermal)
1996 call geothermal_entraining(h, tv, dt, eaml, ebml, g, gv, us, cs%geothermal_CSp, halo=cs%halo_TS_diff)
1997 call cpu_clock_end(id_clock_geothermal)
1998 if (showcalltree)
call calltree_waypoint(
"geothermal (diabatic)")
1999 if (cs%debugConservation)
call mom_state_stats(
'geothermal', u, v, h, tv%T, tv%S, g, gv, us)
2004 call diag_update_remap_grids(cs%diag)
2009 if (
associated(cs%optics)) &
2010 call set_pen_shortwave(cs%optics, fluxes, g, gv, us, cs%diabatic_aux_CSp, cs%opacity_CSp, cs%tracer_flow_CSp)
2012 if (cs%bulkmixedlayer)
then
2013 if (cs%debug)
call mom_forcing_chksum(
"Before mixedlayer", fluxes, g, us, haloshift=0)
2015 if (cs%ML_mix_first > 0.0)
then
2023 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
2025 call cpu_clock_begin(id_clock_mixedlayer)
2026 if (cs%ML_mix_first < 1.0)
then
2028 call bulkmixedlayer(h, u_h, v_h, tv, fluxes, dt*cs%ML_mix_first, &
2029 eaml,ebml, g, gv, us, cs%bulkmixedlayer_CSp, cs%optics, &
2030 hml, cs%aggregate_FW_forcing, dt, last_call=.false.)
2032 call bulkmixedlayer(h, u_h, v_h, tv, fluxes, dt, eaml, ebml, &
2033 g, gv, us, cs%bulkmixedlayer_CSp, cs%optics, &
2034 hml, cs%aggregate_FW_forcing, dt, last_call=.true.)
2043 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
2044 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
2045 call cpu_clock_end(id_clock_mixedlayer)
2047 call mom_state_chksum(
"After mixedlayer ", u, v, h, g, gv, us, haloshift=0)
2048 call mom_forcing_chksum(
"After mixedlayer", fluxes, g, us, haloshift=0)
2050 if (showcalltree)
call calltree_waypoint(
"done with 1st bulkmixedlayer (diabatic)")
2051 if (cs%debugConservation)
call mom_state_stats(
'1st bulkmixedlayer', u, v, h, tv%T, tv%S, g, gv, us)
2056 call mom_state_chksum(
"before find_uv_at_h", u, v, h, g, gv, us, haloshift=0)
2057 if (cs%use_kappa_shear .or. cs%use_CVMix_shear)
then
2058 if ((cs%ML_mix_first > 0.0) .or. cs%use_geothermal)
then
2059 call find_uv_at_h(u, v, h_orig, u_h, v_h, g, gv, us, eaml, ebml)
2061 call hchksum(eaml,
"after find_uv_at_h eaml", g%HI, scale=gv%H_to_m)
2062 call hchksum(ebml,
"after find_uv_at_h ebml", g%HI, scale=gv%H_to_m)
2065 call find_uv_at_h(u, v, h, u_h, v_h, g, gv, us)
2067 if (showcalltree)
call calltree_waypoint(
"done with find_uv_at_h (diabatic)")
2070 call cpu_clock_begin(id_clock_set_diffusivity)
2073 if ((cs%halo_TS_diff > 0) .and. (cs%ML_mix_first > 0.0))
then
2074 if (
associated(tv%T))
call pass_var(tv%T, g%Domain, halo=cs%halo_TS_diff, complete=.false.)
2075 if (
associated(tv%T))
call pass_var(tv%S, g%Domain, halo=cs%halo_TS_diff, complete=.false.)
2076 call pass_var(h, g%domain, halo=cs%halo_TS_diff, complete=.true.)
2079 call mom_state_chksum(
"before set_diffusivity", u, v, h, g, gv, us, haloshift=cs%halo_TS_diff)
2080 if (cs%double_diffuse)
then
2081 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, visc, dt, g, gv, us, cs%set_diff_CSp, &
2082 kd_lay=kd_lay, kd_int=kd_int, kd_extra_t=kd_extra_t, kd_extra_s=kd_extra_s)
2084 call set_diffusivity(u, v, h, u_h, v_h, tv, fluxes, cs%optics, visc, dt, g, gv, us, &
2085 cs%set_diff_CSp, kd_lay=kd_lay, kd_int=kd_int)
2087 call cpu_clock_end(id_clock_set_diffusivity)
2088 if (showcalltree)
call calltree_waypoint(
"done with set_diffusivity (diabatic)")
2091 call mom_state_chksum(
"after set_diffusivity ", u, v, h, g, gv, us, haloshift=0)
2092 call mom_forcing_chksum(
"after set_diffusivity ", fluxes, g, us, haloshift=0)
2093 call mom_thermovar_chksum(
"after set_diffusivity ", tv, g)
2094 call hchksum(kd_lay,
"after set_diffusivity Kd_lay", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2095 call hchksum(kd_int,
"after set_diffusivity Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2100 call cpu_clock_begin(id_clock_kpp)
2106 call calculatebuoyancyflux2d(g, gv, us, fluxes, cs%optics, h, tv%T, tv%S, tv, &
2107 cs%KPP_buoy_flux, cs%KPP_temp_flux, cs%KPP_salt_flux)
2112 if (cs%double_diffuse)
then
2115 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2116 kd_salt(i,j,k) = kd_int(i,j,k) + kd_extra_s(i,j,k)
2117 kd_heat(i,j,k) = kd_int(i,j,k) + kd_extra_t(i,j,k)
2118 enddo ;
enddo ;
enddo
2121 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2122 kd_salt(i,j,k) = kd_int(i,j,k)
2123 kd_heat(i,j,k) = kd_int(i,j,k)
2124 enddo ;
enddo ;
enddo
2127 call kpp_compute_bld(cs%KPP_CSp, g, gv, us, h, tv%T, tv%S, u, v, tv, &
2128 fluxes%ustar, cs%KPP_buoy_flux, waves=waves)
2130 call kpp_calculate(cs%KPP_CSp, g, gv, us, h, fluxes%ustar, cs%KPP_buoy_flux, kd_heat, &
2131 kd_salt, visc%Kv_shear, cs%KPP_NLTheat, cs%KPP_NLTscalar, waves=waves)
2133 if (
associated(hml))
then
2134 call kpp_get_bld(cs%KPP_CSp, hml(:,:), g, us)
2135 call pass_var(hml, g%domain, halo=1)
2137 if (
associated(visc%MLD)) visc%MLD(:,:) = hml(:,:)
2140 if (.not. cs%KPPisPassive)
then
2142 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2143 kd_int(i,j,k) = min( kd_salt(i,j,k), kd_heat(i,j,k) )
2144 enddo ;
enddo ;
enddo
2145 if (cs%double_diffuse)
then
2147 do k=1,nz+1 ;
do j=js,je ;
do i=is,ie
2148 kd_extra_s(i,j,k) = (kd_salt(i,j,k) - kd_int(i,j,k))
2149 kd_extra_t(i,j,k) = (kd_heat(i,j,k) - kd_int(i,j,k))
2150 enddo ;
enddo ;
enddo
2154 call cpu_clock_end(id_clock_kpp)
2155 if (showcalltree)
call calltree_waypoint(
"done with KPP_calculate (diabatic)")
2158 call mom_forcing_chksum(
"after KPP", fluxes, g, us, haloshift=0)
2159 call mom_thermovar_chksum(
"after KPP", tv, g)
2160 call hchksum(kd_lay,
"after KPP Kd_lay", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2161 call hchksum(kd_int,
"after KPP Kd_Int", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
2167 if (cs%use_CVMix_conv)
then
2168 call calculate_cvmix_conv(h, tv, g, gv, us, cs%CVMix_conv_csp, hml, kd=kd_int, kv=visc%Kv_slow)
2172 call cpu_clock_begin(id_clock_kpp)
2174 call hchksum(cs%KPP_temp_flux,
"before KPP_applyNLT netHeat", g%HI, haloshift=0, scale=gv%H_to_m)
2175 call hchksum(cs%KPP_salt_flux,
"before KPP_applyNLT netSalt", g%HI, haloshift=0, scale=gv%H_to_m)
2176 call hchksum(cs%KPP_NLTheat,
"before KPP_applyNLT NLTheat", g%HI, haloshift=0)
2177 call hchksum(cs%KPP_NLTscalar,
"before KPP_applyNLT NLTscalar", g%HI, haloshift=0)
2181 call kpp_nonlocaltransport_temp(cs%KPP_CSp, g, gv, h, cs%KPP_NLTheat, cs%KPP_temp_flux, &
2182 us%T_to_s*dt, tv%T, us%Q_to_J_kg*tv%C_p)
2183 call kpp_nonlocaltransport_saln(cs%KPP_CSp, g, gv, h, cs%KPP_NLTscalar, cs%KPP_salt_flux, &
2185 call cpu_clock_end(id_clock_kpp)
2186 if (showcalltree)
call calltree_waypoint(
"done with KPP_applyNonLocalTransport (diabatic)")
2187 if (cs%debugConservation)
call mom_state_stats(
'KPP_applyNonLocalTransport', u, v, h, tv%T, tv%S, g, gv, us)
2190 call mom_state_chksum(
"after KPP_applyNLT ", u, v, h, g, gv, us, haloshift=0)
2191 call mom_forcing_chksum(
"after KPP_applyNLT ", fluxes, g, us, haloshift=0)
2192 call mom_thermovar_chksum(
"after KPP_applyNLT ", tv, g)
2198 if (cs%double_diffuse .and.
associated(tv%T))
then
2200 call cpu_clock_begin(id_clock_differential_diff)
2201 call differential_diffuse_t_s(h, tv%T, tv%S, kd_extra_t, kd_extra_s, dt, g, gv)
2202 call cpu_clock_end(id_clock_differential_diff)
2203 if (showcalltree)
call calltree_waypoint(
"done with differential_diffuse_T_S (diabatic)")
2204 if (cs%debugConservation)
call mom_state_stats(
'differential_diffuse_T_S', u, v, h, tv%T, tv%S, g, gv, us)
2208 if (.not. cs%useKPP)
then
2210 do k=2,nz ;
do j=js,je ;
do i=is,ie
2211 kd_heat(i,j,k) = kd_heat(i,j,k) + kd_extra_t(i,j,k)
2212 kd_salt(i,j,k) = kd_salt(i,j,k) + kd_extra_s(i,j,k)
2213 enddo ;
enddo ;
enddo
2220 call cpu_clock_begin(id_clock_entrain)
2223 call entrainment_diffusive(h, tv, fluxes, dt, g, gv, us, cs%entrain_diffusive_CSp, &
2224 ea, eb, kb, kd_lay=kd_lay, kd_int=kd_int)
2225 call cpu_clock_end(id_clock_entrain)
2226 if (showcalltree)
call calltree_waypoint(
"done with Entrainment_diffusive (diabatic)")
2229 call mom_forcing_chksum(
"after calc_entrain ", fluxes, g, us, haloshift=0)
2230 call mom_thermovar_chksum(
"after calc_entrain ", tv, g)
2231 call mom_state_chksum(
"after calc_entrain ", u, v, h, g, gv, us, haloshift=0)
2232 call hchksum(ea,
"after calc_entrain ea", g%HI, haloshift=0, scale=gv%H_to_m)
2233 call hchksum(eb,
"after calc_entrain eb", g%HI, haloshift=0, scale=gv%H_to_m)
2237 if (cs%boundary_forcing_tendency_diag)
then
2238 do k=1,nz ;
do j=js,je ;
do i=is,ie
2239 h_diag(i,j,k) = h(i,j,k)
2240 temp_diag(i,j,k) = tv%T(i,j,k)
2241 saln_diag(i,j,k) = tv%S(i,j,k)
2242 enddo ;
enddo ;
enddo
2256 hold(i,j,1) = h(i,j,1)
2257 h(i,j,1) = h(i,j,1) + (eb(i,j,1) - ea(i,j,2))
2258 hold(i,j,nz) = h(i,j,nz)
2259 h(i,j,nz) = h(i,j,nz) + (ea(i,j,nz) - eb(i,j,nz-1))
2260 if (h(i,j,1) <= 0.0)
then
2261 h(i,j,1) = gv%Angstrom_H
2263 if (h(i,j,nz) <= 0.0)
then
2264 h(i,j,nz) = gv%Angstrom_H
2267 do k=2,nz-1 ;
do i=is,ie
2268 hold(i,j,k) = h(i,j,k)
2269 h(i,j,k) = h(i,j,k) + ((ea(i,j,k) - eb(i,j,k-1)) + &
2270 (eb(i,j,k) - ea(i,j,k+1)))
2271 if (h(i,j,k) <= 0.0)
then
2272 h(i,j,k) = gv%Angstrom_H
2277 call diag_update_remap_grids(cs%diag)
2280 call mom_state_chksum(
"after negative check ", u, v, h, g, gv, us, haloshift=0)
2281 call mom_forcing_chksum(
"after negative check ", fluxes, g, us, haloshift=0)
2282 call mom_thermovar_chksum(
"after negative check ", tv, g)
2284 if (showcalltree)
call calltree_waypoint(
"done with h=ea-eb (diabatic)")
2285 if (cs%debugConservation)
call mom_state_stats(
'h=ea-eb', u, v, h, tv%T, tv%S, g, gv, us)
2291 if (cs%bulkmixedlayer)
then
2293 if (
associated(tv%T))
then
2294 call cpu_clock_begin(id_clock_tridiag)
2299 if (cs%massless_match_targets)
then
2303 h_tr = hold(i,j,1) + h_neglect
2304 b1(i) = 1.0 / (h_tr + eb(i,j,1))
2305 d1(i) = h_tr * b1(i)
2306 tv%T(i,j,1) = b1(i) * (h_tr*tv%T(i,j,1))
2307 tv%S(i,j,1) = b1(i) * (h_tr*tv%S(i,j,1))
2309 do k=2,nkmb ;
do i=is,ie
2310 c1(i,k) = eb(i,j,k-1) * b1(i)
2311 h_tr = hold(i,j,k) + h_neglect
2312 b_denom_1 = h_tr + d1(i)*ea(i,j,k)
2313 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
2314 if (k<nkmb) d1(i) = b_denom_1 * b1(i)
2315 tv%T(i,j,k) = b1(i) * (h_tr*tv%T(i,j,k) + ea(i,j,k)*tv%T(i,j,k-1))
2316 tv%S(i,j,k) = b1(i) * (h_tr*tv%S(i,j,k) + ea(i,j,k)*tv%S(i,j,k-1))
2319 do k=nkmb+1,nz ;
do i=is,ie
2320 if (k == kb(i,j))
then
2321 c1(i,k) = eb(i,j,k-1) * b1(i)
2322 d1(i) = (((eb(i,j,nkmb)-eb(i,j,k-1)) + hold(i,j,nkmb) + h_neglect) + &
2323 d1(i)*ea(i,j,nkmb)) * b1(i)
2324 h_tr = hold(i,j,k) + h_neglect
2325 b_denom_1 = h_tr + d1(i)*ea(i,j,k)
2326 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
2327 d1(i) = b_denom_1 * b1(i)
2328 tv%T(i,j,k) = b1(i) * (h_tr*tv%T(i,j,k) + ea(i,j,k)*tv%T(i,j,nkmb))
2329 tv%S(i,j,k) = b1(i) * (h_tr*tv%S(i,j,k) + ea(i,j,k)*tv%S(i,j,nkmb))
2330 elseif (k > kb(i,j))
then
2331 c1(i,k) = eb(i,j,k-1) * b1(i)
2332 h_tr = hold(i,j,k) + h_neglect
2333 b_denom_1 = h_tr + d1(i)*ea(i,j,k)
2334 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
2335 d1(i) = b_denom_1 * b1(i)
2336 tv%T(i,j,k) = b1(i) * (h_tr*tv%T(i,j,k) + ea(i,j,k)*tv%T(i,j,k-1))
2337 tv%S(i,j,k) = b1(i) * (h_tr*tv%S(i,j,k) + ea(i,j,k)*tv%S(i,j,k-1))
2338 elseif (eb(i,j,k) < eb(i,j,k-1))
then
2343 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)
2344 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)
2348 do k=nz-1,nkmb,-1 ;
do i=is,ie
2349 if (k >= kb(i,j))
then
2350 tv%T(i,j,k) = tv%T(i,j,k) + c1(i,k+1)*tv%T(i,j,k+1)
2351 tv%S(i,j,k) = tv%S(i,j,k) + c1(i,k+1)*tv%S(i,j,k+1)
2354 do i=is,ie ;
if (kb(i,j) <= nz)
then
2355 tv%T(i,j,nkmb) = tv%T(i,j,nkmb) + c1(i,kb(i,j))*tv%T(i,j,kb(i,j))
2356 tv%S(i,j,nkmb) = tv%S(i,j,nkmb) + c1(i,kb(i,j))*tv%S(i,j,kb(i,j))
2358 do k=nkmb-1,1,-1 ;
do i=is,ie
2359 tv%T(i,j,k) = tv%T(i,j,k) + c1(i,k+1)*tv%T(i,j,k+1)
2360 tv%S(i,j,k) = tv%S(i,j,k) + c1(i,k+1)*tv%S(i,j,k+1)
2367 if (cs%tracer_tridiag)
then
2368 call tracer_vertdiff(hold, ea, eb, dt, tv%T, g, gv)
2369 call tracer_vertdiff(hold, ea, eb, dt, tv%S, g, gv)
2371 call tridiagts(g, gv, is, ie, js, je, hold, ea, eb, tv%T, tv%S)
2374 call cpu_clock_end(id_clock_tridiag)
2377 if (cs%debugConservation)
call mom_state_stats(
'BML tridiag', u, v, h, tv%T, tv%S, g, gv, us)
2379 if ((cs%ML_mix_first > 0.0) .or. cs%use_geothermal)
then
2383 do k=1,nz ;
do j=js,je ;
do i=is,ie
2384 hold(i,j,k) = h_orig(i,j,k)
2385 ea(i,j,k) = ea(i,j,k) + eaml(i,j,k)
2386 eb(i,j,k) = eb(i,j,k) + ebml(i,j,k)
2387 enddo ;
enddo ;
enddo
2389 call hchksum(ea,
"after ea = ea + eaml", g%HI, haloshift=0, scale=gv%H_to_m)
2390 call hchksum(eb,
"after eb = eb + ebml", g%HI, haloshift=0, scale=gv%H_to_m)
2394 if (cs%ML_mix_first < 1.0)
then
2403 call find_uv_at_h(u, v, hold, u_h, v_h, g, gv, us, ea, eb)
2404 if (cs%debug)
call mom_state_chksum(
"find_uv_at_h1 ", u, v, h, g, gv, us, haloshift=0)
2406 dt_mix = min(dt, dt*(1.0 - cs%ML_mix_first))
2407 call cpu_clock_begin(id_clock_mixedlayer)
2409 call bulkmixedlayer(h, u_h, v_h, tv, fluxes, dt_mix, ea, eb, &
2410 g, gv, us, cs%bulkmixedlayer_CSp, cs%optics, &
2411 hml, cs%aggregate_FW_forcing, dt, last_call=.true.)
2419 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
2420 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
2422 call cpu_clock_end(id_clock_mixedlayer)
2423 if (showcalltree)
call calltree_waypoint(
"done with 2nd bulkmixedlayer (diabatic)")
2424 if (cs%debugConservation)
call mom_state_stats(
'2nd bulkmixedlayer', u, v, h, tv%T, tv%S, g, gv, us)
2430 if (
associated(tv%T))
then
2433 call hchksum(ea,
"before triDiagTS ea ", g%HI, haloshift=0, scale=gv%H_to_m)
2434 call hchksum(eb,
"before triDiagTS eb ", g%HI, haloshift=0, scale=gv%H_to_m)
2436 call cpu_clock_begin(id_clock_tridiag)
2444 if (
associated(tv%S) .and.
associated(tv%salt_deficit)) &
2445 call adjust_salt(h, tv, g, gv, cs%diabatic_aux_CSp)
2447 if (cs%diabatic_diff_tendency_diag)
then
2448 do k=1,nz ;
do j=js,je ;
do i=is,ie
2449 temp_diag(i,j,k) = tv%T(i,j,k)
2450 saln_diag(i,j,k) = tv%S(i,j,k)
2451 enddo ;
enddo ;
enddo
2455 if (cs%tracer_tridiag)
then
2456 call tracer_vertdiff(hold, ea, eb, dt, tv%T, g, gv)
2457 call tracer_vertdiff(hold, ea, eb, dt, tv%S, g, gv)
2459 call tridiagts(g, gv, is, ie, js, je, hold, ea, eb, tv%T, tv%S)
2465 if (cs%diabatic_diff_tendency_diag)
then
2466 call diagnose_diabatic_diff_tendency(tv, hold, temp_diag, saln_diag, dt, g, gv, us, cs)
2467 if (cs%id_diabatic_diff_h > 0)
call post_data(cs%id_diabatic_diff_h, hold, cs%diag, alt_h=hold)
2470 call cpu_clock_end(id_clock_tridiag)
2471 if (showcalltree)
call calltree_waypoint(
"done with triDiagTS (diabatic)")
2474 if (cs%debugConservation)
call mom_state_stats(
'triDiagTS', u, v, h, tv%T, tv%S, g, gv, us)
2479 call mom_state_chksum(
"after mixed layer ", u, v, h, g, gv, us, haloshift=0)
2480 call mom_thermovar_chksum(
"after mixed layer ", tv, g)
2481 call hchksum(ea,
"after mixed layer ea", g%HI, scale=gv%H_to_m)
2482 call hchksum(eb,
"after mixed layer eb", g%HI, scale=gv%H_to_m)
2485 call cpu_clock_begin(id_clock_remap)
2486 call regularize_layers(h, tv, dt, ea, eb, g, gv, us, cs%regularize_layers_CSp)
2487 call cpu_clock_end(id_clock_remap)
2488 if (showcalltree)
call calltree_waypoint(
"done with regularize_layers (diabatic)")
2489 if (cs%debugConservation)
call mom_state_stats(
'regularize_layers', u, v, h, tv%T, tv%S, g, gv, us)
2493 if (
associated(adp%du_dt_dia) .or.
associated(adp%dv_dt_dia)) &
2495 call pass_var(h, g%domain, to_west+to_south+omit_corners, halo=1)
2496 call diag_update_remap_grids(cs%diag)
2500 if ((cs%id_Tdif > 0) .or. (cs%id_Tadv > 0))
then
2501 do j=js,je ;
do i=is,ie
2502 tdif_flx(i,j,1) = 0.0 ; tdif_flx(i,j,nz+1) = 0.0
2503 tadv_flx(i,j,1) = 0.0 ; tadv_flx(i,j,nz+1) = 0.0
2506 do k=2,nz ;
do j=js,je ;
do i=is,ie
2507 tdif_flx(i,j,k) = (idt * 0.5*(ea(i,j,k) + eb(i,j,k-1))) * &
2508 (tv%T(i,j,k-1) - tv%T(i,j,k))
2509 tadv_flx(i,j,k) = (idt * (ea(i,j,k) - eb(i,j,k-1))) * &
2510 0.5*(tv%T(i,j,k-1) + tv%T(i,j,k))
2511 enddo ;
enddo ;
enddo
2513 if ((cs%id_Sdif > 0) .or. (cs%id_Sadv > 0))
then
2514 do j=js,je ;
do i=is,ie
2515 sdif_flx(i,j,1) = 0.0 ; sdif_flx(i,j,nz+1) = 0.0
2516 sadv_flx(i,j,1) = 0.0 ; sadv_flx(i,j,nz+1) = 0.0
2519 do k=2,nz ;
do j=js,je ;
do i=is,ie
2520 sdif_flx(i,j,k) = (idt * 0.5*(ea(i,j,k) + eb(i,j,k-1))) * &
2521 (tv%S(i,j,k-1) - tv%S(i,j,k))
2522 sadv_flx(i,j,k) = (idt * (ea(i,j,k) - eb(i,j,k-1))) * &
2523 0.5*(tv%S(i,j,k-1) + tv%S(i,j,k))
2524 enddo ;
enddo ;
enddo
2528 call cpu_clock_begin(id_clock_tracers)
2529 if (cs%mix_boundary_tracers)
then
2530 tr_ea_bbl = gv%Z_to_H * sqrt(dt*cs%Kd_BBL_tr)
2534 ebtr(i,j,nz) = eb(i,j,nz)
2536 in_boundary(i) = (g%mask2dT(i,j) > 0.0)
2538 do k=nz,2,-1 ;
do i=is,ie
2539 if (in_boundary(i))
then
2540 htot(i) = htot(i) + h(i,j,k)
2549 add_ent = ((dt * cs%Kd_min_tr) * gv%Z_to_H**2) * &
2550 ((h(i,j,k-1)+h(i,j,k)+h_neglect) / &
2551 (h(i,j,k-1)*h(i,j,k)+h_neglect2)) - &
2552 0.5*(ea(i,j,k) + eb(i,j,k-1))
2553 if (htot(i) < tr_ea_bbl)
then
2554 add_ent = max(0.0, add_ent, &
2555 (tr_ea_bbl - htot(i)) - min(ea(i,j,k), eb(i,j,k-1)))
2556 elseif (add_ent < 0.0)
then
2557 add_ent = 0.0 ; in_boundary(i) = .false.
2560 ebtr(i,j,k-1) = eb(i,j,k-1) + add_ent
2561 eatr(i,j,k) = ea(i,j,k) + add_ent
2563 ebtr(i,j,k-1) = eb(i,j,k-1) ; eatr(i,j,k) = ea(i,j,k)
2565 if (cs%double_diffuse)
then ;
if (kd_extra_s(i,j,k) > 0.0)
then
2566 add_ent = ((dt * kd_extra_s(i,j,k)) * gv%Z_to_H**2) / &
2567 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + &
2569 ebtr(i,j,k-1) = ebtr(i,j,k-1) + add_ent
2570 eatr(i,j,k) = eatr(i,j,k) + add_ent
2573 do i=is,ie ; eatr(i,j,1) = ea(i,j,1) ;
enddo
2577 call call_tracer_column_fns(hold, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
2578 cs%optics, cs%tracer_flow_CSp, cs%debug)
2580 elseif (cs%double_diffuse)
then
2582 do j=js,je ;
do i=is,ie
2583 ebtr(i,j,nz) = eb(i,j,nz) ; eatr(i,j,1) = ea(i,j,1)
2586 do k=nz,2,-1 ;
do j=js,je ;
do i=is,ie
2587 if (kd_extra_s(i,j,k) > 0.0)
then
2588 add_ent = ((dt * kd_extra_s(i,j,k)) * gv%Z_to_H**2) / &
2589 (0.25 * ((h(i,j,k-1) + h(i,j,k)) + (hold(i,j,k-1) + hold(i,j,k))) + &
2594 ebtr(i,j,k-1) = eb(i,j,k-1) + add_ent
2595 eatr(i,j,k) = ea(i,j,k) + add_ent
2596 enddo ;
enddo ;
enddo
2598 call call_tracer_column_fns(hold, h, eatr, ebtr, fluxes, hml, dt, g, gv, us, tv, &
2599 cs%optics, cs%tracer_flow_CSp, cs%debug)
2602 call call_tracer_column_fns(hold, h, ea, eb, fluxes, hml, dt, g, gv, us, tv, &
2603 cs%optics, cs%tracer_flow_CSp, cs%debug)
2607 call cpu_clock_end(id_clock_tracers)
2610 if (cs%use_sponge)
then
2611 call cpu_clock_begin(id_clock_sponge)
2613 if (cs%bulkmixedlayer .and.
associated(tv%eqn_of_state))
then
2614 do i=is,ie ; p_ref_cv(i) = tv%P_Ref ;
enddo
2615 eosdom(:) = eos_domain(g%HI)
2619 tv%eqn_of_state, eosdom)
2621 call apply_sponge(h, dt, g, gv, us, ea, eb, cs%sponge_CSp, rcv_ml)
2623 call apply_sponge(h, dt, g, gv, us, ea, eb, cs%sponge_CSp)
2625 call cpu_clock_end(id_clock_sponge)
2628 call mom_thermovar_chksum(
"apply_sponge ", tv, g)
2633 if (
associated(cdp%diapyc_vel))
then
2636 do k=2,nz ;
do i=is,ie
2637 cdp%diapyc_vel(i,j,k) = us%s_to_T*idt * (ea(i,j,k) - eb(i,j,k-1))
2640 cdp%diapyc_vel(i,j,1) = 0.0
2641 cdp%diapyc_vel(i,j,nz+1) = 0.0
2649 if (cs%bulkmixedlayer)
then
2651 call hchksum(ea,
"before net flux rearrangement ea", g%HI, scale=gv%H_to_m)
2652 call hchksum(eb,
"before net flux rearrangement eb", g%HI, scale=gv%H_to_m)
2656 do k=2,gv%nkml ;
do i=is,ie
2657 net_ent = ea(i,j,k) - eb(i,j,k-1)
2658 ea(i,j,k) = max(net_ent, 0.0)
2659 eb(i,j,k-1) = max(-net_ent, 0.0)
2663 call hchksum(ea,
"after net flux rearrangement ea", g%HI, scale=gv%H_to_m)
2664 call hchksum(eb,
"after net flux rearrangement eb", g%HI, scale=gv%H_to_m)
2672 hold(i,js-1,k) = gv%Angstrom_H ; ea(i,js-1,k) = 0.0 ; eb(i,js-1,k) = 0.0
2673 hold(i,je+1,k) = gv%Angstrom_H ; ea(i,je+1,k) = 0.0 ; eb(i,je+1,k) = 0.0
2676 hold(is-1,j,k) = gv%Angstrom_H ; ea(is-1,j,k) = 0.0 ; eb(is-1,j,k) = 0.0
2677 hold(ie+1,j,k) = gv%Angstrom_H ; ea(ie+1,j,k) = 0.0 ; eb(ie+1,j,k) = 0.0
2681 call cpu_clock_begin(id_clock_pass)
2682 if (g%symmetric)
then ; dir_flag = to_all+omit_corners
2683 else ; dir_flag = to_west+to_south+omit_corners ;
endif
2687 call do_group_pass(cs%pass_hold_eb_ea, g%Domain)
2688 call cpu_clock_end(id_clock_pass)
2694 call mom_state_chksum(
"before u/v tridiag ", u, v, h, g, gv, us, haloshift=0)
2695 call hchksum(ea,
"before u/v tridiag ea", g%HI, scale=gv%H_to_m)
2696 call hchksum(eb,
"before u/v tridiag eb", g%HI, scale=gv%H_to_m)
2697 call hchksum(hold,
"before u/v tridiag hold", g%HI, scale=gv%H_to_m)
2699 call cpu_clock_begin(id_clock_tridiag)
2700 idt_accel = 1.0 / dt
2704 if (
associated(adp%du_dt_dia)) adp%du_dt_dia(i,j,1) = u(i,j,1)
2705 hval = (hold(i,j,1) + hold(i+1,j,1)) + (ea(i,j,1) + ea(i+1,j,1)) + h_neglect
2706 b1(i) = 1.0 / (hval + (eb(i,j,1) + eb(i+1,j,1)))
2707 d1(i) = hval * b1(i)
2708 u(i,j,1) = b1(i) * (hval * u(i,j,1))
2710 do k=2,nz ;
do i=isq,ieq
2711 if (
associated(adp%du_dt_dia)) adp%du_dt_dia(i,j,k) = u(i,j,k)
2712 c1(i,k) = (eb(i,j,k-1)+eb(i+1,j,k-1)) * b1(i)
2713 eaval = ea(i,j,k) + ea(i+1,j,k)
2714 hval = hold(i,j,k) + hold(i+1,j,k) + h_neglect
2715 b1(i) = 1.0 / ((eb(i,j,k) + eb(i+1,j,k)) + (hval + d1(i)*eaval))
2716 d1(i) = (hval + d1(i)*eaval) * b1(i)
2717 u(i,j,k) = (hval*u(i,j,k) + eaval*u(i,j,k-1))*b1(i)
2719 do k=nz-1,1,-1 ;
do i=isq,ieq
2720 u(i,j,k) = u(i,j,k) + c1(i,k+1)*u(i,j,k+1)
2721 if (
associated(adp%du_dt_dia)) &
2722 adp%du_dt_dia(i,j,k) = (u(i,j,k) - adp%du_dt_dia(i,j,k)) * idt_accel
2724 if (
associated(adp%du_dt_dia))
then
2726 adp%du_dt_dia(i,j,nz) = (u(i,j,nz)-adp%du_dt_dia(i,j,nz)) * idt_accel
2731 call mom_state_chksum(
"aft 1st loop tridiag ", u, v, h, g, gv, us, haloshift=0)
2736 if (
associated(adp%dv_dt_dia)) adp%dv_dt_dia(i,j,1) = v(i,j,1)
2737 hval = (hold(i,j,1) + hold(i,j+1,1)) + (ea(i,j,1) + ea(i,j+1,1)) + h_neglect
2738 b1(i) = 1.0 / (hval + (eb(i,j,1) + eb(i,j+1,1)))
2739 d1(i) = hval * b1(i)
2740 v(i,j,1) = b1(i) * (hval * v(i,j,1))
2742 do k=2,nz ;
do i=is,ie
2743 if (
associated(adp%dv_dt_dia)) adp%dv_dt_dia(i,j,k) = v(i,j,k)
2744 c1(i,k) = (eb(i,j,k-1)+eb(i,j+1,k-1)) * b1(i)
2745 eaval = ea(i,j,k) + ea(i,j+1,k)
2746 hval = hold(i,j,k) + hold(i,j+1,k) + h_neglect
2747 b1(i) = 1.0 / ((eb(i,j,k) + eb(i,j+1,k)) + (hval + d1(i)*eaval))
2748 d1(i) = (hval + d1(i)*eaval) * b1(i)
2749 v(i,j,k) = (hval*v(i,j,k) + eaval*v(i,j,k-1))*b1(i)
2751 do k=nz-1,1,-1 ;
do i=is,ie
2752 v(i,j,k) = v(i,j,k) + c1(i,k+1)*v(i,j,k+1)
2753 if (
associated(adp%dv_dt_dia)) &
2754 adp%dv_dt_dia(i,j,k) = (v(i,j,k) - adp%dv_dt_dia(i,j,k)) * idt_accel
2756 if (
associated(adp%dv_dt_dia))
then
2758 adp%dv_dt_dia(i,j,nz) = (v(i,j,nz)-adp%dv_dt_dia(i,j,nz)) * idt_accel
2762 call cpu_clock_end(id_clock_tridiag)
2764 call mom_state_chksum(
"after u/v tridiag ", u, v, h, g, gv, us, haloshift=0)
2767 call disable_averaging(cs%diag)
2769 call enable_averages(dt, time_end, cs%diag)
2771 if (cs%id_Kd_interface > 0)
call post_data(cs%id_Kd_interface, kd_int, cs%diag)
2772 if (cs%id_Kd_heat > 0)
call post_data(cs%id_Kd_heat, kd_heat, cs%diag)
2773 if (cs%id_Kd_salt > 0)
call post_data(cs%id_Kd_salt, kd_salt, cs%diag)
2774 if (cs%id_Kd_ePBL > 0)
call post_data(cs%id_Kd_ePBL, kd_epbl, cs%diag)
2776 if (cs%id_ea > 0)
call post_data(cs%id_ea, ea, cs%diag)
2777 if (cs%id_eb > 0)
call post_data(cs%id_eb, eb, cs%diag)
2779 if (cs%id_dudt_dia > 0)
call post_data(cs%id_dudt_dia, adp%du_dt_dia, cs%diag)
2780 if (cs%id_dvdt_dia > 0)
call post_data(cs%id_dvdt_dia, adp%dv_dt_dia, cs%diag)
2781 if (cs%id_wd > 0)
call post_data(cs%id_wd, cdp%diapyc_vel, cs%diag)
2783 if (cs%id_Tdif > 0)
call post_data(cs%id_Tdif, tdif_flx, cs%diag)
2784 if (cs%id_Tadv > 0)
call post_data(cs%id_Tadv, tadv_flx, cs%diag)
2785 if (cs%id_Sdif > 0)
call post_data(cs%id_Sdif, sdif_flx, cs%diag)
2786 if (cs%id_Sadv > 0)
call post_data(cs%id_Sadv, sadv_flx, cs%diag)
2789 if (cs%id_hf_dudt_dia_2d > 0)
then
2790 allocate(hf_dudt_dia_2d(g%IsdB:g%IedB,g%jsd:g%jed))
2791 hf_dudt_dia_2d(:,:) = 0.0
2792 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
2793 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)
2794 enddo ;
enddo ;
enddo
2795 call post_data(cs%id_hf_dudt_dia_2d, hf_dudt_dia_2d, cs%diag)
2796 deallocate(hf_dudt_dia_2d)
2799 if (cs%id_hf_dvdt_dia_2d > 0)
then
2800 allocate(hf_dvdt_dia_2d(g%isd:g%ied,g%JsdB:g%JedB))
2801 hf_dvdt_dia_2d(:,:) = 0.0
2802 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
2803 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)
2804 enddo ;
enddo ;
enddo
2805 call post_data(cs%id_hf_dvdt_dia_2d, hf_dvdt_dia_2d, cs%diag)
2806 deallocate(hf_dvdt_dia_2d)
2809 call disable_averaging(cs%diag)
2811 if (showcalltree)
call calltree_leave(
"layered_diabatic()")
2813 end subroutine layered_diabatic
2817 subroutine extract_diabatic_member(CS, opacity_CSp, optics_CSp, evap_CFL_limit, minimum_forcing_depth, &
2818 KPP_CSp, energetic_PBL_CSp, diabatic_aux_CSp, diabatic_halo)
2821 type(
opacity_cs),
optional,
pointer :: opacity_csp
2822 type(
optics_type),
optional,
pointer :: optics_csp
2823 type(
kpp_cs),
optional,
pointer :: kpp_csp
2825 real,
optional,
intent( out) :: evap_cfl_limit
2827 real,
optional,
intent( out) :: minimum_forcing_depth
2831 integer,
optional,
intent( out) :: diabatic_halo
2835 if (
present(opacity_csp)) opacity_csp => cs%opacity_CSp
2836 if (
present(optics_csp)) optics_csp => cs%optics
2837 if (
present(kpp_csp)) kpp_csp => cs%KPP_CSp
2838 if (
present(energetic_pbl_csp)) energetic_pbl_csp => cs%energetic_PBL_CSp
2841 if (
present(evap_cfl_limit)) evap_cfl_limit = cs%evap_CFL_limit
2842 if (
present(minimum_forcing_depth)) minimum_forcing_depth = cs%minimum_forcing_depth
2843 if (
present(diabatic_halo)) diabatic_halo = cs%halo_TS_diff
2845 end subroutine extract_diabatic_member
2848 subroutine adiabatic(h, tv, fluxes, dt, G, GV, US, CS)
2850 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2853 type(
forcing),
intent(inout) :: fluxes
2854 real,
intent(in) :: dt
2859 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: zeros
2863 call call_tracer_column_fns(h, h, zeros, zeros, fluxes, zeros(:,:,1), dt, g, gv, us, tv, &
2864 cs%optics, cs%tracer_flow_CSp, cs%debug)
2866 end subroutine adiabatic
2872 subroutine diagnose_diabatic_diff_tendency(tv, h, temp_old, saln_old, dt, G, GV, US, CS)
2876 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
2877 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: temp_old
2878 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: saln_old
2879 real,
intent(in) :: dt
2884 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: work_3d
2885 real,
dimension(SZI_(G),SZJ_(G)) :: work_2d
2887 real :: ppt2mks = 0.001
2888 integer :: i, j, k, is, ie, js, je, nz
2889 logical :: do_saln_tend
2891 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
2892 idt = 0.0 ;
if (dt > 0.0) idt = 1. / dt
2893 work_3d(:,:,:) = 0.0
2898 do k=1,nz ;
do j=js,je ;
do i=is,ie
2899 work_3d(i,j,k) = (tv%T(i,j,k)-temp_old(i,j,k))*idt
2900 enddo ;
enddo ;
enddo
2901 if (cs%id_diabatic_diff_temp_tend > 0)
then
2902 call post_data(cs%id_diabatic_diff_temp_tend, work_3d, cs%diag, alt_h=h)
2906 if (cs%id_diabatic_diff_heat_tend > 0 .or. cs%id_diabatic_diff_heat_tend_2d > 0)
then
2907 do k=1,nz ;
do j=js,je ;
do i=is,ie
2908 work_3d(i,j,k) = h(i,j,k)*gv%H_to_RZ * tv%C_p * work_3d(i,j,k)
2909 enddo ;
enddo ;
enddo
2910 if (cs%id_diabatic_diff_heat_tend > 0)
then
2911 call post_data(cs%id_diabatic_diff_heat_tend, work_3d, cs%diag, alt_h=h)
2913 if (cs%id_diabatic_diff_heat_tend_2d > 0)
then
2914 do j=js,je ;
do i=is,ie
2917 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
2920 call post_data(cs%id_diabatic_diff_heat_tend_2d, work_2d, cs%diag)
2925 do_saln_tend = cs%id_diabatic_diff_saln_tend > 0 &
2926 .or. cs%id_diabatic_diff_salt_tend > 0 &
2927 .or. cs%id_diabatic_diff_salt_tend_2d > 0
2929 if (do_saln_tend)
then
2930 do k=1,nz ;
do j=js,je ;
do i=is,ie
2931 work_3d(i,j,k) = (tv%S(i,j,k) - saln_old(i,j,k)) * idt
2932 enddo ;
enddo ;
enddo
2934 if (cs%id_diabatic_diff_saln_tend > 0) &
2935 call post_data(cs%id_diabatic_diff_saln_tend, work_3d, cs%diag, alt_h=h)
2938 if (cs%id_diabatic_diff_salt_tend > 0 .or. cs%id_diabatic_diff_salt_tend_2d > 0)
then
2939 do k=1,nz ;
do j=js,je ;
do i=is,ie
2940 work_3d(i,j,k) = h(i,j,k)*gv%H_to_RZ * ppt2mks * work_3d(i,j,k)
2941 enddo ;
enddo ;
enddo
2942 if (cs%id_diabatic_diff_salt_tend > 0)
then
2943 call post_data(cs%id_diabatic_diff_salt_tend, work_3d, cs%diag, alt_h=h)
2945 if (cs%id_diabatic_diff_salt_tend_2d > 0)
then
2946 do j=js,je ;
do i=is,ie
2949 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
2952 call post_data(cs%id_diabatic_diff_salt_tend_2d, work_2d, cs%diag)
2957 end subroutine diagnose_diabatic_diff_tendency
2964 subroutine diagnose_boundary_forcing_tendency(tv, h, temp_old, saln_old, h_old, &
2969 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2971 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2972 intent(in) :: temp_old
2973 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2974 intent(in) :: saln_old
2975 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
2977 real,
intent(in) :: dt
2982 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: work_3d
2983 real,
dimension(SZI_(G),SZJ_(G)) :: work_2d
2985 real :: ppt2mks = 0.001
2986 integer :: i, j, k, is, ie, js, je, nz
2988 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
2989 idt = 0.0 ;
if (dt > 0.0) idt = 1. / dt
2990 work_3d(:,:,:) = 0.0
2994 if (cs%id_boundary_forcing_h_tendency > 0)
then
2995 do k=1,nz ;
do j=js,je ;
do i=is,ie
2996 work_3d(i,j,k) = (h(i,j,k) - h_old(i,j,k))*idt
2997 enddo ;
enddo ;
enddo
2998 call post_data(cs%id_boundary_forcing_h_tendency, work_3d, cs%diag, alt_h=h_old)
3002 if (cs%id_boundary_forcing_temp_tend > 0)
then
3003 do k=1,nz ;
do j=js,je ;
do i=is,ie
3004 work_3d(i,j,k) = (tv%T(i,j,k)-temp_old(i,j,k))*idt
3005 enddo ;
enddo ;
enddo
3006 call post_data(cs%id_boundary_forcing_temp_tend, work_3d, cs%diag, alt_h=h_old)
3010 if (cs%id_boundary_forcing_heat_tend > 0 .or. cs%id_boundary_forcing_heat_tend_2d > 0)
then
3011 do k=1,nz ;
do j=js,je ;
do i=is,ie
3012 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))
3013 enddo ;
enddo ;
enddo
3014 if (cs%id_boundary_forcing_heat_tend > 0)
then
3015 call post_data(cs%id_boundary_forcing_heat_tend, work_3d, cs%diag, alt_h = h_old)
3017 if (cs%id_boundary_forcing_heat_tend_2d > 0)
then
3018 do j=js,je ;
do i=is,ie
3021 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
3024 call post_data(cs%id_boundary_forcing_heat_tend_2d, work_2d, cs%diag)
3029 if (cs%id_boundary_forcing_saln_tend > 0)
then
3030 do k=1,nz ;
do j=js,je ;
do i=is,ie
3031 work_3d(i,j,k) = (tv%S(i,j,k)-saln_old(i,j,k))*idt
3032 enddo ;
enddo ;
enddo
3033 call post_data(cs%id_boundary_forcing_saln_tend, work_3d, cs%diag, alt_h = h_old)
3037 if (cs%id_boundary_forcing_salt_tend > 0 .or. cs%id_boundary_forcing_salt_tend_2d > 0)
then
3038 do k=1,nz ;
do j=js,je ;
do i=is,ie
3039 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))
3040 enddo ;
enddo ;
enddo
3041 if (cs%id_boundary_forcing_salt_tend > 0)
then
3042 call post_data(cs%id_boundary_forcing_salt_tend, work_3d, cs%diag, alt_h = h_old)
3044 if (cs%id_boundary_forcing_salt_tend_2d > 0)
then
3045 do j=js,je ;
do i=is,ie
3048 work_2d(i,j) = work_2d(i,j) + work_3d(i,j,k)
3051 call post_data(cs%id_boundary_forcing_salt_tend_2d, work_2d, cs%diag)
3055 end subroutine diagnose_boundary_forcing_tendency
3062 subroutine diagnose_frazil_tendency(tv, h, temp_old, dt, G, GV, US, CS)
3066 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
3067 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: temp_old
3068 real,
intent(in) :: dt
3072 real,
dimension(SZI_(G),SZJ_(G)) :: work_2d
3074 integer :: i, j, k, is, ie, js, je, nz
3076 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
3077 idt = 0.0 ;
if (dt > 0.0) idt = 1. / dt
3080 if (cs%id_frazil_temp_tend > 0)
then
3081 do k=1,nz ;
do j=js,je ;
do i=is,ie
3082 cs%frazil_temp_diag(i,j,k) = idt * (tv%T(i,j,k)-temp_old(i,j,k))
3083 enddo ;
enddo ;
enddo
3084 call post_data(cs%id_frazil_temp_tend, cs%frazil_temp_diag(:,:,:), cs%diag)
3088 if (cs%id_frazil_heat_tend > 0 .or. cs%id_frazil_heat_tend_2d > 0)
then
3089 do k=1,nz ;
do j=js,je ;
do i=is,ie
3090 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))
3091 enddo ;
enddo ;
enddo
3092 if (cs%id_frazil_heat_tend > 0)
call post_data(cs%id_frazil_heat_tend, cs%frazil_heat_diag(:,:,:), cs%diag)
3096 if (cs%id_frazil_heat_tend_2d > 0)
then
3097 do j=js,je ;
do i=is,ie
3100 work_2d(i,j) = work_2d(i,j) + cs%frazil_heat_diag(i,j,k)
3103 call post_data(cs%id_frazil_heat_tend_2d, work_2d, cs%diag)
3107 end subroutine diagnose_frazil_tendency
3113 subroutine adiabatic_driver_init(Time, G, param_file, diag, CS, &
3115 type(time_type),
intent(in) :: time
3118 type(
diag_ctrl),
target,
intent(inout) :: diag
3124 #include "version_variable.h"
3125 character(len=40) :: mdl =
"MOM_diabatic_driver"
3127 if (
associated(cs))
then
3128 call mom_error(warning,
"adiabatic_driver_init called with an "// &
3129 "associated control structure.")
3131 else ;
allocate(cs) ;
endif
3134 if (
associated(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
3138 "The following parameters are used for diabatic processes.")
3140 end subroutine adiabatic_driver_init
3144 subroutine diabatic_driver_init(Time, G, GV, US, param_file, useALEalgorithm, diag, &
3145 ADp, CDp, CS, tracer_flow_CSp, sponge_CSp, &
3147 type(time_type),
target :: time
3152 logical,
intent(in) :: usealealgorithm
3153 type(
diag_ctrl),
target,
intent(inout) :: diag
3165 logical :: use_temperature
3166 character(len=20) :: en1, en2, en3
3169 #include "version_variable.h"
3170 character(len=40) :: mdl =
"MOM_diabatic_driver"
3171 character(len=48) :: thickness_units
3172 character(len=40) :: var_name
3173 character(len=160) :: var_descript
3174 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, nz, nbands, m
3175 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
3176 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
3178 if (
associated(cs))
then
3179 call mom_error(warning,
"diabatic_driver_init called with an "// &
3180 "associated control structure.")
3189 if (
associated(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
3190 if (
associated(sponge_csp)) cs%sponge_CSp => sponge_csp
3191 if (
associated(ale_sponge_csp)) cs%ALE_sponge_CSp => ale_sponge_csp
3193 cs%useALEalgorithm = usealealgorithm
3194 cs%bulkmixedlayer = (gv%nkml > 0)
3198 "The following parameters are used for diabatic processes.", &
3199 log_to_all=.true., debugging=.true.)
3200 call get_param(param_file, mdl,
"USE_LEGACY_DIABATIC_DRIVER", cs%use_legacy_diabatic, &
3201 "If true, use a legacy version of the diabatic subroutine. "//&
3202 "This is temporary and is needed to avoid change in answers.", &
3204 call get_param(param_file, mdl,
"SPONGE", cs%use_sponge, &
3205 "If true, sponges may be applied anywhere in the domain. "//&
3206 "The exact location and properties of those sponges are "//&
3207 "specified via calls to initialize_sponge and possibly "//&
3208 "set_up_sponge_field.", default=.false.)
3209 call get_param(param_file, mdl,
"ENABLE_THERMODYNAMICS", use_temperature, &
3210 "If true, temperature and salinity are used as state "//&
3211 "variables.", default=.true.)
3212 call get_param(param_file, mdl,
"ENERGETICS_SFC_PBL", cs%use_energetic_PBL, &
3213 "If true, use an implied energetics planetary boundary "//&
3214 "layer scheme to determine the diffusivity and viscosity "//&
3215 "in the surface boundary layer.", default=.false.)
3216 call get_param(param_file, mdl,
"EPBL_IS_ADDITIVE", cs%ePBL_is_additive, &
3217 "If true, the diffusivity from ePBL is added to all "//&
3218 "other diffusivities. Otherwise, the larger of kappa-shear "//&
3219 "and ePBL diffusivities are used.", default=.true.)
3220 call get_param(param_file, mdl,
"PRANDTL_EPBL", cs%ePBL_Prandtl, &
3221 "The Prandtl number used by ePBL to convert vertical diffusivities into "//&
3222 "viscosities.", default=1.0, units=
"nondim", do_not_log=.not.cs%use_energetic_PBL)
3223 call get_param(param_file, mdl,
"USE_KPP", cs%use_KPP, &
3224 "If true, turns on the [CVMix] KPP scheme of Large et al., 1994, "//&
3225 "to calculate diffusivities and non-local transport in the OBL.", &
3226 default=.false., do_not_log=.true.)
3227 cs%use_CVMix_ddiff = cvmix_ddiff_is_used(param_file)
3229 cs%use_kappa_shear = kappa_shear_is_used(param_file)
3230 cs%use_CVMix_shear = cvmix_shear_is_used(param_file)
3232 if (cs%bulkmixedlayer)
then
3233 call get_param(param_file, mdl,
"ML_MIX_FIRST", cs%ML_mix_first, &
3234 "The fraction of the mixed layer mixing that is applied "//&
3235 "before interior diapycnal mixing. 0 by default.", &
3236 units=
"nondim", default=0.0)
3237 call get_param(param_file, mdl,
"NKBL", cs%nkbl, default=2, do_not_log=.true.)
3239 cs%ML_mix_first = 0.0
3241 if (use_temperature)
then
3242 call get_param(param_file, mdl,
"DO_GEOTHERMAL", cs%use_geothermal, &
3243 "If true, apply geothermal heating.", default=.false.)
3245 cs%use_geothermal = .false.
3247 call get_param(param_file, mdl,
"INTERNAL_TIDES", cs%use_int_tides, &
3248 "If true, use the code that advances a separate set of "//&
3249 "equations for the internal tide energy density.", default=.false.)
3251 if (cs%use_int_tides)
then
3252 call get_param(param_file, mdl,
"INTERNAL_TIDE_MODES", cs%nMode, &
3253 "The number of distinct internal tide modes "//&
3254 "that will be calculated.", default=1, do_not_log=.true.)
3255 call get_param(param_file, mdl,
"UNIFORM_TEST_CG", cs%uniform_test_cg, &
3256 "If positive, a uniform group velocity of internal tide for test case", &
3257 default=-1., units=
"m s-1", scale=us%m_s_to_L_T)
3260 call get_param(param_file, mdl,
"MASSLESS_MATCH_TARGETS", &
3261 cs%massless_match_targets, &
3262 "If true, the temperature and salinity of massless layers "//&
3263 "are kept consistent with their target densities. "//&
3264 "Otherwise the properties of massless layers evolve "//&
3265 "diffusively to match massive neighboring layers.", &
3268 call get_param(param_file, mdl,
"AGGREGATE_FW_FORCING", cs%aggregate_FW_forcing, &
3269 "If true, the net incoming and outgoing fresh water fluxes are combined "//&
3270 "and applied as either incoming or outgoing depending on the sign of the net. "//&
3271 "If false, the net incoming fresh water flux is added to the model and "//&
3272 "thereafter the net outgoing is removed from the topmost non-vanished "//&
3273 "layers of the updated state.", default=.true.)
3275 call get_param(param_file, mdl,
"DEBUG", cs%debug, &
3276 "If true, write out verbose debugging data.", &
3277 default=.false., debuggingparam=.true.)
3278 call get_param(param_file, mdl,
"DEBUG_CONSERVATION", cs%debugConservation, &
3279 "If true, monitor conservation and extrema.", &
3280 default=.false., debuggingparam=.true.)
3282 call get_param(param_file, mdl,
"DEBUG_ENERGY_REQ", cs%debug_energy_req, &
3283 "If true, debug the energy requirements.", default=.false., do_not_log=.true.)
3284 call get_param(param_file, mdl,
"MIX_BOUNDARY_TRACERS", cs%mix_boundary_tracers, &
3285 "If true, mix the passive tracers in massless layers at "//&
3286 "the bottom into the interior as though a diffusivity of "//&
3287 "KD_MIN_TR were operating.", default=.true.)
3289 if (cs%mix_boundary_tracers)
then
3290 call get_param(param_file, mdl,
"KD", kd, default=0.0)
3291 call get_param(param_file, mdl,
"KD_MIN_TR", cs%Kd_min_tr, &
3292 "A minimal diffusivity that should always be applied to "//&
3293 "tracers, especially in massless layers near the bottom. "//&
3294 "The default is 0.1*KD.", units=
"m2 s-1", default=0.1*kd, scale=us%m2_s_to_Z2_T)
3295 call get_param(param_file, mdl,
"KD_BBL_TR", cs%Kd_BBL_tr, &
3296 "A bottom boundary layer tracer diffusivity that will "//&
3297 "allow for explicitly specified bottom fluxes. The "//&
3298 "entrainment at the bottom is at least sqrt(Kd_BBL_tr*dt) "//&
3299 "over the same distance.", units=
"m2 s-1", default=0., scale=us%m2_s_to_Z2_T)
3302 call get_param(param_file, mdl,
"TRACER_TRIDIAG", cs%tracer_tridiag, &
3303 "If true, use the passive tracer tridiagonal solver for T and S", &
3306 call get_param(param_file, mdl,
"MINIMUM_FORCING_DEPTH", cs%minimum_forcing_depth, &
3307 "The smallest depth over which forcing can be applied. This "//&
3308 "only takes effect when near-surface layers become thin "//&
3309 "relative to this scale, in which case the forcing tendencies "//&
3310 "scaled down by distributing the forcing over this depth scale.", &
3311 units=
"m", default=0.001, scale=gv%m_to_H)
3312 call get_param(param_file, mdl,
"EVAP_CFL_LIMIT", cs%evap_CFL_limit, &
3313 "The largest fraction of a layer than can be lost to forcing "//&
3314 "(e.g. evaporation, sea-ice formation) in one time-step. The unused "//&
3315 "mass loss is passed down through the column.", &
3316 units=
"nondim", default=0.8)
3320 thickness_units = get_thickness_units(gv)
3322 cs%id_ea_t = register_diag_field(
'ocean_model',
'ea_t', diag%axesTL, time, &
3323 'Layer (heat) entrainment from above per timestep',
'm', conversion=gv%H_to_m)
3324 cs%id_eb_t = register_diag_field(
'ocean_model',
'eb_t', diag%axesTL, time, &
3325 'Layer (heat) entrainment from below per timestep',
'm', conversion=gv%H_to_m)
3326 cs%id_ea_s = register_diag_field(
'ocean_model',
'ea_s', diag%axesTL, time, &
3327 'Layer (salt) entrainment from above per timestep',
'm', conversion=gv%H_to_m)
3328 cs%id_eb_s = register_diag_field(
'ocean_model',
'eb_s', diag%axesTL, time, &
3329 'Layer (salt) entrainment from below per timestep',
'm', conversion=gv%H_to_m)
3331 cs%id_ea = register_diag_field(
'ocean_model',
'ea', diag%axesTL, time, &
3332 'Layer entrainment from above per timestep',
'm', conversion=gv%H_to_m)
3333 cs%id_eb = register_diag_field(
'ocean_model',
'eb', diag%axesTL, time, &
3334 'Layer entrainment from below per timestep',
'm', conversion=gv%H_to_m)
3335 cs%id_wd = register_diag_field(
'ocean_model',
'wd', diag%axesTi, time, &
3336 'Diapycnal velocity',
'm s-1', conversion=gv%H_to_m)
3337 if (cs%id_wd > 0)
call safe_alloc_ptr(cdp%diapyc_vel,isd,ied,jsd,jed,nz+1)
3339 cs%id_dudt_dia = register_diag_field(
'ocean_model',
'dudt_dia', diag%axesCuL, time, &
3340 'Zonal Acceleration from Diapycnal Mixing',
'm s-2', conversion=us%L_T2_to_m_s2)
3341 cs%id_dvdt_dia = register_diag_field(
'ocean_model',
'dvdt_dia', diag%axesCvL, time, &
3342 'Meridional Acceleration from Diapycnal Mixing',
'm s-2', conversion=us%L_T2_to_m_s2)
3344 cs%id_hf_dudt_dia_2d = register_diag_field(
'ocean_model',
'hf_dudt_dia_2d', diag%axesCu1, time, &
3345 'Depth-sum Fractional Thickness-weighted Zonal Acceleration from Diapycnal Mixing',
'm s-2', &
3346 conversion=us%L_T2_to_m_s2)
3347 if (cs%id_hf_dudt_dia_2d > 0)
then
3348 call safe_alloc_ptr(adp%diag_hfrac_u,isdb,iedb,jsd,jed,nz)
3349 call safe_alloc_ptr(adp%du_dt_dia,isdb,iedb,jsd,jed,nz)
3352 cs%id_hf_dvdt_dia_2d = register_diag_field(
'ocean_model',
'hf_dvdt_dia_2d', diag%axesCv1, time, &
3353 'Depth-sum Fractional Thickness-weighted Meridional Acceleration from Diapycnal Mixing',
'm s-2', &
3354 conversion=us%L_T2_to_m_s2)
3355 if (cs%id_hf_dvdt_dia_2d > 0)
then
3356 call safe_alloc_ptr(adp%diag_hfrac_v,isd,ied,jsd,jedb,nz)
3357 call safe_alloc_ptr(adp%dv_dt_dia,isd,ied,jsdb,jedb,nz)
3360 if (cs%use_int_tides)
then
3361 cs%id_cg1 = register_diag_field(
'ocean_model',
'cn1', diag%axesT1, &
3362 time,
'First baroclinic mode (eigen) speed',
'm s-1', conversion=us%L_T_to_m_s)
3363 allocate(cs%id_cn(cs%nMode)) ; cs%id_cn(:) = -1
3365 write(var_name,
'("cn_mode",i1)') m
3366 write(var_descript,
'("Baroclinic (eigen) speed of mode ",i1)') m
3367 cs%id_cn(m) = register_diag_field(
'ocean_model',var_name, diag%axesT1, &
3368 time, var_descript,
'm s-1', conversion=us%L_T_to_m_s)
3369 call mom_mesg(
"Registering "//trim(var_name)//
", Described as: "//var_descript, 5)
3373 if (use_temperature)
then
3374 cs%id_Tdif = register_diag_field(
'ocean_model',
"Tflx_dia_diff", diag%axesTi, &
3375 time,
"Diffusive diapycnal temperature flux across interfaces", &
3376 "degC m s-1", conversion=gv%H_to_m*us%s_to_T)
3377 cs%id_Tadv = register_diag_field(
'ocean_model',
"Tflx_dia_adv", diag%axesTi, &
3378 time,
"Advective diapycnal temperature flux across interfaces", &
3379 "degC m s-1", conversion=gv%H_to_m*us%s_to_T)
3380 cs%id_Sdif = register_diag_field(
'ocean_model',
"Sflx_dia_diff", diag%axesTi, &
3381 time,
"Diffusive diapycnal salnity flux across interfaces", &
3382 "psu m s-1", conversion=gv%H_to_m*us%s_to_T)
3383 cs%id_Sadv = register_diag_field(
'ocean_model',
"Sflx_dia_adv", diag%axesTi, &
3384 time,
"Advective diapycnal salnity flux across interfaces", &
3385 "psu m s-1", conversion=gv%H_to_m*us%s_to_T)
3386 cs%id_MLD_003 = register_diag_field(
'ocean_model',
'MLD_003', diag%axesT1, time, &
3387 'Mixed layer depth (delta rho = 0.03)',
'm', conversion=us%Z_to_m, &
3388 cmor_field_name=
'mlotst', cmor_long_name=
'Ocean Mixed Layer Thickness Defined by Sigma T', &
3389 cmor_standard_name=
'ocean_mixed_layer_thickness_defined_by_sigma_t')
3390 cs%id_mlotstsq = register_diag_field(
'ocean_model',
'mlotstsq', diag%axesT1, time, &
3391 long_name=
'Square of Ocean Mixed Layer Thickness Defined by Sigma T', &
3392 standard_name=
'square_of_ocean_mixed_layer_thickness_defined_by_sigma_t', &
3393 units=
'm2', conversion=us%Z_to_m**2)
3394 cs%id_MLD_0125 = register_diag_field(
'ocean_model',
'MLD_0125', diag%axesT1, time, &
3395 'Mixed layer depth (delta rho = 0.125)',
'm', conversion=us%Z_to_m)
3396 call get_param(param_file, mdl,
"MLD_EN_VALS", cs%MLD_EN_VALS, &
3397 "The energy values used to compute MLDs. If not set (or all set to 0.), the "//&
3398 "default will overwrite to 25., 2500., 250000.",units=
'J/m2', default=0., &
3399 scale=us%kg_m3_to_R*us%m_to_Z**3*us%T_to_s**2)
3400 if ((cs%MLD_EN_VALS(1)==0.).and.(cs%MLD_EN_VALS(2)==0.).and.(cs%MLD_EN_VALS(3)==0.))
then
3401 cs%MLD_EN_VALS = (/25.*us%kg_m3_to_R*us%m_to_Z*us%m_to_L**2*us%T_to_s**2,&
3402 2500.*us%kg_m3_to_R*us%m_to_Z*us%m_to_L**2*us%T_to_s**2,&
3403 250000.*us%kg_m3_to_R*us%m_to_Z*us%m_to_L**2*us%T_to_s**2/)
3405 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
3406 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
3407 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
3408 cs%id_MLD_EN1 = register_diag_field(
'ocean_model',
'MLD_EN1', diag%axesT1, time, &
3409 'Mixed layer depth for energy value set to '//trim(en1)//
' J/m2 (Energy set by 1st MLD_EN_VALS)', &
3410 'm', conversion=us%Z_to_m)
3411 cs%id_MLD_EN2 = register_diag_field(
'ocean_model',
'MLD_EN2', diag%axesT1, time, &
3412 'Mixed layer depth for energy value set to '//trim(en2)//
' J/m2 (Energy set by 2nd MLD_EN_VALS)', &
3413 'm', conversion=us%Z_to_m)
3414 cs%id_MLD_EN3 = register_diag_field(
'ocean_model',
'MLD_EN3', diag%axesT1, time, &
3415 'Mixed layer depth for energy value set to '//trim(en3)//
' J/m2 (Energy set by 3rd MLD_EN_VALS)', &
3416 'm', conversion=us%Z_to_m)
3417 cs%id_subMLN2 = register_diag_field(
'ocean_model',
'subML_N2', diag%axesT1, time, &
3418 'Squared buoyancy frequency below mixed layer',
's-2', conversion=us%s_to_T**2)
3419 cs%id_MLD_user = register_diag_field(
'ocean_model',
'MLD_user', diag%axesT1, time, &
3420 'Mixed layer depth (used defined)',
'm', conversion=us%Z_to_m)
3422 call get_param(param_file, mdl,
"DIAG_MLD_DENSITY_DIFF", cs%MLDdensityDifference, &
3423 "The density difference used to determine a diagnostic mixed "//&
3424 "layer depth, MLD_user, following the definition of Levitus 1982. "//&
3425 "The MLD is the depth at which the density is larger than the "//&
3426 "surface density by the specified amount.", &
3427 units=
'kg/m3', default=0.1, scale=us%kg_m3_to_R)
3428 call get_param(param_file, mdl,
"DIAG_DEPTH_SUBML_N2", cs%dz_subML_N2, &
3429 "The distance over which to calculate a diagnostic of the "//&
3430 "stratification at the base of the mixed layer.", &
3431 units=
'm', default=50.0, scale=us%m_to_Z)
3433 if (cs%id_dudt_dia > 0)
call safe_alloc_ptr(adp%du_dt_dia,isdb,iedb,jsd,jed,nz)
3434 if (cs%id_dvdt_dia > 0)
call safe_alloc_ptr(adp%dv_dt_dia,isd,ied,jsdb,jedb,nz)
3437 cs%id_u_predia = register_diag_field(
'ocean_model',
'u_predia', diag%axesCuL, time, &
3438 'Zonal velocity before diabatic forcing',
'm s-1', conversion=us%L_T_to_m_s)
3439 cs%id_v_predia = register_diag_field(
'ocean_model',
'v_predia', diag%axesCvL, time, &
3440 'Meridional velocity before diabatic forcing',
'm s-1', conversion=us%L_T_to_m_s)
3441 cs%id_h_predia = register_diag_field(
'ocean_model',
'h_predia', diag%axesTL, time, &
3442 'Layer Thickness before diabatic forcing', &
3443 trim(thickness_units), conversion=gv%H_to_MKS, v_extensive=.true.)
3444 cs%id_e_predia = register_diag_field(
'ocean_model',
'e_predia', diag%axesTi, time, &
3445 'Interface Heights before diabatic forcing',
'm')
3446 if (use_temperature)
then
3447 cs%id_T_predia = register_diag_field(
'ocean_model',
'temp_predia', diag%axesTL, time, &
3448 'Potential Temperature',
'degC')
3449 cs%id_S_predia = register_diag_field(
'ocean_model',
'salt_predia', diag%axesTL, time, &
3455 cs%id_Kd_interface = register_diag_field(
'ocean_model',
'Kd_interface', diag%axesTi, time, &
3456 'Total diapycnal diffusivity at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
3457 if (cs%use_energetic_PBL)
then
3458 cs%id_Kd_ePBL = register_diag_field(
'ocean_model',
'Kd_ePBL', diag%axesTi, time, &
3459 'ePBL diapycnal diffusivity at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s)
3462 cs%id_Kd_heat = register_diag_field(
'ocean_model',
'Kd_heat', diag%axesTi, time, &
3463 'Total diapycnal diffusivity for heat at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s, &
3464 cmor_field_name=
'difvho', &
3465 cmor_standard_name=
'ocean_vertical_heat_diffusivity', &
3466 cmor_long_name=
'Ocean vertical heat diffusivity')
3467 cs%id_Kd_salt = register_diag_field(
'ocean_model',
'Kd_salt', diag%axesTi, time, &
3468 'Total diapycnal diffusivity for salt at interfaces',
'm2 s-1', conversion=us%Z2_T_to_m2_s, &
3469 cmor_field_name=
'difvso', &
3470 cmor_standard_name=
'ocean_vertical_salt_diffusivity', &
3471 cmor_long_name=
'Ocean vertical salt diffusivity')
3475 cs%useKPP = kpp_init(param_file, g, gv, us, diag, time, cs%KPP_CSp, passive=cs%KPPisPassive)
3477 allocate( cs%KPP_NLTheat(isd:ied,jsd:jed,nz+1) ) ; cs%KPP_NLTheat(:,:,:) = 0.
3478 allocate( cs%KPP_NLTscalar(isd:ied,jsd:jed,nz+1) ) ; cs%KPP_NLTscalar(:,:,:) = 0.
3481 allocate( cs%KPP_buoy_flux(isd:ied,jsd:jed,nz+1) ) ; cs%KPP_buoy_flux(:,:,:) = 0.
3482 allocate( cs%KPP_temp_flux(isd:ied,jsd:jed) ) ; cs%KPP_temp_flux(:,:) = 0.
3483 allocate( cs%KPP_salt_flux(isd:ied,jsd:jed) ) ; cs%KPP_salt_flux(:,:) = 0.
3490 cs%id_diabatic_diff_h = register_diag_field(
'ocean_model',
'diabatic_diff_h', diag%axesTL, time, &
3491 long_name=
'Cell thickness used during diabatic diffusion', &
3492 units=
'm', conversion=gv%H_to_m, v_extensive=.true.)
3493 if (cs%useALEalgorithm)
then
3494 cs%id_diabatic_diff_temp_tend = register_diag_field(
'ocean_model', &
3495 'diabatic_diff_temp_tendency', diag%axesTL, time, &
3496 'Diabatic diffusion temperature tendency',
'degC s-1', conversion=us%s_to_T)
3497 if (cs%id_diabatic_diff_temp_tend > 0)
then
3498 cs%diabatic_diff_tendency_diag = .true.
3501 cs%id_diabatic_diff_saln_tend = register_diag_field(
'ocean_model',&
3502 'diabatic_diff_saln_tendency', diag%axesTL, time, &
3503 'Diabatic diffusion salinity tendency',
'psu s-1', conversion=us%s_to_T)
3504 if (cs%id_diabatic_diff_saln_tend > 0)
then
3505 cs%diabatic_diff_tendency_diag = .true.
3508 cs%id_diabatic_diff_heat_tend = register_diag_field(
'ocean_model', &
3509 'diabatic_heat_tendency', diag%axesTL, time, &
3510 'Diabatic diffusion heat tendency', &
3511 'W m-2', conversion=us%QRZ_T_to_W_m2, cmor_field_name=
'opottempdiff', &
3512 cmor_standard_name=
'tendency_of_sea_water_potential_temperature_expressed_as_heat_content_'// &
3513 'due_to_parameterized_dianeutral_mixing', &
3514 cmor_long_name=
'Tendency of sea water potential temperature expressed as heat content '// &
3515 'due to parameterized dianeutral mixing', &
3517 if (cs%id_diabatic_diff_heat_tend > 0)
then
3518 cs%diabatic_diff_tendency_diag = .true.
3521 cs%id_diabatic_diff_salt_tend = register_diag_field(
'ocean_model', &
3522 'diabatic_salt_tendency', diag%axesTL, time, &
3523 'Diabatic diffusion of salt tendency', &
3524 'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s, cmor_field_name=
'osaltdiff', &
3525 cmor_standard_name=
'tendency_of_sea_water_salinity_expressed_as_salt_content_'// &
3526 'due_to_parameterized_dianeutral_mixing', &
3527 cmor_long_name=
'Tendency of sea water salinity expressed as salt content '// &
3528 'due to parameterized dianeutral mixing', &
3530 if (cs%id_diabatic_diff_salt_tend > 0)
then
3531 cs%diabatic_diff_tendency_diag = .true.
3535 cs%id_diabatic_diff_heat_tend_2d = register_diag_field(
'ocean_model', &
3536 'diabatic_heat_tendency_2d', diag%axesT1, time, &
3537 'Depth integrated diabatic diffusion heat tendency', &
3538 'W m-2', conversion=us%QRZ_T_to_W_m2, cmor_field_name=
'opottempdiff_2d', &
3539 cmor_standard_name=
'tendency_of_sea_water_potential_temperature_expressed_as_heat_content_'//&
3540 'due_to_parameterized_dianeutral_mixing_depth_integrated', &
3541 cmor_long_name=
'Tendency of sea water potential temperature expressed as heat content '//&
3542 'due to parameterized dianeutral mixing depth integrated')
3543 if (cs%id_diabatic_diff_heat_tend_2d > 0)
then
3544 cs%diabatic_diff_tendency_diag = .true.
3548 cs%id_diabatic_diff_salt_tend_2d = register_diag_field(
'ocean_model', &
3549 'diabatic_salt_tendency_2d', diag%axesT1, time, &
3550 'Depth integrated diabatic diffusion salt tendency', &
3551 'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s, cmor_field_name=
'osaltdiff_2d', &
3552 cmor_standard_name=
'tendency_of_sea_water_salinity_expressed_as_salt_content_'// &
3553 'due_to_parameterized_dianeutral_mixing_depth_integrated', &
3554 cmor_long_name=
'Tendency of sea water salinity expressed as salt content '// &
3555 'due to parameterized dianeutral mixing depth integrated')
3556 if (cs%id_diabatic_diff_salt_tend_2d > 0)
then
3557 cs%diabatic_diff_tendency_diag = .true.
3563 cs%id_boundary_forcing_h = register_diag_field(
'ocean_model',
'boundary_forcing_h', diag%axesTL, time, &
3564 long_name=
'Cell thickness after applying boundary forcing', &
3565 units=
'm', conversion=gv%H_to_m, v_extensive=.true.)
3566 cs%id_boundary_forcing_h_tendency = register_diag_field(
'ocean_model', &
3567 'boundary_forcing_h_tendency', diag%axesTL, time, &
3568 'Cell thickness tendency due to boundary forcing',
'm s-1', &
3569 conversion=gv%H_to_m*us%s_to_T, v_extensive=.true.)
3570 if (cs%id_boundary_forcing_h_tendency > 0)
then
3571 cs%boundary_forcing_tendency_diag = .true.
3574 cs%id_boundary_forcing_temp_tend = register_diag_field(
'ocean_model',&
3575 'boundary_forcing_temp_tendency', diag%axesTL, time, &
3576 'Boundary forcing temperature tendency',
'degC s-1', conversion=us%s_to_T)
3577 if (cs%id_boundary_forcing_temp_tend > 0)
then
3578 cs%boundary_forcing_tendency_diag = .true.
3581 cs%id_boundary_forcing_saln_tend = register_diag_field(
'ocean_model',&
3582 'boundary_forcing_saln_tendency', diag%axesTL, time, &
3583 'Boundary forcing saln tendency',
'psu s-1', conversion=us%s_to_T)
3584 if (cs%id_boundary_forcing_saln_tend > 0)
then
3585 cs%boundary_forcing_tendency_diag = .true.
3588 cs%id_boundary_forcing_heat_tend = register_diag_field(
'ocean_model',&
3589 'boundary_forcing_heat_tendency', diag%axesTL, time, &
3590 'Boundary forcing heat tendency', &
3591 'W m-2', conversion=us%QRZ_T_to_W_m2, v_extensive = .true.)
3592 if (cs%id_boundary_forcing_heat_tend > 0)
then
3593 cs%boundary_forcing_tendency_diag = .true.
3596 cs%id_boundary_forcing_salt_tend = register_diag_field(
'ocean_model',&
3597 'boundary_forcing_salt_tendency', diag%axesTL, time, &
3598 'Boundary forcing salt tendency',
'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s, &
3599 v_extensive = .true.)
3600 if (cs%id_boundary_forcing_salt_tend > 0)
then
3601 cs%boundary_forcing_tendency_diag = .true.
3605 cs%id_boundary_forcing_heat_tend_2d = register_diag_field(
'ocean_model',&
3606 'boundary_forcing_heat_tendency_2d', diag%axesT1, time, &
3607 'Depth integrated boundary forcing of ocean heat', &
3608 'W m-2', conversion=us%QRZ_T_to_W_m2)
3609 if (cs%id_boundary_forcing_heat_tend_2d > 0)
then
3610 cs%boundary_forcing_tendency_diag = .true.
3614 cs%id_boundary_forcing_salt_tend_2d = register_diag_field(
'ocean_model',&
3615 'boundary_forcing_salt_tendency_2d', diag%axesT1, time, &
3616 'Depth integrated boundary forcing of ocean salt', &
3617 'kg m-2 s-1', conversion=us%RZ_T_to_kg_m2s)
3618 if (cs%id_boundary_forcing_salt_tend_2d > 0)
then
3619 cs%boundary_forcing_tendency_diag = .true.
3624 cs%id_frazil_h = register_diag_field(
'ocean_model',
'frazil_h', diag%axesTL, time, &
3625 long_name=
'Cell Thickness', standard_name=
'cell_thickness', &
3626 units=
'm', conversion=gv%H_to_m, v_extensive=.true.)
3629 cs%id_frazil_temp_tend = register_diag_field(
'ocean_model',&
3630 'frazil_temp_tendency', diag%axesTL, time, &
3631 'Temperature tendency due to frazil formation',
'degC s-1', conversion=us%s_to_T)
3632 if (cs%id_frazil_temp_tend > 0)
then
3633 cs%frazil_tendency_diag = .true.
3637 cs%id_frazil_heat_tend = register_diag_field(
'ocean_model',&
3638 'frazil_heat_tendency', diag%axesTL, time, &
3639 'Heat tendency due to frazil formation', &
3640 'W m-2', conversion=us%QRZ_T_to_W_m2, v_extensive=.true.)
3641 if (cs%id_frazil_heat_tend > 0)
then
3642 cs%frazil_tendency_diag = .true.
3646 cs%id_frazil_heat_tend_2d = register_diag_field(
'ocean_model',&
3647 'frazil_heat_tendency_2d', diag%axesT1, time, &
3648 'Depth integrated heat tendency due to frazil formation', &
3649 'W m-2', conversion=us%QRZ_T_to_W_m2)
3650 if (cs%id_frazil_heat_tend_2d > 0)
then
3651 cs%frazil_tendency_diag = .true.
3654 if (cs%frazil_tendency_diag)
then
3655 allocate(cs%frazil_temp_diag(isd:ied,jsd:jed,nz) ) ; cs%frazil_temp_diag(:,:,:) = 0.
3656 allocate(cs%frazil_heat_diag(isd:ied,jsd:jed,nz) ) ; cs%frazil_heat_diag(:,:,:) = 0.
3660 cs%use_CVMix_conv = cvmix_conv_init(time, g, gv, us, param_file, diag, cs%CVMix_conv_csp)
3662 call entrain_diffusive_init(time, g, gv, us, param_file, diag, cs%entrain_diffusive_CSp, &
3663 just_read_params=cs%useALEalgorithm)
3666 if (cs%use_geothermal) &
3667 call geothermal_init(time, g, gv, us, param_file, diag, cs%geothermal_CSp, usealealgorithm)
3670 if (cs%use_int_tides)
then
3671 call int_tide_input_init(time, g, gv, us, param_file, diag, cs%int_tide_input_CSp, &
3673 call internal_tides_init(time, g, gv, us, param_file, diag, cs%int_tide_CSp)
3677 call set_diffusivity_init(time, g, gv, us, param_file, diag, cs%set_diff_CSp, cs%int_tide_CSp, &
3678 halo_ts=cs%halo_TS_diff, double_diffuse=cs%double_diffuse)
3680 if (cs%useKPP .and. (cs%double_diffuse .and. .not.cs%use_CVMix_ddiff)) &
3681 call mom_error(fatal,
'diabatic_driver_init: DOUBLE_DIFFUSION (old method) does not work '//&
3682 'with KPP. Please set DOUBLE_DIFFUSION=False and USE_CVMIX_DDIFF=True.')
3685 id_clock_entrain = cpu_clock_id(
'(Ocean diabatic entrain)', grain=clock_module)
3686 if (cs%bulkmixedlayer) &
3687 id_clock_mixedlayer = cpu_clock_id(
'(Ocean mixed layer)', grain=clock_module)
3688 id_clock_remap = cpu_clock_id(
'(Ocean vert remap)', grain=clock_module)
3689 if (cs%use_geothermal) &
3690 id_clock_geothermal = cpu_clock_id(
'(Ocean geothermal)', grain=clock_routine)
3691 id_clock_set_diffusivity = cpu_clock_id(
'(Ocean set_diffusivity)', grain=clock_module)
3692 id_clock_kpp = cpu_clock_id(
'(Ocean KPP)', grain=clock_module)
3693 id_clock_tracers = cpu_clock_id(
'(Ocean tracer_columns)', grain=clock_module_driver+5)
3694 if (cs%use_sponge) &
3695 id_clock_sponge = cpu_clock_id(
'(Ocean sponges)', grain=clock_module)
3696 id_clock_tridiag = cpu_clock_id(
'(Ocean diabatic tridiag)', grain=clock_routine)
3697 id_clock_pass = cpu_clock_id(
'(Ocean diabatic message passing)', grain=clock_routine)
3698 id_clock_differential_diff = -1 ;
if (cs%double_diffuse .and. .not.cs%use_CVMix_ddiff) &
3699 id_clock_differential_diff = cpu_clock_id(
'(Ocean differential diffusion)', grain=clock_routine)
3702 call diabatic_aux_init(time, g, gv, us, param_file, diag, cs%diabatic_aux_CSp, &
3703 cs%useALEalgorithm, cs%use_energetic_PBL)
3706 if (cs%bulkmixedlayer) &
3707 call bulkmixedlayer_init(time, g, gv, us, param_file, diag, cs%bulkmixedlayer_CSp)
3708 if (cs%use_energetic_PBL) &
3709 call energetic_pbl_init(time, g, gv, us, param_file, diag, cs%energetic_PBL_CSp)
3711 call regularize_layers_init(time, g, gv, param_file, diag, cs%regularize_layers_CSp)
3713 if (cs%debug_energy_req) &
3714 call diapyc_energy_req_init(time, g, gv, us, param_file, diag, cs%diapyc_en_rec_CSp)
3717 if (use_temperature)
then
3718 call get_param(param_file, mdl,
"PEN_SW_NBANDS", nbands, default=1)
3719 if (nbands > 0)
then
3721 call opacity_init(time, g, gv, us, param_file, diag, cs%opacity_CSp, cs%optics)
3726 call diag_grid_storage_init(cs%diag_grids_prev, g, diag)
3728 end subroutine diabatic_driver_init
3732 subroutine diabatic_driver_end(CS)
3735 if (.not.
associated(cs))
return
3737 call diabatic_aux_end(cs%diabatic_aux_CSp)
3739 call entrain_diffusive_end(cs%entrain_diffusive_CSp)
3740 call set_diffusivity_end(cs%set_diff_CSp)
3743 deallocate( cs%KPP_buoy_flux )
3744 deallocate( cs%KPP_temp_flux )
3745 deallocate( cs%KPP_salt_flux )
3746 deallocate( cs%KPP_NLTheat )
3747 deallocate( cs%KPP_NLTscalar )
3748 call kpp_end(cs%KPP_CSp)
3751 if (cs%use_CVMix_conv)
call cvmix_conv_end(cs%CVMix_conv_csp)
3753 if (cs%use_energetic_PBL) &
3754 call energetic_pbl_end(cs%energetic_PBL_CSp)
3755 if (cs%debug_energy_req) &
3756 call diapyc_energy_req_end(cs%diapyc_en_rec_CSp)
3758 if (
associated(cs%optics))
then
3759 call opacity_end(cs%opacity_CSp, cs%optics)
3760 deallocate(cs%optics)
3771 end subroutine diabatic_driver_end