16 implicit none ;
private 18 #include <MOM_memory.h> 20 public set_opacity, opacity_init, opacity_end, opacity_manizza, opacity_morel
21 public extract_optics_slice, extract_optics_fields, optics_nbands
22 public absorbremainingsw, sumswoverbands
28 real,
pointer,
dimension(:,:,:,:) :: opacity_band => null()
31 real,
pointer,
dimension(:,:,:) :: sw_pen_band => null()
35 real,
pointer,
dimension(:) :: &
36 min_wavelength_band => null(), &
37 max_wavelength_band => null()
39 real :: pensw_flux_absorb
41 real :: pensw_absorb_invlen
43 logical :: answers_2018
53 integer :: opacity_scheme
58 real :: pen_sw_scale_2nd
60 real :: sw_1st_exp_ratio
65 real :: opacity_land_value
71 integer :: id_sw_pen = -1, id_sw_vis_pen = -1
72 integer,
pointer :: id_opacity(:) => null()
77 integer,
parameter :: no_scheme = 0, manizza_05 = 1, morel_88 = 2, single_exp = 3, double_exp = 4
80 character*(10),
parameter :: manizza_05_string =
"MANIZZA_05" 81 character*(10),
parameter :: morel_88_string =
"MOREL_88" 82 character*(10),
parameter :: single_exp_string =
"SINGLE_EXP" 83 character*(10),
parameter :: double_exp_string =
"DOUBLE_EXP" 85 real,
parameter :: op_diag_len = 1e-10
91 subroutine set_opacity(optics, sw_total, sw_vis_dir, sw_vis_dif, sw_nir_dir, sw_nir_dif, &
92 G, GV, US, CS, chl_2d, chl_3d)
95 real,
dimension(:,:),
pointer :: sw_total
96 real,
dimension(:,:),
pointer :: sw_vis_dir
97 real,
dimension(:,:),
pointer :: sw_vis_dif
98 real,
dimension(:,:),
pointer :: sw_nir_dir
99 real,
dimension(:,:),
pointer :: sw_nir_dif
104 real,
dimension(SZI_(G),SZJ_(G)), &
105 optional,
intent(in) :: chl_2d
106 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
107 optional,
intent(in) :: chl_3d
110 integer :: i, j, k, n, is, ie, js, je, nz
111 real :: inv_sw_pen_scale
114 logical :: call_for_surface
115 real :: tmp(szi_(g),szj_(g),szk_(gv))
116 real :: chl(szi_(g),szj_(g),szk_(gv))
117 real :: pen_sw_tot(szi_(g),szj_(g))
119 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
121 if (.not.
associated(cs))
call mom_error(fatal,
"set_opacity: "// &
122 "Module must be initialized via opacity_init before it is used.")
124 if (
present(chl_2d) .or.
present(chl_3d))
then 126 call opacity_from_chl(optics, sw_total, sw_vis_dir, sw_vis_dif, sw_nir_dir, sw_nir_dif, &
127 g, gv, us, cs, chl_2d, chl_3d)
129 if (optics%nbands <= 1)
then ; inv_nbands = 1.0
130 else ; inv_nbands = 1.0 /
real(optics%nbands) ; endif
133 inv_sw_pen_scale = 1.0 / max(cs%pen_sw_scale, 0.1*gv%Angstrom_m, &
134 gv%H_to_m*gv%H_subroundoff)
135 if ( cs%Opacity_scheme == double_exp )
then 137 do k=1,nz ;
do j=js,je ;
do i=is,ie
138 optics%opacity_band(1,i,j,k) = inv_sw_pen_scale
139 optics%opacity_band(2,i,j,k) = 1.0 / max(cs%pen_sw_scale_2nd, &
140 0.1*gv%Angstrom_m,gv%H_to_m*gv%H_subroundoff)
141 enddo ;
enddo ;
enddo 142 if (.not.
associated(sw_total) .or. (cs%pen_SW_scale <= 0.0))
then 144 do j=js,je ;
do i=is,ie ;
do n=1,optics%nbands
145 optics%sw_pen_band(n,i,j) = 0.0
146 enddo ;
enddo ;
enddo 149 do j=js,je ;
do i=is,ie
150 optics%sw_pen_band(1,i,j) = (cs%SW_1st_EXP_RATIO) * sw_total(i,j)
151 optics%sw_pen_band(2,i,j) = (1.-cs%SW_1st_EXP_RATIO) * sw_total(i,j)
155 do k=1,nz ;
do j=js,je ;
do i=is,ie ;
do n=1,optics%nbands
156 optics%opacity_band(n,i,j,k) = inv_sw_pen_scale
157 enddo ;
enddo ;
enddo ;
enddo 158 if (.not.
associated(sw_total) .or. (cs%pen_SW_scale <= 0.0))
then 160 do j=js,je ;
do i=is,ie ;
do n=1,optics%nbands
161 optics%sw_pen_band(n,i,j) = 0.0
162 enddo ;
enddo ;
enddo 165 do j=js,je ;
do i=is,ie ;
do n=1,optics%nbands
166 optics%sw_pen_band(n,i,j) = cs%pen_SW_frac * inv_nbands * sw_total(i,j)
167 enddo ;
enddo ;
enddo 172 if (query_averaging_enabled(cs%diag))
then 173 if (cs%id_sw_pen > 0)
then 175 do j=js,je ;
do i=is,ie
176 pen_sw_tot(i,j) = 0.0
178 pen_sw_tot(i,j) = pen_sw_tot(i,j) + optics%sw_pen_band(n,i,j)
181 call post_data(cs%id_sw_pen, pen_sw_tot, cs%diag)
183 if (cs%id_sw_vis_pen > 0)
then 184 if (cs%opacity_scheme == manizza_05)
then 186 do j=js,je ;
do i=is,ie
187 pen_sw_tot(i,j) = 0.0
188 do n=1,min(optics%nbands,2)
189 pen_sw_tot(i,j) = pen_sw_tot(i,j) + optics%sw_pen_band(n,i,j)
194 do j=js,je ;
do i=is,ie
195 pen_sw_tot(i,j) = 0.0
197 pen_sw_tot(i,j) = pen_sw_tot(i,j) + optics%sw_pen_band(n,i,j)
201 call post_data(cs%id_sw_vis_pen, pen_sw_tot, cs%diag)
203 do n=1,optics%nbands ;
if (cs%id_opacity(n) > 0)
then 205 do k=1,nz ;
do j=js,je ;
do i=is,ie
209 tmp(i,j,k) = tanh(op_diag_len * optics%opacity_band(n,i,j,k)) / op_diag_len
210 enddo ;
enddo ;
enddo 211 call post_data(cs%id_opacity(n), tmp, cs%diag)
215 end subroutine set_opacity
220 subroutine opacity_from_chl(optics, sw_total, sw_vis_dir, sw_vis_dif, sw_nir_dir, sw_nir_dif, &
221 G, GV, US, CS, chl_2d, chl_3d)
222 type(optics_type),
intent(inout) :: optics
224 real,
dimension(:,:),
pointer :: sw_total
225 real,
dimension(:,:),
pointer :: sw_vis_dir
226 real,
dimension(:,:),
pointer :: sw_vis_dif
227 real,
dimension(:,:),
pointer :: sw_nir_dir
228 real,
dimension(:,:),
pointer :: sw_nir_dif
229 type(ocean_grid_type),
intent(in) :: G
230 type(verticalGrid_type),
intent(in) :: GV
231 type(unit_scale_type),
intent(in) :: US
232 type(opacity_CS),
pointer :: CS
233 real,
dimension(SZI_(G),SZJ_(G)), &
234 optional,
intent(in) :: chl_2d
235 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
236 optional,
intent(in) :: chl_3d
238 real :: chl_data(SZI_(G),SZJ_(G))
241 real :: Inv_nbands_nir
249 type(time_type) :: day
250 character(len=128) :: mesg
251 integer :: i, j, k, n, is, ie, js, je, nz, nbands
252 logical :: multiband_vis_input, multiband_nir_input
254 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
268 nbands = optics%nbands
270 if (nbands <= 1)
then ; inv_nbands = 1.0
271 else ; inv_nbands = 1.0 /
real(nbands) ; endif
273 if (nbands <= 2)
then ; inv_nbands_nir = 0.0
274 else ; inv_nbands_nir = 1.0 /
real(nbands - 2.0) ; endif
276 multiband_vis_input = (
associated(sw_vis_dir) .and. &
277 associated(sw_vis_dif))
278 multiband_nir_input = (
associated(sw_nir_dir) .and. &
279 associated(sw_nir_dif))
282 if (
present(chl_3d))
then 283 do j=js,je ;
do i=is,ie ; chl_data(i,j) = chl_3d(i,j,1) ;
enddo ;
enddo 284 do k=1,nz;
do j=js,je ;
do i=is,ie
285 if ((g%mask2dT(i,j) > 0.5) .and. (chl_3d(i,j,k) < 0.0))
then 286 write(mesg,
'(" Negative chl_3d of ",(1pe12.4)," found at i,j,k = ", & 287 & 3(1x,i3), " lon/lat = ",(1pe12.4)," E ", (1pe12.4), " N.")') &
288 chl_3d(i,j,k), i, j, k, g%geoLonT(i,j), g%geoLatT(i,j)
289 call mom_error(fatal,
"MOM_opacity opacity_from_chl: "//trim(mesg))
291 enddo ;
enddo ;
enddo 292 elseif (
present(chl_2d))
then 293 do j=js,je ;
do i=is,ie ; chl_data(i,j) = chl_2d(i,j) ;
enddo ;
enddo 294 do j=js,je ;
do i=is,ie
295 if ((g%mask2dT(i,j) > 0.5) .and. (chl_2d(i,j) < 0.0))
then 296 write(mesg,
'(" Negative chl_2d of ",(1pe12.4)," at i,j = ", & 297 & 2(i3), "lon/lat = ",(1pe12.4)," E ", (1pe12.4), " N.")') &
298 chl_data(i,j), i, j, g%geoLonT(i,j), g%geoLatT(i,j)
299 call mom_error(fatal,
"MOM_opacity opacity_from_chl: "//trim(mesg))
303 call mom_error(fatal,
"Either chl_2d or chl_3d must be preesnt in a call to opacity_form_chl.")
306 select case (cs%opacity_scheme)
309 do j=js,je ;
do i=is,ie
310 sw_vis_tot = 0.0 ; sw_nir_tot = 0.0
311 if (g%mask2dT(i,j) > 0.5)
then 312 if (multiband_vis_input)
then 313 sw_vis_tot = sw_vis_dir(i,j) + sw_vis_dif(i,j)
315 sw_vis_tot = 0.42 * sw_total(i,j)
317 if (multiband_nir_input)
then 318 sw_nir_tot = sw_nir_dir(i,j) + sw_nir_dif(i,j)
320 sw_nir_tot = sw_total(i,j) - sw_vis_tot
325 optics%sw_pen_band(1,i,j) = cs%blue_frac*sw_vis_tot
328 optics%sw_pen_band(2,i,j) = (1.0-cs%blue_frac)*sw_vis_tot
331 optics%sw_pen_band(n,i,j) = inv_nbands_nir * sw_nir_tot
336 do j=js,je ;
do i=is,ie
338 if (g%mask2dT(i,j) > 0.5)
then ;
if (multiband_vis_input)
then 339 sw_pen_tot = sw_pen_frac_morel(chl_data(i,j)) * (sw_vis_dir(i,j) + sw_vis_dif(i,j))
341 sw_pen_tot = sw_pen_frac_morel(chl_data(i,j)) * 0.5*sw_total(i,j)
345 optics%sw_pen_band(n,i,j) = inv_nbands*sw_pen_tot
349 call mom_error(fatal,
"opacity_from_chl: CS%opacity_scheme is not valid.")
354 if (
present(chl_3d))
then 355 do j=js,je ;
do i=is,ie ; chl_data(i,j) = chl_3d(i,j,k) ;
enddo ;
enddo 358 select case (cs%opacity_scheme)
360 do j=js,je ;
do i=is,ie
361 if (g%mask2dT(i,j) <= 0.5)
then 363 optics%opacity_band(n,i,j,k) = cs%opacity_land_value
367 optics%opacity_band(1,i,j,k) = 0.0232 + 0.074*chl_data(i,j)**0.674
369 optics%opacity_band(2,i,j,k) = 0.225 + 0.037*chl_data(i,j)**0.629
371 do n=3,nbands ; optics%opacity_band(n,i,j,k) = 2.86 ;
enddo 375 do j=js,je ;
do i=is,ie
376 optics%opacity_band(1,i,j,k) = cs%opacity_land_value
377 if (g%mask2dT(i,j) > 0.5) &
378 optics%opacity_band(1,i,j,k) = opacity_morel(chl_data(i,j))
381 optics%opacity_band(n,i,j,k) = optics%opacity_band(1,i,j,k)
386 call mom_error(fatal,
"opacity_from_chl: CS%opacity_scheme is not valid.")
391 end subroutine opacity_from_chl
395 function opacity_morel(chl_data)
396 real,
intent(in) :: chl_data
397 real :: opacity_morel
404 real,
dimension(6),
parameter :: &
405 z2_coef=(/7.925, -6.644, 3.662, -1.815, -0.218, 0.502/)
408 chl = log10(min(max(chl_data,0.02),60.0)) ; chl2 = chl*chl
409 opacity_morel = 1.0 / ( (z2_coef(1) + z2_coef(2)*chl) + chl2 * &
410 ((z2_coef(3) + chl*z2_coef(4)) + chl2*(z2_coef(5) + chl*z2_coef(6))) )
415 function sw_pen_frac_morel(chl_data)
416 real,
intent(in) :: chl_data
417 real :: sw_pen_frac_morel
425 real,
dimension(6),
parameter :: &
426 v1_coef=(/0.321, 0.008, 0.132, 0.038, -0.017, -0.007/)
428 chl = log10(min(max(chl_data,0.02),60.0)) ; chl2 = chl*chl
429 sw_pen_frac_morel = 1.0 - ( (v1_coef(1) + v1_coef(2)*chl) + chl2 * &
430 ((v1_coef(3) + chl*v1_coef(4)) + chl2*(v1_coef(5) + chl*v1_coef(6))) )
431 end function sw_pen_frac_morel
435 function opacity_manizza(chl_data)
436 real,
intent(in) :: chl_data
437 real :: opacity_manizza
440 opacity_manizza = 0.0232 + 0.074*chl_data**0.674
445 subroutine extract_optics_slice(optics, j, G, GV, opacity, opacity_scale, penSW_top, penSW_scale)
448 integer,
intent(in) :: j
451 real,
dimension(max(optics%nbands,1),SZI_(G),SZK_(GV)), &
452 optional,
intent(out) :: opacity
453 real,
optional,
intent(in) :: opacity_scale
454 real,
dimension(max(optics%nbands,1),SZI_(G)), &
455 optional,
intent(out) :: pensw_top
458 real,
optional,
intent(in) :: pensw_scale
461 real :: scale_opacity, scale_pensw
462 integer :: i, is, ie, k, nz, n
463 is = g%isc ; ie = g%iec ; nz = g%ke
465 scale_opacity = 1.0 ;
if (
present(opacity_scale)) scale_opacity = opacity_scale
466 scale_pensw = 1.0 ;
if (
present(pensw_scale)) scale_pensw = pensw_scale
468 if (
present(opacity))
then ;
do k=1,nz ;
do i=is,ie
470 opacity(n,i,k) = scale_opacity * optics%opacity_band(n,i,j,k)
472 enddo ;
enddo ;
endif 474 if (
present(pensw_top))
then ;
do k=1,nz ;
do i=is,ie
476 pensw_top(n,i) = scale_pensw * optics%sw_pen_band(n,i,j)
478 enddo ;
enddo ;
endif 480 end subroutine extract_optics_slice
483 subroutine extract_optics_fields(optics, nbands)
486 integer,
optional,
intent(out) :: nbands
488 if (
present(nbands)) nbands = optics%nbands
490 end subroutine extract_optics_fields
493 function optics_nbands(optics)
496 integer :: optics_nbands
498 optics_nbands = optics%nbands
499 end function optics_nbands
508 subroutine absorbremainingsw(G, GV, US, h, opacity_band, nsw, optics, j, dt, H_limit_fluxes, &
509 adjustAbsorptionProfile, absorbAllSW, T, Pen_SW_bnd, &
510 eps, ksort, htot, Ttot, TKE, dSV_dT)
515 integer,
intent(in) :: nsw
517 real,
dimension(SZI_(G),SZK_(GV)),
intent(in) :: h
518 real,
dimension(max(1,nsw),SZI_(G),SZK_(GV)),
intent(in) :: opacity_band
523 integer,
intent(in) :: j
524 real,
intent(in) :: dt
525 real,
intent(in) :: h_limit_fluxes
531 logical,
intent(in) :: adjustabsorptionprofile
537 logical,
intent(in) :: absorballsw
543 real,
dimension(SZI_(G),SZK_(GV)),
intent(inout) :: t
545 real,
dimension(max(1,nsw),SZI_(G)),
intent(inout) :: pen_sw_bnd
550 real,
dimension(SZI_(G),SZK_(GV)),
optional,
intent(in) :: eps
553 integer,
dimension(SZI_(G),SZK_(GV)),
optional,
intent(in) :: ksort
554 real,
dimension(SZI_(G)),
optional,
intent(in) :: htot
555 real,
dimension(SZI_(G)),
optional,
intent(inout) :: ttot
557 real,
dimension(SZI_(G),SZK_(GV)),
optional,
intent(in) :: dsv_dt
559 real,
dimension(SZI_(G),SZK_(GV)),
optional,
intent(inout) :: tke
563 real,
dimension(SZI_(G),SZK_(GV)) :: &
571 real,
dimension(SZI_(G)) :: &
602 logical :: sw_remains
607 integer :: is, ie, nz, i, k, ks, n
610 min_sw_heat = optics%PenSW_flux_absorb * dt
611 i_habs = optics%PenSW_absorb_Invlen
613 h_min_heat = 2.0*gv%Angstrom_H + gv%H_subroundoff
614 is = g%isc ; ie = g%iec ; nz = g%ke
615 c1_6 = 1.0 / 6.0 ; c1_60 = 1.0 / 60.0
617 tke_calc = (
present(tke) .and.
present(dsv_dt))
619 if (optics%answers_2018)
then 620 g_hconv2 = (us%L_to_Z**2*gv%g_Earth * gv%H_to_RZ) * gv%H_to_RZ
622 g_hconv2 = us%L_to_Z**2*gv%g_Earth * gv%H_to_RZ**2
626 if (
present(htot))
then ;
do i=is,ie ; h_heat(i) = htot(i) ;
enddo ;
endif 630 do ks=1,nz ;
do i=is,ie
632 if (
present(ksort))
then 633 if (ksort(i,ks) <= 0) cycle
636 epsilon = 0.0 ;
if (
present(eps)) epsilon = eps(i,k)
638 t_chg_above(i,k) = 0.0
640 if (h(i,k) > 1.5*epsilon)
then 641 do n=1,nsw ;
if (pen_sw_bnd(n,i) > 0.0)
then 643 opt_depth = h(i,k) * opacity_band(n,i,k)
644 exp_od = exp(-opt_depth)
649 if (optics%answers_2018)
then 650 if (nsw*pen_sw_bnd(n,i)*sw_trans < min_sw_heat*min(1.0, i_habs*h(i,k)) ) sw_trans = 0.0
651 elseif ((nsw*pen_sw_bnd(n,i)*sw_trans < min_sw_heat) .and. (h(i,k) > h_min_heat))
then 652 if (nsw*pen_sw_bnd(n,i) <= min_sw_heat * (i_habs*(h(i,k) - h_min_heat)))
then 655 sw_trans = min(sw_trans, &
656 1.0 - (min_sw_heat*(i_habs*(h(i,k) - h_min_heat))) / (nsw*pen_sw_bnd(n,i)))
660 heat_bnd = pen_sw_bnd(n,i) * (1.0 - sw_trans)
661 if (adjustabsorptionprofile .and. (h_heat(i) > 0.0))
then 672 if (opt_depth > 1e-5)
then 673 swa = ((opt_depth + (opt_depth + 2.0)*exp_od) - 2.0) / &
674 ((opt_depth + opacity_band(n,i,k) * h_heat(i)) * &
679 swa = h(i,k) * (opt_depth * (1.0 - opt_depth)) / &
680 ((h_heat(i) + h(i,k)) * (6.0 - 3.0*opt_depth))
683 if (swa*(h_heat(i) + h(i,k)) > h_heat(i))
then 684 coswa_frac = (swa*(h_heat(i) + h(i,k)) - h_heat(i) ) / &
685 (swa*(h_heat(i) + h(i,k)))
686 swa = h_heat(i) / (h_heat(i) + h(i,k))
689 t_chg_above(i,k) = t_chg_above(i,k) + (swa * heat_bnd) / h_heat(i)
690 t(i,k) = t(i,k) + ((1.0 - swa) * heat_bnd) / h(i,k)
693 t(i,k) = t(i,k) + pen_sw_bnd(n,i) * (1.0 - sw_trans) / h(i,k)
697 if (opt_depth > 1e-2)
then 698 tke(i,k) = tke(i,k) - coswa_frac*heat_bnd*dsv_dt(i,k)* &
699 (0.5*h(i,k)*g_hconv2) * &
700 (opt_depth*(1.0+exp_od) - 2.0*(1.0-exp_od)) / (opt_depth*(1.0-exp_od))
704 tke(i,k) = tke(i,k) - coswa_frac*heat_bnd*dsv_dt(i,k)* &
705 (0.5*h(i,k)*g_hconv2) * &
706 (c1_6*opt_depth * (1.0 - c1_60*opt_depth**2))
710 pen_sw_bnd(n,i) = pen_sw_bnd(n,i) * sw_trans
716 if (h(i,k) >= 2.0*h_min_heat)
then 717 h_heat(i) = h_heat(i) + h(i,k)
718 elseif (h(i,k) > h_min_heat)
then 719 h_heat(i) = h_heat(i) + (2.0*h(i,k) - 2.0*h_min_heat)
726 do i=is,ie ; t_chg(i) = 0.0 ;
enddo 728 if (absorballsw)
then 733 pen_sw_rem(i) = pen_sw_bnd(1,i)
734 do n=2,nsw ; pen_sw_rem(i) = pen_sw_rem(i) + pen_sw_bnd(n,i) ;
enddo 736 do i=is,ie ;
if (pen_sw_rem(i) > 0.0) sw_remains = .true. ;
enddo 738 ih_limit = 1.0 / h_limit_fluxes
739 do i=is,ie ;
if ((pen_sw_rem(i) > 0.0) .and. (h_heat(i) > 0.0))
then 740 if (h_heat(i)*ih_limit >= 1.0)
then 741 t_chg(i) = pen_sw_rem(i) / h_heat(i) ; unabsorbed = 0.0
743 t_chg(i) = pen_sw_rem(i) * ih_limit
744 unabsorbed = 1.0 - h_heat(i)*ih_limit
746 do n=1,nsw ; pen_sw_bnd(n,i) = unabsorbed * pen_sw_bnd(n,i) ;
enddo 750 if (absorballsw .or. adjustabsorptionprofile)
then 751 do ks=nz,1,-1 ;
do i=is,ie
753 if (
present(ksort))
then 754 if (ksort(i,ks) <= 0) cycle
758 if (t_chg(i) > 0.0)
then 760 if (h(i,k) >= 2.0*h_min_heat)
then ; t(i,k) = t(i,k) + t_chg(i)
761 elseif (h(i,k) > h_min_heat)
then 762 t(i,k) = t(i,k) + t_chg(i) * (2.0 - 2.0*h_min_heat/h(i,k))
766 t_chg(i) = t_chg(i) + t_chg_above(i,k)
768 if (
present(htot) .and.
present(ttot))
then 769 do i=is,ie ; ttot(i) = ttot(i) + t_chg(i) * htot(i) ;
enddo 773 end subroutine absorbremainingsw
779 subroutine sumswoverbands(G, GV, US, h, nsw, optics, j, dt, &
780 H_limit_fluxes, absorbAllSW, iPen_SW_bnd, netPen)
784 real,
dimension(SZI_(G),SZK_(GV)), &
786 integer,
intent(in) :: nsw
790 integer,
intent(in) :: j
791 real,
intent(in) :: dt
792 real,
intent(in) :: h_limit_fluxes
795 logical,
intent(in) :: absorballsw
797 real,
dimension(max(nsw,1),SZI_(G)),
intent(in) :: ipen_sw_bnd
801 real,
dimension(SZI_(G),SZK_(GV)+1), &
802 intent(inout) :: netpen
806 real :: h_heat(szi_(g))
808 real :: pen_sw_rem(szi_(g))
813 real,
dimension(max(nsw,1),SZI_(G)) :: pen_sw_bnd
828 logical :: sw_remains
831 integer :: is, ie, nz, i, k, ks, n
834 min_sw_heat = optics%PenSW_flux_absorb*dt
835 i_habs = 1e3*gv%H_to_m
837 h_min_heat = 2.0*gv%Angstrom_H + gv%H_subroundoff
838 is = g%isc ; ie = g%iec ; nz = g%ke
840 pen_sw_bnd(:,:) = ipen_sw_bnd(:,:)
841 do i=is,ie ; h_heat(i) = 0.0 ;
enddo 845 netpen(i,1) = netpen(i,1) + pen_sw_bnd(n,i)
856 if (h(i,k) > 0.0)
then 857 do n=1,nsw ;
if (pen_sw_bnd(n,i) > 0.0)
then 859 opt_depth = h(i,k)*gv%H_to_m * optics%opacity_band(n,i,j,k)
860 exp_od = exp(-opt_depth)
865 if (optics%answers_2018)
then 866 if (nsw*pen_sw_bnd(n,i)*sw_trans < min_sw_heat*min(1.0, i_habs*h(i,k)) ) sw_trans = 0.0
867 elseif ((nsw*pen_sw_bnd(n,i)*sw_trans < min_sw_heat) .and. (h(i,k) > h_min_heat))
then 868 if (nsw*pen_sw_bnd(n,i) <= min_sw_heat * (i_habs*(h(i,k) - h_min_heat)))
then 871 sw_trans = min(sw_trans, &
872 1.0 - (min_sw_heat*(i_habs*(h(i,k) - h_min_heat))) / (nsw*pen_sw_bnd(n,i)))
876 pen_sw_bnd(n,i) = pen_sw_bnd(n,i) * sw_trans
877 netpen(i,k+1) = netpen(i,k+1) + pen_sw_bnd(n,i)
883 if (h(i,k) >= 2.0*h_min_heat)
then 884 h_heat(i) = h_heat(i) + h(i,k)
885 elseif (h(i,k) > h_min_heat)
then 886 h_heat(i) = h_heat(i) + (2.0*h(i,k) - 2.0*h_min_heat)
891 if (absorballsw)
then 897 pen_sw_rem(i) = pen_sw_bnd(1,i)
898 do n=2,nsw ; pen_sw_rem(i) = pen_sw_rem(i) + pen_sw_bnd(n,i) ;
enddo 900 do i=is,ie ;
if (pen_sw_rem(i) > 0.0) sw_remains = .true. ;
enddo 902 ih_limit = 1.0 / h_limit_fluxes
903 do i=is,ie ;
if ((pen_sw_rem(i) > 0.0) .and. (h_heat(i) > 0.0))
then 904 if (h_heat(i)*ih_limit < 1.0)
then 905 unabsorbed = 1.0 - h_heat(i)*ih_limit
909 do n=1,nsw ; pen_sw_bnd(n,i) = unabsorbed * pen_sw_bnd(n,i) ;
enddo 914 end subroutine sumswoverbands
919 subroutine opacity_init(Time, G, GV, US, param_file, diag, CS, optics)
920 type(time_type),
target,
intent(in) :: time
926 type(
diag_ctrl),
target,
intent(inout) :: diag
933 character(len=200) :: tmpstr
934 character(len=40) :: mdl =
"MOM_opacity" 935 character(len=40) :: bandnum, shortname
936 character(len=200) :: longname
937 character(len=40) :: scheme_string
939 # include "version_variable.h" 940 real :: pensw_absorb_minthick
942 real :: pensw_minthick_dflt
943 logical :: default_2018_answers
944 logical :: use_scheme
945 integer :: isd, ied, jsd, jed, nz, n
946 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
948 if (
associated(cs))
then 949 call mom_error(warning,
"opacity_init called with an associated"// &
950 "associated control structure.")
952 else ;
allocate(cs) ;
endif 960 call get_param(param_file, mdl,
"VAR_PEN_SW", cs%var_pen_sw, &
961 "If true, use one of the CHL_A schemes specified by "//&
962 "OPACITY_SCHEME to determine the e-folding depth of "//&
963 "incoming short wave radiation.", default=.false.)
965 cs%opacity_scheme = no_scheme ; scheme_string =
'' 966 if (cs%var_pen_sw)
then 967 call get_param(param_file, mdl,
"OPACITY_SCHEME", tmpstr, &
968 "This character string specifies how chlorophyll "//&
969 "concentrations are translated into opacities. Currently "//&
970 "valid options include:\n"//&
971 " \t\t MANIZZA_05 - Use Manizza et al., GRL, 2005. \n"//&
972 " \t\t MOREL_88 - Use Morel, JGR, 1988.", &
973 default=manizza_05_string)
974 if (len_trim(tmpstr) > 0)
then 975 tmpstr = uppercase(tmpstr)
977 case (manizza_05_string)
978 cs%opacity_scheme = manizza_05 ; scheme_string = manizza_05_string
979 case (morel_88_string)
980 cs%opacity_scheme = morel_88 ; scheme_string = morel_88_string
982 call mom_error(fatal,
"opacity_init: #DEFINE OPACITY_SCHEME "//&
983 trim(tmpstr) //
"in input file is invalid.")
985 call mom_mesg(
'opacity_init: opacity scheme set to "'//trim(tmpstr)//
'".', 5)
987 if (cs%opacity_scheme == no_scheme)
then 988 call mom_error(warning,
"opacity_init: No scheme has successfully "//&
989 "been specified for the opacity. Using the default MANIZZA_05.")
990 cs%opacity_scheme = manizza_05 ; scheme_string = manizza_05_string
993 call get_param(param_file, mdl,
"BLUE_FRAC_SW", cs%blue_frac, &
994 "The fraction of the penetrating shortwave radiation "//&
995 "that is in the blue band.", default=0.5, units=
"nondim")
997 call get_param(param_file, mdl,
"EXP_OPACITY_SCHEME", tmpstr, &
998 "This character string specifies which exponential "//&
999 "opacity scheme to utilize. Currently "//&
1000 "valid options include:\n"//&
1001 " \t\t SINGLE_EXP - Single Exponent decay. \n"//&
1002 " \t\t DOUBLE_EXP - Double Exponent decay.", &
1003 default=single_exp_string)
1004 if (len_trim(tmpstr) > 0)
then 1005 tmpstr = uppercase(tmpstr)
1006 select case (tmpstr)
1007 case (single_exp_string)
1008 cs%opacity_scheme = single_exp ; scheme_string = single_exp_string
1009 case (double_exp_string)
1010 cs%opacity_scheme = double_exp ; scheme_string = double_exp_string
1012 call mom_mesg(
'opacity_init: opacity scheme set to "'//trim(tmpstr)//
'".', 5)
1014 call get_param(param_file, mdl,
"PEN_SW_SCALE", cs%pen_sw_scale, &
1015 "The vertical absorption e-folding depth of the "//&
1016 "penetrating shortwave radiation.", units=
"m", default=0.0)
1018 if (cs%Opacity_scheme == double_exp )
then 1019 call get_param(param_file, mdl,
"PEN_SW_SCALE_2ND", cs%pen_sw_scale_2nd, &
1020 "The (2nd) vertical absorption e-folding depth of the "//&
1021 "penetrating shortwave radiation "//&
1022 "(use if SW_EXP_MODE==double.)",&
1023 units=
"m", default=0.0)
1024 call get_param(param_file, mdl,
"SW_1ST_EXP_RATIO", cs%sw_1st_exp_ratio, &
1025 "The fraction of 1st vertical absorption e-folding depth "//&
1026 "penetrating shortwave radiation if SW_EXP_MODE==double.",&
1027 units=
"m", default=0.0)
1028 elseif (cs%OPACITY_SCHEME == single_exp)
then 1030 cs%pen_sw_scale_2nd = 0.0
1031 cs%sw_1st_exp_ratio = 1.0
1033 call get_param(param_file, mdl,
"PEN_SW_FRAC", cs%pen_sw_frac, &
1034 "The fraction of the shortwave radiation that penetrates "//&
1035 "below the surface.", units=
"nondim", default=0.0)
1038 call get_param(param_file, mdl,
"PEN_SW_NBANDS", optics%nbands, &
1039 "The number of bands of penetrating shortwave radiation.", &
1041 if (cs%Opacity_scheme == double_exp )
then 1042 if (optics%nbands /= 2)
call mom_error(fatal, &
1043 "set_opacity: \Cannot use a double_exp opacity scheme with nbands!=2.")
1044 elseif (cs%Opacity_scheme == single_exp )
then 1045 if (optics%nbands /= 1)
call mom_error(fatal, &
1046 "set_opacity: \Cannot use a single_exp opacity scheme with nbands!=1.")
1049 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
1050 "This sets the default value for the various _2018_ANSWERS parameters.", &
1052 call get_param(param_file, mdl,
"OPTICS_2018_ANSWERS", optics%answers_2018, &
1053 "If true, use the order of arithmetic and expressions that recover the "//&
1054 "answers from the end of 2018. Otherwise, use updated expressions for "//&
1055 "handling the absorption of small remaining shortwave fluxes.", &
1056 default=default_2018_answers)
1058 call get_param(param_file, mdl,
"PEN_SW_FLUX_ABSORB", optics%PenSW_flux_absorb, &
1059 "A minimum remaining shortwave heating rate that will be simply absorbed in "//&
1060 "the next sufficiently thick layers for computational efficiency, instead of "//&
1061 "continuing to penetrate. The default, 2.5e-11 degC m s-1, is about 1e-4 W m-2 "//&
1062 "or 0.08 degC m century-1, but 0 is also a valid value.", &
1063 default=2.5e-11, units=
"degC m s-1", scale=gv%m_to_H*us%T_to_s)
1065 if (optics%answers_2018)
then ; pensw_minthick_dflt = 0.001 ;
else ; pensw_minthick_dflt = 1.0 ;
endif 1066 call get_param(param_file, mdl,
"PEN_SW_ABSORB_MINTHICK", pensw_absorb_minthick, &
1067 "A thickness that is used to absorb the remaining penetrating shortwave heat "//&
1068 "flux when it drops below PEN_SW_FLUX_ABSORB.", &
1069 default=pensw_minthick_dflt, units=
"m", scale=gv%m_to_H)
1070 optics%PenSW_absorb_Invlen = 1.0 / (pensw_absorb_minthick + gv%H_subroundoff)
1072 if (.not.
associated(optics%min_wavelength_band)) &
1073 allocate(optics%min_wavelength_band(optics%nbands))
1074 if (.not.
associated(optics%max_wavelength_band)) &
1075 allocate(optics%max_wavelength_band(optics%nbands))
1077 if (cs%opacity_scheme == manizza_05)
then 1078 optics%min_wavelength_band(1) =0
1079 optics%max_wavelength_band(1) =550
1080 if (optics%nbands >= 2)
then 1081 optics%min_wavelength_band(2)=550
1082 optics%max_wavelength_band(2)=700
1084 if (optics%nbands > 2)
then 1085 do n=3,optics%nbands
1086 optics%min_wavelength_band(n) =700
1087 optics%max_wavelength_band(n) =2800
1092 call get_param(param_file, mdl,
"OPACITY_LAND_VALUE", cs%opacity_land_value, &
1093 "The value to use for opacity over land. The default is "//&
1094 "10 m-1 - a value for muddy water.", units=
"m-1", default=10.0)
1096 if (.not.
associated(optics%opacity_band)) &
1097 allocate(optics%opacity_band(optics%nbands,isd:ied,jsd:jed,nz))
1098 if (.not.
associated(optics%sw_pen_band)) &
1099 allocate(optics%sw_pen_band(optics%nbands,isd:ied,jsd:jed))
1100 allocate(cs%id_opacity(optics%nbands)) ; cs%id_opacity(:) = -1
1102 cs%id_sw_pen = register_diag_field(
'ocean_model',
'SW_pen', diag%axesT1, time, &
1103 'Penetrating shortwave radiation flux into ocean',
'W m-2', conversion=us%QRZ_T_to_W_m2)
1104 cs%id_sw_vis_pen = register_diag_field(
'ocean_model',
'SW_vis_pen', diag%axesT1, time, &
1105 'Visible penetrating shortwave radiation flux into ocean',
'W m-2', conversion=us%QRZ_T_to_W_m2)
1106 do n=1,optics%nbands
1107 write(bandnum,
'(i3)') n
1108 shortname =
'opac_'//trim(adjustl(bandnum))
1109 longname =
'Opacity for shortwave radiation in band '//trim(adjustl(bandnum)) &
1110 //
', saved as L^-1 tanh(Opacity * L) for L = 10^-10 m' 1111 cs%id_opacity(n) = register_diag_field(
'ocean_model', shortname, diag%axesTL, time, &
1115 end subroutine opacity_init
1118 subroutine opacity_end(CS, optics)
1122 if (
associated(cs%id_opacity))
deallocate(cs%id_opacity)
1123 if (
associated(cs))
deallocate(cs)
1125 if (
present(optics))
then ;
if (
associated(optics))
then 1126 if (
associated(optics%opacity_band))
deallocate(optics%opacity_band)
1127 if (
associated(optics%sw_pen_band))
deallocate(optics%sw_pen_band)
1130 end subroutine opacity_end
Ocean grid type. See mom_grid for details.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Make a diagnostic available for averaging or output.
The control structure with paramters for the MOM_opacity module.
Describes various unit conversion factors.
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.
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.
Handy functions for manipulating strings.
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.