7 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
8 use mom_cpu_clock, only : clock_module_driver, clock_module, clock_routine
25 use time_interp_external_mod
, only : init_external_field, time_interp_external
26 use time_interp_external_mod
, only : time_interp_external_init
28 implicit none ;
private 30 #include <MOM_memory.h> 32 public diabatic_aux_init, diabatic_aux_end
33 public make_frazil, adjust_salt, differential_diffuse_t_s, tridiagts
34 public find_uv_at_h, diagnosemldbydensitydifference, applyboundaryfluxesinout, set_pen_shortwave
35 public diagnosemldbyenergy
44 logical :: do_rivermix = .false.
46 real :: rivermix_depth = 0.0
47 logical :: reclaim_frazil
50 logical :: pressure_dependent_frazil
54 logical :: ignore_fluxes_over_land
58 logical :: use_river_heat_content
61 logical :: use_calving_heat_content
68 logical :: chl_from_file
70 type(time_type),
pointer :: time => null()
74 integer :: id_createdh = -1
75 integer :: id_brine_lay = -1
76 integer :: id_pensw_diag = -1
77 integer :: id_penswflux_diag = -1
78 integer :: id_nonpensw_diag = -1
79 integer :: id_chl = -1
82 real,
allocatable,
dimension(:,:) :: createdh
84 real,
allocatable,
dimension(:,:,:) :: pensw_diag
86 real,
allocatable,
dimension(:,:,:) :: penswflux_diag
88 real,
allocatable,
dimension(:,:) :: nonpensw_diag
94 integer :: id_clock_uv_at_h, id_clock_frazil
103 subroutine make_frazil(h, tv, G, GV, US, CS, p_surf, halo)
106 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
113 real,
dimension(SZI_(G),SZJ_(G)), &
114 optional,
intent(in) :: p_surf
115 integer,
optional,
intent(in) :: halo
118 real,
dimension(SZI_(G)) :: &
122 real,
dimension(SZI_(G),SZK_(G)) :: &
128 integer :: i, j, k, is, ie, js, je, nz
130 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
131 if (
present(halo))
then 132 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
135 call cpu_clock_begin(id_clock_frazil)
137 if (.not.cs%pressure_dependent_frazil)
then 138 do k=1,nz ;
do i=is,ie ; pressure(i,k) = 0.0 ;
enddo ;
enddo 140 h_to_rl2_t2 = gv%H_to_RZ * gv%g_Earth
146 if (
PRESENT(p_surf))
then ;
do i=is,ie
150 do i=is,ie ; fraz_col(i) = 0.0 ;
enddo 152 if (cs%pressure_dependent_frazil)
then 154 pressure(i,1) = ps(i) + (0.5*h_to_rl2_t2)*h(i,j,1)
156 do k=2,nz ;
do i=is,ie
157 pressure(i,k) = pressure(i,k-1) + &
158 (0.5*h_to_rl2_t2) * (h(i,j,k) + h(i,j,k-1))
162 if (cs%reclaim_frazil)
then 164 do i=is,ie ;
if (tv%frazil(i,j) > 0.0)
then 165 if (.not.t_fr_set)
then 167 1, ie-i+1, tv%eqn_of_state, pres_scale=us%RL2_T2_to_Pa)
171 if (tv%T(i,j,1) > t_freeze(i))
then 174 hc = (tv%C_p*gv%H_to_RZ) * h(i,j,1)
175 if (tv%frazil(i,j) - hc * (tv%T(i,j,1) - t_freeze(i)) <= 0.0)
then 176 tv%T(i,j,1) = tv%T(i,j,1) - tv%frazil(i,j) / hc
179 tv%frazil(i,j) = tv%frazil(i,j) - hc * (tv%T(i,j,1) - t_freeze(i))
180 tv%T(i,j,1) = t_freeze(i)
189 if ((g%mask2dT(i,j) > 0.0) .and. &
190 ((tv%T(i,j,k) < 0.0) .or. (fraz_col(i) > 0.0)))
then 191 if (.not.t_fr_set)
then 193 1, ie-i+1, tv%eqn_of_state, pres_scale=us%RL2_T2_to_Pa)
197 hc = (tv%C_p*gv%H_to_RZ) * h(i,j,k)
198 if (h(i,j,k) <= 10.0*(gv%Angstrom_H + gv%H_subroundoff))
then 200 if (tv%T(i,j,k) < t_freeze(i))
then 201 fraz_col(i) = fraz_col(i) + hc * (t_freeze(i) - tv%T(i,j,k))
202 tv%T(i,j,k) = t_freeze(i)
204 elseif ((fraz_col(i) > 0.0) .or. (tv%T(i,j,k) < t_freeze(i)))
then 205 if (fraz_col(i) + hc * (t_freeze(i) - tv%T(i,j,k)) < 0.0)
then 206 tv%T(i,j,k) = tv%T(i,j,k) - fraz_col(i) / hc
209 fraz_col(i) = fraz_col(i) + hc * (t_freeze(i) - tv%T(i,j,k))
210 tv%T(i,j,k) = t_freeze(i)
217 tv%frazil(i,j) = tv%frazil(i,j) + fraz_col(i)
220 call cpu_clock_end(id_clock_frazil)
222 end subroutine make_frazil
226 subroutine differential_diffuse_t_s(h, tv, visc, dt, G, GV)
229 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
235 real,
intent(in) :: dt
238 real,
dimension(SZI_(G)) :: &
241 real,
dimension(SZI_(G),SZK_(G)) :: &
243 real,
dimension(SZI_(G),SZK_(G)+1) :: &
253 real,
dimension(:,:,:),
pointer :: T=>null(), s=>null()
254 real,
dimension(:,:,:),
pointer :: Kd_T=>null(), kd_s=>null()
255 integer :: i, j, k, is, ie, js, je, nz
257 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
258 h_neglect = gv%H_subroundoff
260 if (.not.
associated(tv%T))
call mom_error(fatal, &
261 "differential_diffuse_T_S: Called with an unassociated tv%T")
262 if (.not.
associated(tv%S))
call mom_error(fatal, &
263 "differential_diffuse_T_S: Called with an unassociated tv%S")
264 if (.not.
associated(visc%Kd_extra_T))
call mom_error(fatal, &
265 "differential_diffuse_T_S: Called with an unassociated visc%Kd_extra_T")
266 if (.not.
associated(visc%Kd_extra_S))
call mom_error(fatal, &
267 "differential_diffuse_T_S: Called with an unassociated visc%Kd_extra_S")
269 t => tv%T ; s => tv%S
270 kd_t => visc%Kd_extra_T ; kd_s => visc%Kd_extra_S
276 i_h_int = 1.0 / (0.5 * (h(i,j,1) + h(i,j,2)) + h_neglect)
277 mix_t(i,2) = ((dt * kd_t(i,j,2)) * gv%Z_to_H**2) * i_h_int
278 mix_s(i,2) = ((dt * kd_s(i,j,2)) * gv%Z_to_H**2) * i_h_int
280 h_tr = h(i,j,1) + h_neglect
281 b1_t(i) = 1.0 / (h_tr + mix_t(i,2))
282 b1_s(i) = 1.0 / (h_tr + mix_s(i,2))
283 d1_t(i) = h_tr * b1_t(i)
284 d1_s(i) = h_tr * b1_s(i)
285 t(i,j,1) = (b1_t(i)*h_tr)*t(i,j,1)
286 s(i,j,1) = (b1_s(i)*h_tr)*s(i,j,1)
288 do k=2,nz-1 ;
do i=is,ie
290 i_h_int = 1.0 / (0.5 * (h(i,j,k) + h(i,j,k+1)) + h_neglect)
291 mix_t(i,k+1) = ((dt * kd_t(i,j,k+1)) * gv%Z_to_H**2) * i_h_int
292 mix_s(i,k+1) = ((dt * kd_s(i,j,k+1)) * gv%Z_to_H**2) * i_h_int
294 c1_t(i,k) = mix_t(i,k) * b1_t(i)
295 c1_s(i,k) = mix_s(i,k) * b1_s(i)
297 h_tr = h(i,j,k) + h_neglect
298 b_denom_t = h_tr + d1_t(i)*mix_t(i,k)
299 b_denom_s = h_tr + d1_s(i)*mix_s(i,k)
300 b1_t(i) = 1.0 / (b_denom_t + mix_t(i,k+1))
301 b1_s(i) = 1.0 / (b_denom_s + mix_s(i,k+1))
302 d1_t(i) = b_denom_t * b1_t(i)
303 d1_s(i) = b_denom_s * b1_s(i)
305 t(i,j,k) = b1_t(i) * (h_tr*t(i,j,k) + mix_t(i,k)*t(i,j,k-1))
306 s(i,j,k) = b1_s(i) * (h_tr*s(i,j,k) + mix_s(i,k)*s(i,j,k-1))
309 c1_t(i,nz) = mix_t(i,nz) * b1_t(i)
310 c1_s(i,nz) = mix_s(i,nz) * b1_s(i)
312 h_tr = h(i,j,nz) + h_neglect
313 b1_t(i) = 1.0 / (h_tr + d1_t(i)*mix_t(i,nz))
314 b1_s(i) = 1.0 / (h_tr + d1_s(i)*mix_s(i,nz))
316 t(i,j,nz) = b1_t(i) * (h_tr*t(i,j,nz) + mix_t(i,nz)*t(i,j,nz-1))
317 s(i,j,nz) = b1_s(i) * (h_tr*s(i,j,nz) + mix_s(i,nz)*s(i,j,nz-1))
319 do k=nz-1,1,-1 ;
do i=is,ie
320 t(i,j,k) = t(i,j,k) + c1_t(i,k+1)*t(i,j,k+1)
321 s(i,j,k) = s(i,j,k) + c1_s(i,k+1)*s(i,j,k+1)
324 end subroutine differential_diffuse_t_s
329 subroutine adjust_salt(h, tv, G, GV, CS, halo)
332 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
338 integer,
optional,
intent(in) :: halo
341 real :: salt_add_col(szi_(g),szj_(g))
344 integer :: i, j, k, is, ie, js, je, nz
346 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
347 if (
present(halo))
then 348 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
353 s_min = tv%min_salinity
355 salt_add_col(:,:) = 0.0
359 do k=nz,1,-1 ;
do i=is,ie
360 if ( (g%mask2dT(i,j) > 0.0) .and. &
361 ((tv%S(i,j,k) < s_min) .or. (salt_add_col(i,j) > 0.0)) )
then 362 mc = gv%H_to_RZ * h(i,j,k)
363 if (h(i,j,k) <= 10.0*gv%Angstrom_H)
then 365 if (tv%S(i,j,k) < s_min)
then 366 salt_add_col(i,j) = salt_add_col(i,j) + mc * (s_min - tv%S(i,j,k))
369 elseif (salt_add_col(i,j) + mc * (s_min - tv%S(i,j,k)) <= 0.0)
then 370 tv%S(i,j,k) = tv%S(i,j,k) - salt_add_col(i,j) / mc
371 salt_add_col(i,j) = 0.0
373 salt_add_col(i,j) = salt_add_col(i,j) + mc * (s_min - tv%S(i,j,k))
379 tv%salt_deficit(i,j) = tv%salt_deficit(i,j) + salt_add_col(i,j)
384 end subroutine adjust_salt
388 subroutine tridiagts(G, GV, is, ie, js, je, hold, ea, eb, T, S)
391 integer,
intent(in) :: is
392 integer,
intent(in) :: ie
393 integer,
intent(in) :: js
394 integer,
intent(in) :: je
395 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: hold
397 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: ea
399 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: eb
401 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: T
402 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: S
405 real :: b1(szib_(g)), d1(szib_(g))
406 real :: c1(szib_(g),szk_(g))
407 real :: h_tr, b_denom_1
413 h_tr = hold(i,j,1) + gv%H_subroundoff
414 b1(i) = 1.0 / (h_tr + eb(i,j,1))
416 t(i,j,1) = (b1(i)*h_tr)*t(i,j,1)
417 s(i,j,1) = (b1(i)*h_tr)*s(i,j,1)
419 do k=2,g%ke ;
do i=is,ie
420 c1(i,k) = eb(i,j,k-1) * b1(i)
421 h_tr = hold(i,j,k) + gv%H_subroundoff
422 b_denom_1 = h_tr + d1(i)*ea(i,j,k)
423 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
424 d1(i) = b_denom_1 * b1(i)
425 t(i,j,k) = b1(i) * (h_tr*t(i,j,k) + ea(i,j,k)*t(i,j,k-1))
426 s(i,j,k) = b1(i) * (h_tr*s(i,j,k) + ea(i,j,k)*s(i,j,k-1))
428 do k=g%ke-1,1,-1 ;
do i=is,ie
429 t(i,j,k) = t(i,j,k) + c1(i,k+1)*t(i,j,k+1)
430 s(i,j,k) = s(i,j,k) + c1(i,k+1)*s(i,j,k+1)
433 end subroutine tridiagts
437 subroutine find_uv_at_h(u, v, h, u_h, v_h, G, GV, US, ea, eb)
441 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
443 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
445 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
447 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
449 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
451 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
452 optional,
intent(in) :: ea
455 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
456 optional,
intent(in) :: eb
464 real :: b1(szi_(g)), d1(szi_(g)), c1(szi_(g),szk_(g))
465 real :: a_n(szi_(g)), a_s(szi_(g))
466 real :: a_e(szi_(g)), a_w(szi_(g))
468 real :: sum_area, Idenom
469 logical :: mix_vertically
470 integer :: i, j, k, is, ie, js, je, nz
471 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
472 call cpu_clock_begin(id_clock_uv_at_h)
473 h_neglect = gv%H_subroundoff
475 mix_vertically =
present(ea)
476 if (
present(ea) .neqv.
present(eb))
call mom_error(fatal, &
477 "find_uv_at_h: Either both ea and eb or neither one must be present "// &
478 "in call to find_uv_at_h.")
484 sum_area = g%areaCu(i-1,j) + g%areaCu(i,j)
485 if (sum_area>0.0)
then 486 idenom = sqrt(0.5*g%IareaT(i,j) / sum_area)
487 a_w(i) = g%areaCu(i-1,j) * idenom
488 a_e(i) = g%areaCu(i,j) * idenom
490 a_w(i) = 0.0 ; a_e(i) = 0.0
493 sum_area = g%areaCv(i,j-1) + g%areaCv(i,j)
494 if (sum_area>0.0)
then 495 idenom = sqrt(0.5*g%IareaT(i,j) / sum_area)
496 a_s(i) = g%areaCv(i,j-1) * idenom
497 a_n(i) = g%areaCv(i,j) * idenom
499 a_s(i) = 0.0 ; a_n(i) = 0.0
503 if (mix_vertically)
then 505 b_denom_1 = h(i,j,1) + h_neglect
506 b1(i) = 1.0 / (b_denom_1 + eb(i,j,1))
507 d1(i) = b_denom_1 * b1(i)
508 u_h(i,j,1) = (h(i,j,1)*b1(i)) * (a_e(i)*u(i,j,1) + a_w(i)*u(i-1,j,1))
509 v_h(i,j,1) = (h(i,j,1)*b1(i)) * (a_n(i)*v(i,j,1) + a_s(i)*v(i,j-1,1))
511 do k=2,nz ;
do i=is,ie
512 c1(i,k) = eb(i,j,k-1) * b1(i)
513 b_denom_1 = h(i,j,k) + d1(i)*ea(i,j,k) + h_neglect
514 b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
515 d1(i) = b_denom_1 * b1(i)
516 u_h(i,j,k) = (h(i,j,k) * (a_e(i)*u(i,j,k) + a_w(i)*u(i-1,j,k)) + &
517 ea(i,j,k)*u_h(i,j,k-1))*b1(i)
518 v_h(i,j,k) = (h(i,j,k) * (a_n(i)*v(i,j,k) + a_s(i)*v(i,j-1,k)) + &
519 ea(i,j,k)*v_h(i,j,k-1))*b1(i)
521 do k=nz-1,1,-1 ;
do i=is,ie
522 u_h(i,j,k) = u_h(i,j,k) + c1(i,k+1)*u_h(i,j,k+1)
523 v_h(i,j,k) = v_h(i,j,k) + c1(i,k+1)*v_h(i,j,k+1)
526 do k=1,nz ;
do i=is,ie
527 u_h(i,j,k) = a_e(i)*u(i,j,k) + a_w(i)*u(i-1,j,k)
528 v_h(i,j,k) = a_n(i)*v(i,j,k) + a_s(i)*v(i,j-1,k)
533 call cpu_clock_end(id_clock_uv_at_h)
534 end subroutine find_uv_at_h
537 subroutine set_pen_shortwave(optics, fluxes, G, GV, US, CS, opacity_CSp, tracer_flow_CSp)
540 type(
forcing),
intent(inout) :: fluxes
551 real,
dimension(SZI_(G),SZJ_(G)) :: chl_2d
552 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: chl_3d
553 character(len=128) :: mesg
554 integer :: i, j, k, is, ie, js, je
555 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
557 if (.not.
associated(optics))
return 559 if (cs%var_pen_sw)
then 560 if (cs%chl_from_file)
then 563 call time_interp_external(cs%sbc_chl, cs%Time, chl_2d)
564 do j=js,je ;
do i=is,ie
565 if ((g%mask2dT(i,j) > 0.5) .and. (chl_2d(i,j) < 0.0))
then 566 write(mesg,
'(" Time_interp negative chl of ",(1pe12.4)," at i,j = ",& 567 & 2(i3), "lon/lat = ",(1pe12.4)," E ", (1pe12.4), " N.")') &
568 chl_2d(i,j), i, j, g%geoLonT(i,j), g%geoLatT(i,j)
569 call mom_error(fatal,
"MOM_diabatic_aux set_pen_shortwave: "//trim(mesg))
573 if (cs%id_chl > 0)
call post_data(cs%id_chl, chl_2d, cs%diag)
575 call set_opacity(optics, fluxes%sw, fluxes%sw_vis_dir, fluxes%sw_vis_dif, &
576 fluxes%sw_nir_dir, fluxes%sw_nir_dif, g, gv, us, opacity_csp, chl_2d=chl_2d)
578 if (.not.
associated(tracer_flow_csp))
call mom_error(fatal, &
579 "The tracer flow control structure must be associated when the model sets "//&
580 "the chlorophyll internally in set_pen_shortwave.")
581 call get_chl_from_model(chl_3d, g, tracer_flow_csp)
583 if (cs%id_chl > 0)
call post_data(cs%id_chl, chl_3d(:,:,1), cs%diag)
585 call set_opacity(optics, fluxes%sw, fluxes%sw_vis_dir, fluxes%sw_vis_dif, &
586 fluxes%sw_nir_dir, fluxes%sw_nir_dif, g, gv, us, opacity_csp, chl_3d=chl_3d)
589 call set_opacity(optics, fluxes%sw, fluxes%sw_vis_dir, fluxes%sw_vis_dif, &
590 fluxes%sw_nir_dir, fluxes%sw_nir_dif, g, gv, us, opacity_csp)
593 end subroutine set_pen_shortwave
598 subroutine diagnosemldbydensitydifference(id_MLD, h, tv, densityDiff, G, GV, US, diagPtr, &
599 id_N2subML, id_MLDsq, dz_subML)
603 integer,
intent(in) :: id_MLD
604 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
608 real,
intent(in) :: densityDiff
610 integer,
optional,
intent(in) :: id_N2subML
611 integer,
optional,
intent(in) :: id_MLDsq
612 real,
optional,
intent(in) :: dz_subML
616 real,
dimension(SZI_(G)) :: deltaRhoAtKm1, deltaRhoAtK
617 real,
dimension(SZI_(G)) :: pRef_MLD, pRef_N2
618 real,
dimension(SZI_(G)) :: H_subML, dH_N2
619 real,
dimension(SZI_(G)) :: T_subML, T_deeper
620 real,
dimension(SZI_(G)) :: S_subML, S_deeper
621 real,
dimension(SZI_(G)) :: rho_subML, rho_deeper
622 real,
dimension(SZI_(G)) :: dK, dKm1
623 real,
dimension(SZI_(G)) :: rhoSurf
624 real,
dimension(SZI_(G), SZJ_(G)) :: MLD
625 real,
dimension(SZI_(G), SZJ_(G)) :: subMLN2
626 real,
dimension(SZI_(G), SZJ_(G)) :: MLD2
627 logical,
dimension(SZI_(G)) :: N2_region_set
633 integer,
dimension(2) :: EOSdom
634 integer :: i, j, is, ie, js, je, k, nz, id_N2, id_SQ
636 id_n2 = -1 ;
if (
PRESENT(id_n2subml)) id_n2 = id_n2subml
638 id_sq = -1 ;
if (
PRESENT(id_mldsq)) id_sq = id_mldsq
640 ge_rho0 = us%L_to_Z**2*gv%g_Earth / (gv%Rho0)
641 dh_subml = 50.*gv%m_to_H ;
if (
present(dz_subml)) dh_subml = gv%Z_to_H*dz_subml
643 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
646 eosdom(:) = eos_domain(g%HI)
648 do i=is,ie ; dk(i) = 0.5 * h(i,j,1) * gv%H_to_Z ;
enddo 649 call calculate_density(tv%T(:,j,1), tv%S(:,j,1), pref_mld, rhosurf, tv%eqn_of_state, eosdom)
655 h_subml(i) = h(i,j,1) ; dh_n2(i) = 0.0
656 t_subml(i) = 0.0 ; s_subml(i) = 0.0 ; t_deeper(i) = 0.0 ; s_deeper(i) = 0.0
657 n2_region_set(i) = (g%mask2dT(i,j)<0.5)
663 dk(i) = dk(i) + 0.5 * ( h(i,j,k) + h(i,j,k-1) ) * gv%H_to_Z
670 if (mld(i,j)==0.0)
then 671 h_subml(i) = h_subml(i) + h(i,j,k)
672 elseif (.not.n2_region_set(i))
then 673 if (dh_n2(i)==0.0)
then 674 t_subml(i) = tv%T(i,j,k) ; s_subml(i) = tv%S(i,j,k)
675 h_subml(i) = h_subml(i) + 0.5 * h(i,j,k)
676 dh_n2(i) = 0.5 * h(i,j,k)
677 elseif (dh_n2(i) + h(i,j,k) < dh_subml)
then 678 dh_n2(i) = dh_n2(i) + h(i,j,k)
680 t_deeper(i) = tv%T(i,j,k) ; s_deeper(i) = tv%S(i,j,k)
681 dh_n2(i) = dh_n2(i) + 0.5 * h(i,j,k)
682 n2_region_set(i) = .true.
689 do i=is,ie ; deltarhoatkm1(i) = deltarhoatk(i) ;
enddo 690 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pref_mld, deltarhoatk, tv%eqn_of_state, eosdom)
692 deltarhoatk(i) = deltarhoatk(i) - rhosurf(i)
693 ddrho = deltarhoatk(i) - deltarhoatkm1(i)
694 if ((mld(i,j)==0.) .and. (ddrho>0.) .and. &
695 (deltarhoatkm1(i)<densitydiff) .and. (deltarhoatk(i)>=densitydiff))
then 696 afac = ( densitydiff - deltarhoatkm1(i) ) / ddrho
697 mld(i,j) = dk(i) * afac + dkm1(i) * (1. - afac)
699 if (id_sq > 0) mld2(i,j) = mld(i,j)**2
703 if ((mld(i,j)==0.) .and. (deltarhoatk(i)<densitydiff)) mld(i,j) = dk(i)
707 do i=is,ie ; pref_n2(i) = (gv%g_Earth * gv%H_to_RZ) * (h_subml(i) + 0.5*dh_n2(i)) ;
enddo 713 call calculate_density(t_subml, s_subml, pref_n2, rho_subml, tv%eqn_of_state, eosdom)
714 call calculate_density(t_deeper, s_deeper, pref_n2, rho_deeper, tv%eqn_of_state, eosdom)
715 do i=is,ie ;
if ((g%mask2dT(i,j)>0.5) .and. n2_region_set(i))
then 716 submln2(i,j) = ge_rho0 * (rho_deeper(i) - rho_subml(i)) / (gv%H_to_z * dh_n2(i))
721 if (id_mld > 0)
call post_data(id_mld, mld, diagptr)
722 if (id_n2 > 0)
call post_data(id_n2, submln2, diagptr)
723 if (id_sq > 0)
call post_data(id_sq, mld2, diagptr)
725 end subroutine diagnosemldbydensitydifference
729 subroutine diagnosemldbyenergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, diagPtr)
751 integer,
dimension(3),
intent(in) :: id_MLD
755 real,
dimension(3),
intent(in) :: Mixing_Energy
756 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
763 real,
dimension(SZI_(G), SZJ_(G),3) :: MLD
764 real,
dimension(SZK_(G)) :: Z_L, Z_U, dZ, Rho_c, pRef_MLD
765 real,
dimension(3) :: PE_threshold
768 real :: PE_Threshold_fraction, PE, PE_Mixed, PE_Mixed_TST
769 real :: RhoDZ_ML, H_ML, RhoDZ_ML_TST, H_ML_TST
772 real :: R1, D1, R2, D2
773 real :: Ca, Cb,D ,Cc, Cd, Ca2, Cb2, C, Cc2
774 real :: Gx, Gpx, Hx, iHx, Hpx, Ix, Ipx, Fgx, Fpx, X, X2
777 integer :: i, j, is, ie, js, je, k, nz
779 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
783 pe_threshold_fraction = 1.e-4
786 pe_threshold(im) = mixing_energy(im)/gv%g_earth
789 do j=js,je;
do i=is,ie
790 if (g%mask2dT(i,j) > 0.0)
then 792 call calculate_density(tv%T(i,j,:), tv%S(i,j,:), pref_mld, rho_c, 1, nz, &
793 tv%eqn_of_state, scale=us%kg_m3_to_R)
796 dz(k) = h(i,j,k) * gv%H_to_Z
802 z_l(k) = z_l(k-1)-dz(k)
818 pe = pe + 0.5*rho_c(k)*(z_u(k)**2-z_l(k)**2)
821 h_ml_tst = h_ml + dz(k)
822 rhodz_ml_tst = rhodz_ml + rho_c(k)*dz(k)
825 rho_ml = rhodz_ml_tst/h_ml_tst
830 pe_mixed_tst = 0.5*rho_ml*(0.**2-(0.-h_ml_tst)**2)
833 if (pe_mixed_tst-pe<=pe_threshold(im))
then 835 rhodz_ml = rhodz_ml_tst
860 cb = -( r1*d1 + r2*(2.*d1) )
862 cc = -( r1*d1*(2*d1) + r2*d )
866 c = d2**2. + d1**2. + 2.*d1*d2
881 gx = 0.5*(ca*x**3+cb*x**2+cc*x+cd)
882 gpx = 0.5*(3*ca*x**2+2*cb*x+cc)
888 ix = 0.5*(ca2*x**2. + cb2*x + cc2)
889 ipx = 0.5*(2.*ca2*x+cb2)
893 fgx = (pe_mixed-(pe+pe_threshold(im)))
894 fpx = (gpx*hx-hpx*gx)*ihx**2+ipx
899 if (abs(fgx)>pe_threshold(im)*pe_threshold_fraction)
then 902 if (x2<0. .or. x2>dz(k))
then 926 if (id_mld(1) > 0)
call post_data(id_mld(1), mld(:,:,1), diagptr)
927 if (id_mld(2) > 0)
call post_data(id_mld(2), mld(:,:,2), diagptr)
928 if (id_mld(3) > 0)
call post_data(id_mld(3), mld(:,:,3), diagptr)
930 end subroutine diagnosemldbyenergy
935 subroutine applyboundaryfluxesinout(CS, G, GV, US, dt, fluxes, optics, nsw, h, tv, &
936 aggregate_FW_forcing, evap_CFL_limit, &
937 minimum_forcing_depth, cTKE, dSV_dT, dSV_dS, &
943 real,
intent(in) :: dt
944 type(
forcing),
intent(inout) :: fluxes
946 integer,
intent(in) :: nsw
948 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
952 logical,
intent(in) :: aggregate_FW_forcing
953 real,
intent(in) :: evap_CFL_limit
955 real,
intent(in) :: minimum_forcing_depth
957 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
958 optional,
intent(out) :: cTKE
960 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
961 optional,
intent(out) :: dSV_dT
963 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
964 optional,
intent(out) :: dSV_dS
966 real,
dimension(SZI_(G),SZJ_(G)), &
967 optional,
intent(out) :: SkinBuoyFlux
970 integer,
parameter :: maxGroundings = 5
971 integer :: numberOfGroundings, iGround(maxgroundings), jGround(maxgroundings)
972 real :: H_limit_fluxes
973 real :: IforcingDepthScale
975 real :: dThickness, dTemp, dSalt
976 real :: fractionOfForcing, hOld, Ithickness
977 real :: RivermixConst
979 real,
dimension(SZI_(G)) :: &
999 real,
dimension(SZI_(G), SZK_(G)) :: &
1005 real,
dimension(SZI_(G)) :: &
1008 real,
dimension(max(nsw,1),SZI_(G)) :: &
1013 real,
dimension(max(nsw,1),SZI_(G),SZK_(G)) :: &
1016 real,
dimension(maxGroundings) :: hGrounding
1017 real :: Temp_in, Salin_in
1022 logical :: calculate_energetics
1023 logical :: calculate_buoyancy
1024 integer,
dimension(2) :: EOSdom
1025 integer :: i, j, is, ie, js, je, k, nz, n, nb
1026 character(len=45) :: mesg
1028 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1032 calculate_energetics = (
present(ctke) .and.
present(dsv_dt) .and.
present(dsv_ds))
1033 calculate_buoyancy =
present(skinbuoyflux)
1034 if (calculate_buoyancy) skinbuoyflux(:,:) = 0.0
1035 if (
present(ctke)) ctke(:,:,:) = 0.0
1036 g_hconv2 = (us%L_to_Z**2*gv%g_Earth * gv%H_to_RZ) * gv%H_to_RZ
1037 eosdom(:) = eos_domain(g%HI)
1040 if (.not.
associated(fluxes%sw) .and. .not.calculate_energetics)
return 1042 if (calculate_buoyancy)
then 1043 surfpressure(:) = 0.0
1044 gorho = us%L_to_Z**2*gv%g_Earth / gv%Rho0
1052 h_limit_fluxes = max(gv%Angstrom_H, 1.e-30*gv%m_to_H)
1055 if (cs%id_createdH>0) cs%createdH(:,:) = 0.
1056 numberofgroundings = 0
1077 do k=1,nz ;
do i=is,ie
1079 t2d(i,k) = tv%T(i,j,k)
1082 if (calculate_energetics)
then 1086 if (
associated(tv%p_surf))
then 1087 do i=is,ie ; pres(i) = tv%p_surf(i,j) ;
enddo 1089 do i=is,ie ; pres(i) = 0.0 ;
enddo 1093 d_pres(i) = (gv%g_Earth * gv%H_to_RZ) * h2d(i,k)
1094 p_lay(i) = pres(i) + 0.5*d_pres(i)
1095 pres(i) = pres(i) + d_pres(i)
1098 dsv_dt(:,j,k), dsv_ds(:,j,k), tv%eqn_of_state, eosdom)
1099 do i=is,ie ; dsv_dt_2d(i,k) = dsv_dt(i,j,k) ;
enddo 1101 pen_tke_2d(:,:) = 0.0
1105 if (.not.
associated(fluxes%sw)) cycle
1107 if (nsw>0)
call extract_optics_slice(optics, j, g, gv, opacity=opacityband, opacity_scale=(1.0/gv%m_to_H))
1152 if (calculate_buoyancy)
then 1153 call extractfluxes1d(g, gv, us, fluxes, optics, nsw, j, dt, &
1154 h_limit_fluxes, cs%use_river_heat_content, cs%use_calving_heat_content, &
1155 h2d, t2d, netmassinout, netmassout, netheat, netsalt, &
1156 pen_sw_bnd, tv, aggregate_fw_forcing, nonpensw=nonpensw, &
1157 net_heat_rate=netheat_rate, net_salt_rate=netsalt_rate, &
1158 netmassinout_rate=netmassinout_rate, pen_sw_bnd_rate=pen_sw_bnd_rate)
1160 call extractfluxes1d(g, gv, us, fluxes, optics, nsw, j, dt, &
1161 h_limit_fluxes, cs%use_river_heat_content, cs%use_calving_heat_content, &
1162 h2d, t2d, netmassinout, netmassout, netheat, netsalt, &
1163 pen_sw_bnd, tv, aggregate_fw_forcing, nonpensw=nonpensw)
1168 if (aggregate_fw_forcing)
then 1169 netmassout(i) = netmassinout(i)
1172 netmassin(i) = netmassinout(i) - netmassout(i)
1174 if (g%mask2dT(i,j)>0.0)
then 1175 fluxes%netMassOut(i,j) = netmassout(i)
1176 fluxes%netMassIn(i,j) = netmassin(i)
1178 fluxes%netMassOut(i,j) = 0.0
1179 fluxes%netMassIn(i,j) = 0.0
1189 if (g%mask2dT(i,j)>0.)
then 1195 dthickness = netmassin(i)
1201 netmassin(i) = netmassin(i) - dthickness
1205 dtemp = dtemp + dthickness*temp_in
1208 if (
associated(fluxes%heat_content_massin)) &
1209 fluxes%heat_content_massin(i,j) = fluxes%heat_content_massin(i,j) + &
1210 t2d(i,k) * max(0.,dthickness) * gv%H_to_RZ * fluxes%C_p * idt
1211 if (
associated(fluxes%heat_content_massout)) &
1212 fluxes%heat_content_massout(i,j) = fluxes%heat_content_massout(i,j) + &
1213 t2d(i,k) * min(0.,dthickness) * gv%H_to_RZ * fluxes%C_p * idt
1214 if (
associated(tv%TempxPmE)) tv%TempxPmE(i,j) = tv%TempxPmE(i,j) + &
1215 t2d(i,k) * dthickness * gv%H_to_RZ
1218 if (calculate_energetics .and.
associated(fluxes%lrunoff) .and. cs%do_rivermix)
then 1229 if (gv%Boussinesq)
then 1230 rivermixconst = -0.5*(cs%rivermix_depth*dt) * ( us%L_to_Z**2*gv%g_Earth ) * gv%Rho0
1232 rivermixconst = -0.5*(cs%rivermix_depth*dt) * gv%Rho0 * ( us%L_to_Z**2*gv%g_Earth )
1234 ctke(i,j,k) = ctke(i,j,k) + max(0.0, rivermixconst*dsv_ds(i,j,1) * &
1235 (fluxes%lrunoff(i,j) + fluxes%frunoff(i,j)) * tv%S(i,j,1))
1240 h2d(i,k) = h2d(i,k) + dthickness
1241 if (h2d(i,k) > 0.0)
then 1242 if (calculate_energetics .and. (dthickness > 0.))
then 1245 ctke(i,j,k) = ctke(i,j,k) + 0.5*g_hconv2*(hold*dthickness) * &
1246 ((t2d(i,k) - temp_in) * dsv_dt(i,j,k) + (tv%S(i,j,k) - salin_in) * dsv_ds(i,j,k))
1248 ithickness = 1.0/h2d(i,k)
1250 if (dthickness /= 0. .or. dtemp /= 0.) t2d(i,k) = (hold*t2d(i,k) + dtemp)*ithickness
1251 if (dthickness /= 0. .or. dsalt /= 0.) tv%S(i,j,k) = (hold*tv%S(i,j,k) + dsalt)*ithickness
1263 iforcingdepthscale = 1. / max(gv%H_subroundoff, minimum_forcing_depth - netmassout(i) )
1265 fractionofforcing = min(1.0, h2d(i,k)*iforcingdepthscale)
1270 if (-fractionofforcing*netmassout(i) > evap_cfl_limit*h2d(i,k))
then 1271 fractionofforcing = -evap_cfl_limit*h2d(i,k)/netmassout(i)
1276 dthickness = max( fractionofforcing*netmassout(i), -h2d(i,k) )
1277 dtemp = fractionofforcing*netheat(i)
1279 dsalt = max( fractionofforcing*netsalt(i), -0.9999*h2d(i,k)*tv%S(i,j,k))
1283 netmassout(i) = netmassout(i) - dthickness
1284 netheat(i) = netheat(i) - dtemp
1285 netsalt(i) = netsalt(i) - dsalt
1288 dtemp = dtemp + dthickness*t2d(i,k)
1291 if (
associated(fluxes%heat_content_massin)) &
1292 fluxes%heat_content_massin(i,j) = fluxes%heat_content_massin(i,j) + &
1293 t2d(i,k) * max(0.,dthickness) * gv%H_to_RZ * fluxes%C_p * idt
1294 if (
associated(fluxes%heat_content_massout)) &
1295 fluxes%heat_content_massout(i,j) = fluxes%heat_content_massout(i,j) + &
1296 t2d(i,k) * min(0.,dthickness) * gv%H_to_RZ * fluxes%C_p * idt
1297 if (
associated(tv%TempxPmE)) tv%TempxPmE(i,j) = tv%TempxPmE(i,j) + &
1298 t2d(i,k) * dthickness * gv%H_to_RZ
1302 h2d(i,k) = h2d(i,k) + dthickness
1304 if (h2d(i,k) > 0.)
then 1305 if (calculate_energetics)
then 1311 ctke(i,j,k) = ctke(i,j,k) - (0.5*h2d(i,k)*g_hconv2) * &
1312 ((dtemp - dthickness*t2d(i,k)) * dsv_dt(i,j,k) + &
1313 (dsalt - dthickness*tv%S(i,j,k)) * dsv_ds(i,j,k))
1315 ithickness = 1.0/h2d(i,k)
1316 t2d(i,k) = (hold*t2d(i,k) + dtemp)*ithickness
1317 tv%S(i,j,k) = (hold*tv%S(i,j,k) + dsalt)*ithickness
1318 elseif (h2d(i,k) < 0.0)
then 1319 call forcing_singlepointprint(fluxes,g,i,j,
'applyBoundaryFluxesInOut (h<0)')
1320 write(0,*)
'applyBoundaryFluxesInOut(): lon,lat=',g%geoLonT(i,j),g%geoLatT(i,j)
1321 write(0,*)
'applyBoundaryFluxesInOut(): netT,netS,netH=',netheat(i),netsalt(i),netmassinout(i)
1322 write(0,*)
'applyBoundaryFluxesInOut(): dT,dS,dH=',dtemp,dsalt,dthickness
1323 write(0,*)
'applyBoundaryFluxesInOut(): h(n),h(n+1),k=',hold,h2d(i,k),k
1324 call mom_error(fatal,
"MOM_diabatic_driver.F90, applyBoundaryFluxesInOut(): "//&
1325 "Complete mass loss in column!")
1331 elseif ((abs(netheat(i))+abs(netsalt(i))+abs(netmassin(i))+abs(netmassout(i)))>0.)
then 1333 if (.not. cs%ignore_fluxes_over_land)
then 1334 call forcing_singlepointprint(fluxes,g,i,j,
'applyBoundaryFluxesInOut (land)')
1335 write(0,*)
'applyBoundaryFluxesInOut(): lon,lat=',g%geoLonT(i,j),g%geoLatT(i,j)
1336 write(0,*)
'applyBoundaryFluxesInOut(): netHeat,netSalt,netMassIn,netMassOut=',&
1337 netheat(i),netsalt(i),netmassin(i),netmassout(i)
1339 call mom_error(fatal,
"MOM_diabatic_driver.F90, applyBoundaryFluxesInOut(): "//&
1340 "Mass loss over land?")
1346 if (netmassin(i)+netmassout(i) /= 0.0)
then 1348 numberofgroundings = numberofgroundings +1
1349 if (numberofgroundings<=maxgroundings)
then 1350 iground(numberofgroundings) = i
1351 jground(numberofgroundings) = j
1352 hgrounding(numberofgroundings) = netmassin(i)+netmassout(i)
1355 if (cs%id_createdH>0) cs%createdH(i,j) = cs%createdH(i,j) - (netmassin(i)+netmassout(i))/dt
1366 if (cs%id_penSW_diag > 0 .or. cs%id_penSWflux_diag > 0)
then 1367 do k=1,nz ;
do i=is,ie
1368 cs%penSW_diag(i,j,k) = t2d(i,k)
1369 cs%penSWflux_diag(i,j,k) = 0.0
1372 cs%penSWflux_diag(i,j,k) = 0.0
1376 if (calculate_energetics)
then 1377 call absorbremainingsw(g, gv, us, h2d, opacityband, nsw, optics, j, dt, h_limit_fluxes, &
1378 .false., .true., t2d, pen_sw_bnd, tke=pen_tke_2d, dsv_dt=dsv_dt_2d)
1380 do k=1,nz ;
do i=is,ie
1381 ctke(i,j,k) = ctke(i,j,k) + pen_tke_2d(i,k)
1384 call absorbremainingsw(g, gv, us, h2d, opacityband, nsw, optics, j, dt, h_limit_fluxes, &
1385 .false., .true., t2d, pen_sw_bnd)
1391 do k=1,nz ;
do i=is,ie
1393 tv%T(i,j,k) = t2d(i,k)
1398 if (cs%id_penSW_diag > 0 .or. cs%id_penSWflux_diag > 0)
then 1401 do k=1,nz ;
do i=is,ie
1402 cs%penSW_diag(i,j,k) = (t2d(i,k)-cs%penSW_diag(i,j,k))*h(i,j,k) * idt * tv%C_p * gv%H_to_RZ
1411 if (cs%id_penSWflux_diag > 0)
then 1412 do k=nz,1,-1 ;
do i=is,ie
1413 cs%penSWflux_diag(i,j,k) = cs%penSW_diag(i,j,k) + cs%penSWflux_diag(i,j,k+1)
1420 if (cs%id_nonpenSW_diag > 0)
then 1422 cs%nonpenSW_diag(i,j) = nonpensw(i) * idt * tv%C_p * gv%H_to_RZ
1431 if (calculate_buoyancy)
then 1434 netpen_rate(:) = 0.0
1441 do i=is,ie ;
do nb=1,nsw ; netpen_rate(i) = netpen_rate(i) + pen_sw_bnd_rate(nb,i) ;
enddo ;
enddo 1444 if (
associated(tv%p_surf))
then ;
do i=is,ie ; surfpressure(i) = tv%p_surf(i,j) ;
enddo ;
endif 1446 tv%eqn_of_state, eosdom)
1453 skinbuoyflux(i,j) = - gorho * gv%H_to_Z * &
1454 (drhods(i) * (netsalt_rate(i) - tv%S(i,j,1)*netmassinout_rate(i)) + &
1455 drhodt(i) * ( netheat_rate(i) + netpen_rate(i)) )
1462 if (cs%id_createdH > 0)
call post_data(cs%id_createdH , cs%createdH , cs%diag)
1463 if (cs%id_penSW_diag > 0)
call post_data(cs%id_penSW_diag , cs%penSW_diag , cs%diag)
1464 if (cs%id_penSWflux_diag > 0)
call post_data(cs%id_penSWflux_diag, cs%penSWflux_diag, cs%diag)
1465 if (cs%id_nonpenSW_diag > 0)
call post_data(cs%id_nonpenSW_diag , cs%nonpenSW_diag , cs%diag)
1468 if (numberofgroundings>0 .and. .not. cs%ignore_fluxes_over_land)
then 1469 do i = 1, min(numberofgroundings, maxgroundings)
1470 call forcing_singlepointprint(fluxes,g,iground(i),jground(i),
'applyBoundaryFluxesInOut (grounding)')
1471 write(mesg(1:45),
'(3es15.3)') g%geoLonT( iground(i), jground(i) ), &
1472 g%geoLatT( iground(i), jground(i)), hgrounding(i)*gv%H_to_m
1473 call mom_error(warning,
"MOM_diabatic_driver.F90, applyBoundaryFluxesInOut(): "//&
1474 "Mass created. x,y,dh= "//trim(mesg), all_print=.true.)
1477 if (numberofgroundings - maxgroundings > 0)
then 1478 write(mesg,
'(i4)') numberofgroundings - maxgroundings
1479 call mom_error(warning,
"MOM_diabatic_driver:F90, applyBoundaryFluxesInOut(): "//&
1480 trim(mesg) //
" groundings remaining")
1484 end subroutine applyboundaryfluxesinout
1487 subroutine diabatic_aux_init(Time, G, GV, US, param_file, diag, CS, useALEalgorithm, use_ePBL)
1488 type(time_type),
target,
intent(in) :: Time
1493 type(
diag_ctrl),
target,
intent(inout) :: diag
1496 logical,
intent(in) :: useALEalgorithm
1498 logical,
intent(in) :: use_ePBL
1503 #include "version_variable.h" 1504 character(len=40) :: mdl =
"MOM_diabatic_aux" 1505 character(len=48) :: thickness_units
1506 character(len=200) :: inputdir
1507 character(len=240) :: chl_filename
1508 character(len=128) :: chl_file
1510 character(len=32) :: chl_varname
1511 logical :: use_temperature
1512 integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, nz, nbands
1513 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
1514 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1516 if (
associated(cs))
then 1517 call mom_error(warning,
"diabatic_aux_init called with an "// &
1518 "associated control structure.")
1529 "The following parameters are used for auxiliary diabatic processes.")
1531 call get_param(param_file, mdl,
"ENABLE_THERMODYNAMICS", use_temperature, &
1532 "If true, temperature and salinity are used as state "//&
1533 "variables.", default=.true.)
1535 call get_param(param_file, mdl,
"RECLAIM_FRAZIL", cs%reclaim_frazil, &
1536 "If true, try to use any frazil heat deficit to cool any "//&
1537 "overlying layers down to the freezing point, thereby "//&
1538 "avoiding the creation of thin ice when the SST is above "//&
1539 "the freezing point.", default=.true.)
1540 call get_param(param_file, mdl,
"PRESSURE_DEPENDENT_FRAZIL", &
1541 cs%pressure_dependent_frazil, &
1542 "If true, use a pressure dependent freezing temperature "//&
1543 "when making frazil. The default is false, which will be "//&
1544 "faster but is inappropriate with ice-shelf cavities.", &
1548 call get_param(param_file, mdl,
"IGNORE_FLUXES_OVER_LAND", cs%ignore_fluxes_over_land,&
1549 "If true, the model does not check if fluxes are being applied "//&
1550 "over land points. This is needed when the ocean is coupled "//&
1551 "with ice shelves and sea ice, since the sea ice mask needs to "//&
1552 "be different than the ocean mask to avoid sea ice formation "//&
1553 "under ice shelves. This flag only works when use_ePBL = True.", default=.false.)
1554 call get_param(param_file, mdl,
"DO_RIVERMIX", cs%do_rivermix, &
1555 "If true, apply additional mixing wherever there is "//&
1556 "runoff, so that it is mixed down to RIVERMIX_DEPTH "//&
1557 "if the ocean is that deep.", default=.false.)
1558 if (cs%do_rivermix) &
1559 call get_param(param_file, mdl,
"RIVERMIX_DEPTH", cs%rivermix_depth, &
1560 "The depth to which rivers are mixed if DO_RIVERMIX is "//&
1561 "defined.", units=
"m", default=0.0, scale=us%m_to_Z)
1563 cs%do_rivermix = .false. ; cs%rivermix_depth = 0.0 ; cs%ignore_fluxes_over_land = .false.
1566 if (gv%nkml == 0)
then 1567 call get_param(param_file, mdl,
"USE_RIVER_HEAT_CONTENT", cs%use_river_heat_content, &
1568 "If true, use the fluxes%runoff_Hflx field to set the "//&
1569 "heat carried by runoff, instead of using SST*CP*liq_runoff.", &
1571 call get_param(param_file, mdl,
"USE_CALVING_HEAT_CONTENT", cs%use_calving_heat_content, &
1572 "If true, use the fluxes%calving_Hflx field to set the "//&
1573 "heat carried by runoff, instead of using SST*CP*froz_runoff.", &
1576 cs%use_river_heat_content = .false.
1577 cs%use_calving_heat_content = .false.
1580 if (usealealgorithm)
then 1581 cs%id_createdH = register_diag_field(
'ocean_model',
"created_H",diag%axesT1, &
1582 time,
"The volume flux added to stop the ocean from drying out and becoming negative in depth", &
1583 "m s-1", conversion=gv%H_to_m*us%s_to_T)
1584 if (cs%id_createdH>0)
allocate(cs%createdH(isd:ied,jsd:jed))
1587 cs%id_penSW_diag = register_diag_field(
'ocean_model',
'rsdoabsorb', &
1588 diag%axesTL, time,
'Convergence of Penetrative Shortwave Flux in Sea Water Layer',&
1589 'W m-2', conversion=us%QRZ_T_to_W_m2, &
1590 standard_name=
'net_rate_of_absorption_of_shortwave_energy_in_ocean_layer', v_extensive=.true.)
1594 cs%id_penSWflux_diag = register_diag_field(
'ocean_model',
'rsdo', &
1595 diag%axesTi, time,
'Downwelling Shortwave Flux in Sea Water at Grid Cell Upper Interface',&
1596 'W m-2', conversion=us%QRZ_T_to_W_m2, standard_name=
'downwelling_shortwave_flux_in_sea_water')
1599 if (cs%id_penSW_diag>0 .or. cs%id_penSWflux_diag>0)
then 1600 allocate(cs%penSW_diag(isd:ied,jsd:jed,nz)) ; cs%penSW_diag(:,:,:) = 0.0
1601 allocate(cs%penSWflux_diag(isd:ied,jsd:jed,nz+1)) ; cs%penSWflux_diag(:,:,:) = 0.0
1605 cs%id_nonpenSW_diag = register_diag_field(
'ocean_model',
'nonpenSW', &
1606 diag%axesT1, time, &
1607 'Non-downwelling SW radiation (i.e., SW absorbed in ocean surface with LW,SENS,LAT)',&
1608 'W m-2', conversion=us%QRZ_T_to_W_m2, &
1609 standard_name=
'nondownwelling_shortwave_flux_in_sea_water')
1610 if (cs%id_nonpenSW_diag > 0)
then 1611 allocate(cs%nonpenSW_diag(isd:ied,jsd:jed)) ; cs%nonpenSW_diag(:,:) = 0.0
1615 if (use_temperature)
then 1616 call get_param(param_file, mdl,
"VAR_PEN_SW", cs%var_pen_sw, &
1617 "If true, use one of the CHL_A schemes specified by "//&
1618 "OPACITY_SCHEME to determine the e-folding depth of "//&
1619 "incoming short wave radiation.", default=.false.)
1620 if (cs%var_pen_sw)
then 1622 call get_param(param_file, mdl,
"CHL_FROM_FILE", cs%chl_from_file, &
1623 "If true, chl_a is read from a file.", default=.true.)
1624 if (cs%chl_from_file)
then 1625 call time_interp_external_init()
1627 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
1628 call get_param(param_file, mdl,
"CHL_FILE", chl_file, &
1629 "CHL_FILE is the file containing chl_a concentrations in "//&
1630 "the variable CHL_A. It is used when VAR_PEN_SW and "//&
1631 "CHL_FROM_FILE are true.", fail_if_missing=.true.)
1632 chl_filename = trim(slasher(inputdir))//trim(chl_file)
1633 call log_param(param_file, mdl,
"INPUTDIR/CHL_FILE", chl_filename)
1634 call get_param(param_file, mdl,
"CHL_VARNAME", chl_varname, &
1635 "Name of CHL_A variable in CHL_FILE.", default=
'CHL_A')
1636 cs%sbc_chl = init_external_field(chl_filename, trim(chl_varname), domain=g%Domain%mpp_domain)
1639 cs%id_chl = register_diag_field(
'ocean_model',
'Chl_opac', diag%axesT1, time, &
1640 'Surface chlorophyll A concentration used to find opacity',
'mg m-3')
1644 id_clock_uv_at_h = cpu_clock_id(
'(Ocean find_uv_at_h)', grain=clock_routine)
1645 id_clock_frazil = cpu_clock_id(
'(Ocean frazil)', grain=clock_routine)
1647 end subroutine diabatic_aux_init
1651 subroutine diabatic_aux_end(CS)
1655 if (.not.
associated(cs))
return 1657 if (cs%id_createdH >0)
deallocate(cs%createdH)
1658 if (cs%id_penSW_diag >0)
deallocate(cs%penSW_diag)
1659 if (cs%id_penSWflux_diag >0)
deallocate(cs%penSWflux_diag)
1660 if (cs%id_nonpenSW_diag >0)
deallocate(cs%nonpenSW_diag)
1662 if (
associated(cs))
deallocate(cs)
1664 end subroutine diabatic_aux_end
This module implements boundary forcing for MOM6.
Calculate the derivatives of specific volume with temperature and salinity from T, S, and P.
Provides functions for some diabatic processes such as fraxil, brine rejection, tendency due to surfa...
Ocean grid type. See mom_grid for details.
Calculates density of sea water from T, S and P.
Calculates the freezing point of sea water from T, S and P.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Vertical viscosities, drag coefficients, and related fields.
Wraps the MPP cpu clock functions.
This module contains I/O framework code.
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Orchestrates the registration and calling of tracer packages.
Make a diagnostic available for averaging or output.
The control structure with paramters for the MOM_opacity module.
Describes various unit conversion factors.
Control structure for diabatic_aux.
Routines used to calculate the opacity of the ocean.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
The control structure for orchestrating the calling of tracer packages.
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Provides subroutines for quantities specific to the equation of state.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
This type is used to store information about ocean optical properties.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.