21 use mom_ale,
only : ts_plm_edge_values, ts_ppm_edge_values,
ale_cs
23 implicit none ;
private
25 #include <MOM_memory.h>
27 public pressureforce_fv_init, pressureforce_fv_end
28 public pressureforce_fv_bouss, pressureforce_fv_nonbouss
42 type(time_type),
pointer :: time
45 logical :: usemasswghtinterp
46 logical :: boundary_extrap
49 logical :: reconstruct
54 integer :: recon_scheme
57 real :: stanley_t2_det_coeff
60 integer :: id_e_tidal = -1
61 integer :: id_tvar_sgs = -1
76 subroutine pressureforce_fv_nonbouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm, pbce, eta)
80 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
82 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(out) :: pfu
83 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(out) :: pfv
85 type(
ale_cs),
pointer :: ale_csp
86 real,
dimension(:,:),
optional,
pointer :: p_atm
88 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
optional,
intent(out) :: pbce
91 real,
dimension(SZI_(G),SZJ_(G)),
optional,
intent(out) :: eta
95 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: p
96 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
target :: &
101 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
106 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
111 real,
dimension(SZI_(G),SZJ_(G)) :: &
121 real,
dimension(SZI_(G)) :: rho_cv_bl
123 real,
dimension(SZIB_(G),SZJ_(G)) :: &
126 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)) :: &
128 real,
dimension(SZI_(G),SZJB_(G)) :: &
131 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)) :: &
133 real :: p_ref(szi_(g))
148 real :: rho_in_situ(szi_(g))
152 real,
parameter :: c1_6 = 1.0/6.0
153 integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz, nkmb
154 integer,
dimension(2) :: eosdom
157 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
158 nkmb=gv%nk_rho_varies
159 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
160 eosdom(1) = isq - (g%isd-1) ; eosdom(2) = g%iec+1 - (g%isd-1)
162 if (.not.
associated(cs))
call mom_error(fatal, &
163 "MOM_PressureForce_FV_nonBouss: Module must be initialized before it is used.")
164 if (cs%Stanley_T2_det_coeff>=0.)
call mom_error(fatal, &
165 "MOM_PressureForce_FV_nonBouss: The Stanley parameterization is not yet"//&
166 "implemented in non-Boussinesq mode.")
169 if (
present(p_atm))
then ;
if (
associated(p_atm)) use_p_atm = .true. ;
endif
170 use_eos =
associated(tv%eqn_of_state)
172 if (
associated(ale_csp)) use_ale = cs%reconstruct .and. use_eos
174 h_to_rl2_t2 = gv%g_Earth*gv%H_to_RZ
175 dp_neglect = gv%g_Earth*gv%H_to_RZ * gv%H_subroundoff
176 alpha_ref = 1.0 / cs%Rho0
177 i_gearth = 1.0 / gv%g_Earth
181 do j=jsq,jeq+1 ;
do i=isq,ieq+1
182 p(i,j,1) = p_atm(i,j)
186 do j=jsq,jeq+1 ;
do i=isq,ieq+1
191 do j=jsq,jeq+1 ;
do k=2,nz+1 ;
do i=isq,ieq+1
192 p(i,j,k) = p(i,j,k-1) + h_to_rl2_t2 * h(i,j,k-1)
193 enddo ;
enddo ;
enddo
201 tv_tmp%T => t_tmp ; tv_tmp%S => s_tmp
202 tv_tmp%eqn_of_state => tv%eqn_of_state
203 do i=isq,ieq+1 ; p_ref(i) = tv%P_Ref ;
enddo
206 do k=1,nkmb ;
do i=isq,ieq+1
207 tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
210 tv%eqn_of_state, eosdom)
211 do k=nkmb+1,nz ;
do i=isq,ieq+1
212 if (gv%Rlay(k) < rho_cv_bl(i))
then
213 tv_tmp%T(i,j,k) = tv%T(i,j,nkmb) ; tv_tmp%S(i,j,k) = tv%S(i,j,nkmb)
215 tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
220 tv_tmp%T => tv%T ; tv_tmp%S => tv%S
221 tv_tmp%eqn_of_state => tv%eqn_of_state
230 if ( cs%Recon_Scheme == 1 )
then
231 call ts_plm_edge_values(ale_csp, s_t, s_b, t_t, t_b, g, gv, tv, h, cs%boundary_extrap)
232 elseif ( cs%Recon_Scheme == 2)
then
233 call ts_ppm_edge_values(ale_csp, s_t, s_b, t_t, t_b, g, gv, tv, h, cs%boundary_extrap)
243 if ( cs%Recon_Scheme == 1 )
then
244 call int_spec_vol_dp_generic_plm( t_t(:,:,k), t_b(:,:,k), s_t(:,:,k), s_b(:,:,k), &
245 p(:,:,k), p(:,:,k+1), alpha_ref, dp_neglect, p(:,:,nz+1), g%HI, &
246 tv%eqn_of_state, us, dza(:,:,k), intp_dza(:,:,k), intx_dza(:,:,k), inty_dza(:,:,k), &
247 usemasswghtinterp=cs%useMassWghtInterp)
248 elseif ( cs%Recon_Scheme == 2 )
then
249 call mom_error(fatal,
"PressureForce_FV_nonBouss: "//&
250 "int_spec_vol_dp_generic_ppm does not exist yet.")
257 call int_specific_vol_dp(tv_tmp%T(:,:,k), tv_tmp%S(:,:,k), p(:,:,k), &
258 p(:,:,k+1), alpha_ref, g%HI, tv%eqn_of_state, &
259 us, dza(:,:,k), intp_dza(:,:,k), intx_dza(:,:,k), &
260 inty_dza(:,:,k), bathyp=p(:,:,nz+1), dp_tiny=dp_neglect, &
261 usemasswghtinterp=cs%useMassWghtInterp)
264 alpha_anom = 1.0 / gv%Rlay(k) - alpha_ref
265 do j=jsq,jeq+1 ;
do i=isq,ieq+1
266 dp(i,j) = h_to_rl2_t2 * h(i,j,k)
267 dza(i,j,k) = alpha_anom * dp(i,j)
268 intp_dza(i,j,k) = 0.5 * alpha_anom * dp(i,j)**2
270 do j=js,je ;
do i=isq,ieq
271 intx_dza(i,j,k) = 0.5 * alpha_anom * (dp(i,j)+dp(i+1,j))
273 do j=jsq,jeq ;
do i=is,ie
274 inty_dza(i,j,k) = 0.5 * alpha_anom * (dp(i,j)+dp(i,j+1))
290 za(i,j) = alpha_ref*p(i,j,nz+1) - gv%g_Earth*g%bathyT(i,j)
292 do k=nz,1,-1 ;
do i=isq,ieq+1
293 za(i,j) = za(i,j) + dza(i,j,k)
300 do j=jsq,jeq+1 ;
do i=isq,ieq+1
301 ssh(i,j) = (za(i,j) - alpha_ref*p(i,j,1)) * i_gearth
303 call calc_tidal_forcing(cs%Time, ssh, e_tidal, g, cs%tides_CSp, m_to_z=us%m_to_Z)
305 do j=jsq,jeq+1 ;
do i=isq,ieq+1
306 za(i,j) = za(i,j) - gv%g_Earth * e_tidal(i,j)
310 if (cs%GFS_scale < 1.0)
then
315 call calculate_density(tv_tmp%T(:,j,1), tv_tmp%S(:,j,1), p(:,j,1), rho_in_situ, &
316 tv%eqn_of_state, eosdom)
319 dm(i,j) = (cs%GFS_scale - 1.0) * (p(i,j,1)*(1.0/rho_in_situ(i) - alpha_ref) + za(i,j))
324 do j=jsq,jeq+1 ;
do i=isq,ieq+1
325 dm(i,j) = (cs%GFS_scale - 1.0) * (p(i,j,1)*(1.0/gv%Rlay(1) - alpha_ref) + za(i,j))
338 do j=js,je ;
do i=isq,ieq
339 intx_za(i,j) = 0.5*(za(i,j) + za(i+1,j))
342 do j=jsq,jeq ;
do i=is,ie
343 inty_za(i,j) = 0.5*(za(i,j) + za(i,j+1))
349 do j=jsq,jeq+1 ;
do i=isq,ieq+1
350 dp(i,j) = h_to_rl2_t2 * h(i,j,k)
351 za(i,j) = za(i,j) - dza(i,j,k)
354 do j=js,je ;
do i=isq,ieq
355 intx_za(i,j) = intx_za(i,j) - intx_dza(i,j,k)
356 pfu(i,j,k) = ( ((za(i,j)*dp(i,j) + intp_dza(i,j,k)) - &
357 (za(i+1,j)*dp(i+1,j) + intp_dza(i+1,j,k))) + &
358 ((dp(i+1,j) - dp(i,j)) * intx_za(i,j) - &
359 (p(i+1,j,k) - p(i,j,k)) * intx_dza(i,j,k)) ) * &
360 (2.0*g%IdxCu(i,j) / ((dp(i,j) + dp(i+1,j)) + dp_neglect))
363 do j=jsq,jeq ;
do i=is,ie
364 inty_za(i,j) = inty_za(i,j) - inty_dza(i,j,k)
365 pfv(i,j,k) = (((za(i,j)*dp(i,j) + intp_dza(i,j,k)) - &
366 (za(i,j+1)*dp(i,j+1) + intp_dza(i,j+1,k))) + &
367 ((dp(i,j+1) - dp(i,j)) * inty_za(i,j) - &
368 (p(i,j+1,k) - p(i,j,k)) * inty_dza(i,j,k))) * &
369 (2.0*g%IdyCv(i,j) / ((dp(i,j) + dp(i,j+1)) + dp_neglect))
372 if (cs%GFS_scale < 1.0)
then
375 do j=js,je ;
do i=isq,ieq
376 pfu(i,j,k) = pfu(i,j,k) - (dm(i+1,j) - dm(i,j)) * g%IdxCu(i,j)
379 do j=jsq,jeq ;
do i=is,ie
380 pfv(i,j,k) = pfv(i,j,k) - (dm(i,j+1) - dm(i,j)) * g%IdyCv(i,j)
385 if (
present(pbce))
then
386 call set_pbce_nonbouss(p, tv_tmp, g, gv, us, cs%GFS_scale, pbce)
389 if (
present(eta))
then
390 pa_to_h = 1.0 / (gv%g_Earth * gv%H_to_RZ)
393 do j=jsq,jeq+1 ;
do i=isq,ieq+1
394 eta(i,j) = (p(i,j,nz+1) - p_atm(i,j))*pa_to_h
398 do j=jsq,jeq+1 ;
do i=isq,ieq+1
399 eta(i,j) = p(i,j,nz+1)*pa_to_h
404 if (cs%id_e_tidal>0)
call post_data(cs%id_e_tidal, e_tidal, cs%diag)
406 end subroutine pressureforce_fv_nonbouss
416 subroutine pressureforce_fv_bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm, pbce, eta)
420 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
422 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(out) :: pfu
423 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(out) :: pfv
425 type(
ale_cs),
pointer :: ale_csp
426 real,
dimension(:,:),
optional,
pointer :: p_atm
428 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
optional,
intent(out) :: pbce
431 real,
dimension(SZI_(G),SZJ_(G)),
optional,
intent(out) :: eta
435 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: e
436 real,
dimension(SZI_(G),SZJ_(G)) :: &
441 real,
dimension(SZI_(G)) :: &
444 real,
dimension(SZI_(G),SZJ_(G)) :: &
452 real,
dimension(SZIB_(G),SZJ_(G)) :: &
456 real,
dimension(SZI_(G),SZJB_(G)) :: &
461 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
target :: &
466 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
469 real :: rho_in_situ(szi_(g))
470 real :: p_ref(szi_(g))
488 real,
parameter :: c1_6 = 1.0/6.0
489 integer,
dimension(2) :: eosdom
490 integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz, nkmb
493 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
494 nkmb=gv%nk_rho_varies
495 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
496 eosdom(1) = isq - (g%isd-1) ; eosdom(2) = g%iec+1 - (g%isd-1)
498 if (.not.
associated(cs))
call mom_error(fatal, &
499 "MOM_PressureForce_FV_Bouss: Module must be initialized before it is used.")
502 if (
present(p_atm))
then ;
if (
associated(p_atm)) use_p_atm = .true. ;
endif
503 use_eos =
associated(tv%eqn_of_state)
504 do i=isq,ieq+1 ; p0(i) = 0.0 ;
enddo
506 if (
associated(ale_csp)) use_ale = cs%reconstruct .and. use_eos
508 if (cs%Stanley_T2_det_coeff>=0.)
then
509 if (.not.
associated(tv%varT))
call safe_alloc_ptr(tv%varT, g%isd, g%ied, g%jsd, g%jed, gv%ke)
510 do k=1, nz ;
do j=js-1,je+1 ;
do i=is-1,ie+1
530 hl(1) = h(i,j,k) * g%mask2dT(i,j)
531 hl(2) = h(i-1,j,k) * g%mask2dCu(i-1,j)
532 hl(3) = h(i+1,j,k) * g%mask2dCu(i,j)
533 hl(4) = h(i,j-1,k) * g%mask2dCv(i,j-1)
534 hl(5) = h(i,j+1,k) * g%mask2dCv(i,j)
535 r_sm_h = 1. / ( ( hl(1) + ( ( hl(2) + hl(3) ) + ( hl(4) + hl(5) ) ) ) + gv%H_subroundoff )
537 tl(1) = tv%T(i,j,k) ; tl(2) = tv%T(i-1,j,k) ; tl(3) = tv%T(i+1,j,k)
538 tl(4) = tv%T(i,j-1,k) ; tl(5) = tv%T(i,j+1,k)
539 mn_t = ( hl(1)*tl(1) + ( ( hl(2)*tl(2) + hl(3)*tl(3) ) + ( hl(4)*tl(4) + hl(5)*tl(5) ) ) ) * r_sm_h
541 tl(:) = tl(:) - mn_t ; mn_t = 0.
543 mn_t2 = ( hl(1)*tl(1)*tl(1) + ( ( hl(2)*tl(2)*tl(2) + hl(3)*tl(3)*tl(3) ) &
544 + ( hl(4)*tl(4)*tl(4) + hl(5)*tl(5)*tl(5) ) ) ) * r_sm_h
547 tv%varT(i,j,k) = cs%Stanley_T2_det_coeff * max(0., mn_t2)
548 enddo ;
enddo ;
enddo
551 h_neglect = gv%H_subroundoff
552 dz_neglect = gv%H_subroundoff * gv%H_to_Z
553 i_rho0 = 1.0 / gv%Rho0
554 g_rho0 = gv%g_Earth / gv%Rho0
565 e(i,j,1) = -g%bathyT(i,j)
567 do k=1,nz ;
do i=isq,ieq+1
568 e(i,j,1) = e(i,j,1) + h(i,j,k)*gv%H_to_Z
571 call calc_tidal_forcing(cs%Time, e(:,:,1), e_tidal, g, cs%tides_CSp, m_to_z=us%m_to_Z)
577 do j=jsq,jeq+1 ;
do i=isq,ieq+1
578 e(i,j,nz+1) = -(g%bathyT(i,j) + e_tidal(i,j))
582 do j=jsq,jeq+1 ;
do i=isq,ieq+1
583 e(i,j,nz+1) = -g%bathyT(i,j)
587 do j=jsq,jeq+1;
do k=nz,1,-1 ;
do i=isq,ieq+1
588 e(i,j,k) = e(i,j,k+1) + h(i,j,k)*gv%H_to_Z
589 enddo ;
enddo ;
enddo
598 tv_tmp%T => t_tmp ; tv_tmp%S => s_tmp
599 tv_tmp%eqn_of_state => tv%eqn_of_state
601 do i=isq,ieq+1 ; p_ref(i) = tv%P_Ref ;
enddo
604 do k=1,nkmb ;
do i=isq,ieq+1
605 tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
608 tv%eqn_of_state, eosdom)
610 do k=nkmb+1,nz ;
do i=isq,ieq+1
611 if (gv%Rlay(k) < rho_cv_bl(i))
then
612 tv_tmp%T(i,j,k) = tv%T(i,j,nkmb) ; tv_tmp%S(i,j,k) = tv%S(i,j,nkmb)
614 tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
619 tv_tmp%T => tv%T ; tv_tmp%S => tv%S
620 tv_tmp%eqn_of_state => tv%eqn_of_state
624 if (cs%GFS_scale < 1.0)
then
630 call calculate_density(tv_tmp%T(:,j,1), tv_tmp%S(:,j,1), p_atm(:,j), rho_in_situ, &
631 tv%eqn_of_state, eosdom)
634 tv%eqn_of_state, eosdom)
637 dm(i,j) = (cs%GFS_scale - 1.0) * (g_rho0 * rho_in_situ(i)) * e(i,j,1)
642 do j=jsq,jeq+1 ;
do i=isq,ieq+1
643 dm(i,j) = (cs%GFS_scale - 1.0) * (g_rho0 * gv%Rlay(1)) * e(i,j,1)
654 if ( cs%Recon_Scheme == 1 )
then
655 call ts_plm_edge_values(ale_csp, s_t, s_b, t_t, t_b, g, gv, tv, h, cs%boundary_extrap)
656 elseif ( cs%Recon_Scheme == 2 )
then
657 call ts_ppm_edge_values(ale_csp, s_t, s_b, t_t, t_b, g, gv, tv, h, cs%boundary_extrap)
666 do j=jsq,jeq+1 ;
do i=isq,ieq+1
667 pa(i,j) = (rho_ref*gv%g_Earth)*e(i,j,1) + p_atm(i,j)
671 do j=jsq,jeq+1 ;
do i=isq,ieq+1
672 pa(i,j) = (rho_ref*gv%g_Earth)*e(i,j,1)
676 do j=js,je ;
do i=isq,ieq
677 intx_pa(i,j) = 0.5*(pa(i,j) + pa(i+1,j))
680 do j=jsq,jeq ;
do i=is,ie
681 inty_pa(i,j) = 0.5*(pa(i,j) + pa(i,j+1))
695 if ( cs%Recon_Scheme == 1 )
then
696 call int_density_dz_generic_plm(k, tv, t_t, t_b, s_t, s_b, e, &
697 rho_ref, cs%Rho0, gv%g_Earth, dz_neglect, g%bathyT, &
698 g%HI, gv, tv%eqn_of_state, us, dpa, intz_dpa, intx_dpa, inty_dpa, &
699 usemasswghtinterp=cs%useMassWghtInterp)
700 elseif ( cs%Recon_Scheme == 2 )
then
701 call int_density_dz_generic_ppm(k, tv, t_t, t_b, s_t, s_b, e, &
702 rho_ref, cs%Rho0, gv%g_Earth, dz_neglect, g%bathyT, &
703 g%HI, gv, tv%eqn_of_state, us, dpa, intz_dpa, intx_dpa, inty_dpa, &
704 usemasswghtinterp=cs%useMassWghtInterp)
707 call int_density_dz(tv_tmp%T(:,:,k), tv_tmp%S(:,:,k), e(:,:,k), e(:,:,k+1), &
708 rho_ref, cs%Rho0, gv%g_Earth, g%HI, tv%eqn_of_state, us, dpa, &
709 intz_dpa, intx_dpa, inty_dpa, g%bathyT, dz_neglect, cs%useMassWghtInterp)
712 do j=jsq,jeq+1 ;
do i=isq,ieq+1
713 intz_dpa(i,j) = intz_dpa(i,j)*gv%Z_to_H
717 do j=jsq,jeq+1 ;
do i=isq,ieq+1
718 dz_geo(i,j) = gv%g_Earth * gv%H_to_Z*h(i,j,k)
719 dpa(i,j) = (gv%Rlay(k) - rho_ref) * dz_geo(i,j)
720 intz_dpa(i,j) = 0.5*(gv%Rlay(k) - rho_ref) * dz_geo(i,j)*h(i,j,k)
723 do j=js,je ;
do i=isq,ieq
724 intx_dpa(i,j) = 0.5*(gv%Rlay(k) - rho_ref) * (dz_geo(i,j) + dz_geo(i+1,j))
727 do j=jsq,jeq ;
do i=is,ie
728 inty_dpa(i,j) = 0.5*(gv%Rlay(k) - rho_ref) * (dz_geo(i,j) + dz_geo(i,j+1))
734 do j=js,je ;
do i=isq,ieq
735 pfu(i,j,k) = (((pa(i,j)*h(i,j,k) + intz_dpa(i,j)) - &
736 (pa(i+1,j)*h(i+1,j,k) + intz_dpa(i+1,j))) + &
737 ((h(i+1,j,k) - h(i,j,k)) * intx_pa(i,j) - &
738 (e(i+1,j,k+1) - e(i,j,k+1)) * intx_dpa(i,j) * gv%Z_to_H)) * &
739 ((2.0*i_rho0*g%IdxCu(i,j)) / &
740 ((h(i,j,k) + h(i+1,j,k)) + h_neglect))
741 intx_pa(i,j) = intx_pa(i,j) + intx_dpa(i,j)
745 do j=jsq,jeq ;
do i=is,ie
746 pfv(i,j,k) = (((pa(i,j)*h(i,j,k) + intz_dpa(i,j)) - &
747 (pa(i,j+1)*h(i,j+1,k) + intz_dpa(i,j+1))) + &
748 ((h(i,j+1,k) - h(i,j,k)) * inty_pa(i,j) - &
749 (e(i,j+1,k+1) - e(i,j,k+1)) * inty_dpa(i,j) * gv%Z_to_H)) * &
750 ((2.0*i_rho0*g%IdyCv(i,j)) / &
751 ((h(i,j,k) + h(i,j+1,k)) + h_neglect))
752 inty_pa(i,j) = inty_pa(i,j) + inty_dpa(i,j)
755 do j=jsq,jeq+1 ;
do i=isq,ieq+1
756 pa(i,j) = pa(i,j) + dpa(i,j)
760 if (cs%GFS_scale < 1.0)
then
763 do j=js,je ;
do i=isq,ieq
764 pfu(i,j,k) = pfu(i,j,k) - (dm(i+1,j) - dm(i,j)) * g%IdxCu(i,j)
767 do j=jsq,jeq ;
do i=is,ie
768 pfv(i,j,k) = pfv(i,j,k) - (dm(i,j+1) - dm(i,j)) * g%IdyCv(i,j)
773 if (
present(pbce))
then
774 call set_pbce_bouss(e, tv_tmp, g, gv, us, cs%Rho0, cs%GFS_scale, pbce)
777 if (
present(eta))
then
782 do j=jsq,jeq+1 ;
do i=isq,ieq+1
783 eta(i,j) = e(i,j,1)*gv%Z_to_H + e_tidal(i,j)*gv%Z_to_H
787 do j=jsq,jeq+1 ;
do i=isq,ieq+1
788 eta(i,j) = e(i,j,1)*gv%Z_to_H
793 if (cs%id_e_tidal>0)
call post_data(cs%id_e_tidal, e_tidal, cs%diag)
794 if (cs%id_tvar_sgs>0)
call post_data(cs%id_tvar_sgs, tv%varT, cs%diag)
796 end subroutine pressureforce_fv_bouss
799 subroutine pressureforce_fv_init(Time, G, GV, US, param_file, diag, CS, tides_CSp)
800 type(time_type),
target,
intent(in) :: time
805 type(
diag_ctrl),
target,
intent(inout) :: diag
809 # include "version_variable.h"
810 character(len=40) :: mdl
813 if (
associated(cs))
then
814 call mom_error(warning,
"PressureForce_init called with an associated "// &
815 "control structure.")
817 else ;
allocate(cs) ;
endif
819 cs%diag => diag ; cs%Time => time
820 if (
present(tides_csp))
then
821 if (
associated(tides_csp)) cs%tides_CSp => tides_csp
824 mdl =
"MOM_PressureForce_FV"
826 call get_param(param_file, mdl,
"RHO_0", cs%Rho0, &
827 "The mean ocean density used with BOUSSINESQ true to "//&
828 "calculate accelerations and the mass for conservation "//&
829 "properties, or with BOUSSINSEQ false to convert some "//&
830 "parameters from vertical units of m to kg m-2.", &
831 units=
"kg m-3", default=1035.0, scale=us%kg_m3_to_R)
832 call get_param(param_file, mdl,
"TIDES", cs%tides, &
833 "If true, apply tidal momentum forcing.", default=.false.)
834 call get_param(param_file,
"MOM",
"USE_REGRIDDING", use_ale, &
835 "If True, use the ALE algorithm (regridding/remapping). "//&
836 "If False, use the layered isopycnal algorithm.", default=.false. )
837 call get_param(param_file, mdl,
"MASS_WEIGHT_IN_PRESSURE_GRADIENT", cs%useMassWghtInterp, &
838 "If true, use mass weighting when interpolating T/S for "//&
839 "integrals near the bathymetry in FV pressure gradient "//&
840 "calculations.", default=.false.)
841 call get_param(param_file, mdl,
"RECONSTRUCT_FOR_PRESSURE", cs%reconstruct, &
842 "If True, use vertical reconstruction of T & S within "//&
843 "the integrals of the FV pressure gradient calculation. "//&
844 "If False, use the constant-by-layer algorithm. "//&
845 "The default is set by USE_REGRIDDING.", &
847 call get_param(param_file, mdl,
"PRESSURE_RECONSTRUCTION_SCHEME", cs%Recon_Scheme, &
848 "Order of vertical reconstruction of T/S to use in the "//&
849 "integrals within the FV pressure gradient calculation.\n"//&
850 " 0: PCM or no reconstruction.\n"//&
851 " 1: PLM reconstruction.\n"//&
852 " 2: PPM reconstruction.", default=1)
853 call get_param(param_file, mdl,
"BOUNDARY_EXTRAPOLATION_PRESSURE", cs%boundary_extrap, &
854 "If true, the reconstruction of T & S for pressure in "//&
855 "boundary cells is extrapolated, rather than using PCM "//&
856 "in these cells. If true, the same order polynomial is "//&
857 "used as is used for the interior cells.", default=.true.)
858 call get_param(param_file, mdl,
"PGF_STANLEY_T2_DET_COEFF", cs%Stanley_T2_det_coeff, &
859 "The coefficient correlating SGS temperature variance with "// &
860 "the mean temperature gradient in the deterministic part of "// &
861 "the Stanley form of the Brankart correction. "// &
862 "Negative values disable the scheme.", units=
"nondim", default=-1.0)
863 if (cs%Stanley_T2_det_coeff>=0.)
then
864 cs%id_tvar_sgs = register_diag_field(
'ocean_model',
'tvar_sgs_pgf', diag%axesTL, &
865 time,
'SGS temperature variance used in PGF',
'degC2')
868 cs%id_e_tidal = register_diag_field(
'ocean_model',
'e_tidal', diag%axesT1, &
869 time,
'Tidal Forcing Astronomical and SAL Height Anomaly',
'meter', conversion=us%Z_to_m)
873 if (gv%g_prime(1) /= gv%g_Earth) cs%GFS_scale = gv%g_prime(1) / gv%g_Earth
875 call log_param(param_file, mdl,
"GFS / G_EARTH", cs%GFS_scale)
877 end subroutine pressureforce_fv_init
880 subroutine pressureforce_fv_end(CS)
883 if (
associated(cs))
deallocate(cs)
884 end subroutine pressureforce_fv_end