Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
25 changes: 16 additions & 9 deletions physics/GFS_PBL_generic_post.F90
Original file line number Diff line number Diff line change
Expand Up @@ -19,8 +19,8 @@ subroutine GFS_PBL_generic_post_run (im, levs, nvdiff, ntrac,
index_of_process_pbl, dqsfc_cpl, dusfci_cpl, dvsfci_cpl, dtsfci_cpl, dqsfci_cpl, dusfc_diag, dvsfc_diag, dtsfc_diag, &
dqsfc_diag, dusfci_diag, dvsfci_diag, dtsfci_diag, dqsfci_diag, &
rd, cp, fvirt, hvap, t1, q1, prsl, hflx, ushfsfci, oceanfrac, kdt, dusfc_cice, dvsfc_cice, &
dtsfc_cice, dqsfc_cice, wet, dry, icy, wind, stress_wat, hflx_wat, evap_wat, ugrs1, vgrs1, hffac, &
ugrs, vgrs, tgrs, qgrs, save_u, save_v, save_t, save_q, huge, errmsg, errflg)
dtsfc_cice, dqsfc_cice, use_med_flux, dtsfc_med, dqsfc_med, dusfc_med, dvsfc_med, wet, dry, icy, wind, stress_wat, &
hflx_wat, evap_wat, ugrs1, vgrs1, hffac, ugrs, vgrs, tgrs, qgrs, save_u, save_v, save_t, save_q, huge, errmsg, errflg)

use machine, only : kind_phys
use GFS_PBL_generic_common, only : set_aerosol_tracer_index
Expand All @@ -36,7 +36,7 @@ subroutine GFS_PBL_generic_post_run (im, levs, nvdiff, ntrac,
integer, intent(in) :: imp_physics_zhao_carr, imp_physics_mg, imp_physics_fer_hires
integer, intent(in) :: imp_physics_nssl
logical, intent(in) :: nssl_ccn_on, nssl_hail_on
logical, intent(in) :: ltaerosol, cplflx, cplaqm, cplchm, lssav, ldiag3d, lsidea
logical, intent(in) :: ltaerosol, cplflx, cplaqm, cplchm, lssav, ldiag3d, lsidea, use_med_flux
logical, intent(in) :: hybedmf, do_shoc, satmedmf, shinhong, do_ysu

logical, intent(in) :: flag_for_pbl_generic_tend
Expand All @@ -48,7 +48,7 @@ subroutine GFS_PBL_generic_post_run (im, levs, nvdiff, ntrac,
real(kind=kind_phys), dimension(:), intent(in) :: t1, q1, hflx, oceanfrac
real(kind=kind_phys), dimension(:,:), intent(in) :: prsl
real(kind=kind_phys), dimension(:), intent(in) :: dusfc_cice, dvsfc_cice, dtsfc_cice, dqsfc_cice, &
wind, stress_wat, hflx_wat, evap_wat, ugrs1, vgrs1
dtsfc_med, dqsfc_med, dusfc_med, dvsfc_med, wind, stress_wat, hflx_wat, evap_wat, ugrs1, vgrs1

real(kind=kind_phys), dimension(:,:, :), intent(in) :: qgrs
real(kind=kind_phys), dimension(:,:), intent(in) :: ugrs, vgrs, tgrs
Expand Down Expand Up @@ -318,11 +318,18 @@ subroutine GFS_PBL_generic_post_run (im, levs, nvdiff, ntrac,
endif
dtsfci_cpl(i) = cp * rho * hflx_wat(i) ! sensible heat flux over open ocean
dqsfci_cpl(i) = hvap * rho * evap_wat(i) ! latent heat flux over open ocean
else ! use results from PBL scheme for 100% open ocean
dusfci_cpl(i) = dusfc1(i)
dvsfci_cpl(i) = dvsfc1(i)
dtsfci_cpl(i) = dtsfc1(i)*hffac(i)
dqsfci_cpl(i) = dqsfc1(i)
else ! 100% open ocean
if (use_med_flux .and. kdt > 1) then ! use results from CMEPS mediator
dusfci_cpl(i) = dusfc_med(i)
dvsfci_cpl(i) = dvsfc_med(i)
dtsfci_cpl(i) = dtsfc_med(i)
dqsfci_cpl(i) = dqsfc_med(i)
else ! use results from PBL scheme
dusfci_cpl(i) = dusfc1(i)
dvsfci_cpl(i) = dvsfc1(i)
dtsfci_cpl(i) = dtsfc1(i)*hffac(i)
dqsfci_cpl(i) = dqsfc1(i)
end if
endif
!
dusfc_cpl (i) = dusfc_cpl(i) + dusfci_cpl(i) * dtf
Expand Down
39 changes: 39 additions & 0 deletions physics/GFS_PBL_generic_post.meta
Original file line number Diff line number Diff line change
Expand Up @@ -766,6 +766,45 @@
type = real
kind = kind_phys
intent = in
[use_med_flux]
standard_name = do_mediator_atmosphere_ocean_fluxes
long_name = flag for using atmosphere-ocean fluxes from mediator
units = flag
dimensions = ()
type = logical
intent = in
[dqsfc_med]
standard_name = surface_upward_latent_heat_flux_over_ocean_from_mediator
long_name = sfc latent heat flux input over ocean for coupling
units = W m-2
dimensions = (horizontal_loop_extent)
type = real
kind = kind_phys
intent = in
[dtsfc_med]
standard_name = surface_upward_sensible_heat_flux_over_ocean_from_mediator
long_name = sfc sensible heat flux input over ocean for coupling
units = W m-2
dimensions = (horizontal_loop_extent)
type = real
kind = kind_phys
intent = in
[dusfc_med]
standard_name = surface_x_momentum_flux_over_ocean_from_mediator
long_name = sfc x momentum flux over ocean for coupling
units = Pa
dimensions = (horizontal_loop_extent)
type = real
kind = kind_phys
intent = in
[dvsfc_med]
standard_name = surface_y_momentum_flux_over_ocean_from_mediator
long_name = sfc y momentum flux over ocean for coupling
units = Pa
dimensions = (horizontal_loop_extent)
type = real
kind = kind_phys
intent = in
[wet]
standard_name = flag_nonzero_wet_surface_fraction
long_name = flag indicating presence of some ocean or lake surface area fraction
Expand Down
7 changes: 7 additions & 0 deletions physics/GFS_debug.F90
Original file line number Diff line number Diff line change
Expand Up @@ -870,6 +870,13 @@ subroutine GFS_diagtoscreen_run (Model, Statein, Stateout, Sfcprop, Coupling,
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%v10mi_cpl ', Coupling%v10mi_cpl )
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%tsfci_cpl ', Coupling%tsfci_cpl )
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%psurfi_cpl ', Coupling%psurfi_cpl )
if (Model%use_med_flux) then
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%dusfcino_cpl ', Coupling%dusfcino_cpl )
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%dvsfcino_cpl ', Coupling%dvsfcino_cpl )
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%dtsfcino_cpl ', Coupling%dtsfcino_cpl )
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%dqsfcino_cpl ', Coupling%dqsfcino_cpl )
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%ulwsfcino_cpl', Coupling%ulwsfcino_cpl )
end if
Comment thread
uturuncoglu marked this conversation as resolved.
end if
if (Model%cplchm) then
call print_var(mpirank, omprank, blkno, Grid%xlat_d, Grid%xlon_d, 'Coupling%rainc_cpl', Coupling%rainc_cpl)
Expand Down
19 changes: 13 additions & 6 deletions physics/dcyc2t3.f
Original file line number Diff line number Diff line change
Expand Up @@ -180,9 +180,9 @@ subroutine dcyc2t3_run &
& sfcnirbmd,sfcnirdfd,sfcvisbmd,sfcvisdfd, &
& im, levs, deltim, fhswr, &
& dry, icy, wet, damp_LW_fluxadj, lfnc_k, lfnc_p0, &
& use_LW_jacobian, sfculw, fluxlwUP_jac, &
& t_lay, p_lay, p_lev, flux2D_lwUP, flux2D_lwDOWN, &
& pert_radtend, do_sppt,ca_global, tsfc_radtime, &
& use_LW_jacobian, sfculw, use_med_flux, sfculw_med, &
& fluxlwUP_jac, t_lay, p_lay, p_lev, flux2D_lwUP, &
& flux2D_lwDOWN,pert_radtend,do_sppt,ca_global,tsfc_radtime, &
! & dry, icy, wet, lprnt, ipr, &
! --- input/output:
& dtdt,dtdtnp,htrlw, &
Expand Down Expand Up @@ -213,14 +213,14 @@ subroutine dcyc2t3_run &
! logical lprnt
logical, dimension(:), intent(in) :: dry, icy, wet
logical, intent(in) :: use_LW_jacobian, damp_LW_fluxadj, &
& pert_radtend
& pert_radtend, use_med_flux
logical, intent(in) :: do_sppt,ca_global
real(kind=kind_phys), intent(in) :: solhr, slag, cdec, sdec, &
& deltim, fhswr, lfnc_k, lfnc_p0

real(kind=kind_phys), dimension(:), intent(in) :: &
& sinlat, coslat, xlon, coszen, tf, tsflw, sfcdlw, &
& sfcdsw, sfcnsw, sfculw, tsfc, tsfc_radtime
& sfcdsw, sfcnsw, sfculw, sfculw_med, tsfc, tsfc_radtime

real(kind=kind_phys), dimension(:), intent(in) :: &
& tsfc_lnd, tsfc_ice, tsfc_wat, &
Expand Down Expand Up @@ -345,8 +345,15 @@ subroutine dcyc2t3_run &
endif
if (wet(i)) then
tem2 = tsfc_wat(i) * tsfc_wat(i)
adjsfculw_wat(i) = sfcemis_wat(i) * con_sbc * tem2 * tem2
adjsfculw_wat(i) = sfcemis_wat(i) * con_sbc *
& tem2 * tem2
& + (one - sfcemis_wat(i)) * adjsfcdlw(i)
!> - replace upward longwave flux provided by the mediator (zero over lakes)
if (use_med_flux) then
if (sfculw_med(i) > f_eps) then
adjsfculw_wat(i) = sfculw_med(i)
end if
end if
endif

! if (lprnt .and. i == ipr) write(0,*)' in dcyc3: dry==',dry(i)
Expand Down
15 changes: 15 additions & 0 deletions physics/dcyc2t3.meta
Original file line number Diff line number Diff line change
Expand Up @@ -384,6 +384,21 @@
type = real
kind = kind_phys
intent = in
[use_med_flux]
standard_name = do_mediator_atmosphere_ocean_fluxes
long_name = flag for using atmosphere-ocean fluxes from mediator
units = flag
dimensions = ()
type = logical
intent = in
[sfculw_med]
standard_name = surface_upwelling_longwave_flux_over_ocean_from_mediator
long_name = surface upwelling LW flux over ocean for coupling
units = W m-2
dimensions = (horizontal_loop_extent)
type = real
kind = kind_phys
intent = in
[fluxlwUP_jac]
standard_name = RRTMGP_jacobian_of_lw_flux_upward
long_name = RRTMGP Jacobian upward longwave flux profile
Expand Down
2 changes: 1 addition & 1 deletion physics/sfc_diff.f
Original file line number Diff line number Diff line change
Expand Up @@ -513,7 +513,7 @@ subroutine stability &

dtv = thv1 - tvs
adtv = max(abs(dtv),0.001_kp)
dtv = sign(1.,dtv) * adtv
dtv = sign(1.0_kp,dtv) * adtv

if(thsfc_loc) then ! Use local potential temperature
rb = max(-5000.0_kp, (grav+grav) * dtv * z1
Expand Down
54 changes: 40 additions & 14 deletions physics/sfc_ocean.F
Original file line number Diff line number Diff line change
Expand Up @@ -29,7 +29,7 @@ subroutine sfc_ocean_run &
& ( im, hvap, cp, rd, eps, epsm1, rvrdm1, ps, u1, v1, t1, q1, &
& tskin, cm, ch, lseaspray, fm, fm10, &
& prsl1, prslki, wet, use_flake, wind, &, ! --- inputs
& flag_iter, &
& flag_iter, use_med_flux, dqsfc_med, dtsfc_med, &
& qsurf, cmm, chh, gflux, evap, hflx, ep, & ! --- outputs
& errmsg, errflg &
& )
Expand All @@ -43,6 +43,7 @@ subroutine sfc_ocean_run &
! inputs: !
! ( im, ps, u1, v1, t1, q1, tskin, cm, ch, lseaspray, fm, fm10, !
! prsl1, prslki, wet, use_flake, wind, flag_iter, !
! use_med_flux, !
! outputs: !
! qsurf, cmm, chh, gflux, evap, hflx, ep ) !
! !
Expand All @@ -58,6 +59,8 @@ subroutine sfc_ocean_run &
! reformatted the code and added program documentation !
! sep 2009 -- s. moorthi removed rcl and made pa as pressure unit !
! and furthur reformatted the code !
! dec 2021 -- u. turuncoglu added support for receiving fluxes !
! from mediator !
! !
! !
! ==================== defination of variables ==================== !
Expand All @@ -80,6 +83,9 @@ subroutine sfc_ocean_run &
! wet - logical, =T if any ocean/lak, =F otherwise im !
! wind - real, wind speed (m/s) im !
! flag_iter- logical, im !
! use_med_flux - logical, =T to use fluxes coming from med 1 !
! dqsfc_med- real, latent heat flux im !
! dtsfc_med- real, sensible heat flux im !
! !
! outputs: !
! qsurf - real, specific humidity at sfc im !
Expand Down Expand Up @@ -112,6 +118,12 @@ subroutine sfc_ocean_run &
logical, intent(in) :: lseaspray
!
logical, dimension(:), intent(in) :: flag_iter, wet, use_flake
!
logical, intent(in) :: use_med_flux

! To receive fluxes from mediator
real (kind=kind_phys), dimension(:), intent(in) :: &
& dqsfc_med, dtsfc_med

! --- outputs:
real (kind=kind_phys), dimension(:), intent(inout) :: qsurf, &
Expand Down Expand Up @@ -160,27 +172,42 @@ subroutine sfc_ocean_run &
! rho is density, qss is sat. hum. at surface

if ( flag(i) ) then
q0(i) = max( q1(i), qmin )
rho(i) = prsl1(i) / (rd*t1(i)*(one + rvrdm1*q0(i)))
if (use_med_flux) then
q0(i) = max( q1(i), qmin )
rho(i) = prsl1(i) / (rd*t1(i)*(one + rvrdm1*q0(i)))

tem = ch(i) * wind(i)
cmm(i) = cm(i) * wind(i)
chh(i) = rho(i) * tem

qss = fpvs( tskin(i) )
qss = eps*qss / (ps(i) + epsm1*qss)
hflx(i) = dtsfc_med(i)
evap(i) = dqsfc_med(i)

qsurf(i) = q1(i) + dqsfc_med(i) / (hvap*chh(i))
gflux(i) = zero
else
q0(i) = max( q1(i), qmin )
rho(i) = prsl1(i) / (rd*t1(i)*(one + rvrdm1*q0(i)))

qss = fpvs( tskin(i) )
qss = eps*qss / (ps(i) + epsm1*qss)

! --- ... rcp = rho cp ch v

rch = rho(i) * cp * ch(i) * wind(i)
tem = ch(i) * wind(i)
cmm(i) = cm(i) * wind(i)
chh(i) = rho(i) * tem
rch = rho(i) * cp * ch(i) * wind(i)
tem = ch(i) * wind(i)
cmm(i) = cm(i) * wind(i)
chh(i) = rho(i) * tem

! --- ... sensible and latent heat flux over open water

hflx(i) = rch * (tskin(i) - t1(i) * prslki(i))
hflx(i) = rch * (tskin(i) - t1(i) * prslki(i))

evap(i) = elocp * rch * (qss - q0(i))
evap(i) = elocp * rch * (qss - q0(i))

qsurf(i) = qss
gflux(i) = zero
qsurf(i) = qss
gflux(i) = zero
endif
endif
enddo
!
Expand Down Expand Up @@ -220,7 +247,6 @@ subroutine sfc_ocean_run &
endif
enddo
!

return
!...................................
end subroutine sfc_ocean_run
Expand Down
23 changes: 23 additions & 0 deletions physics/sfc_ocean.meta
Original file line number Diff line number Diff line change
Expand Up @@ -194,6 +194,29 @@
dimensions = (horizontal_loop_extent)
type = logical
intent = in
[use_med_flux]
standard_name = do_mediator_atmosphere_ocean_fluxes
long_name = flag for using atmosphere-ocean fluxes from mediator
units = flag
dimensions = ()
type = logical
intent = in
[dqsfc_med]
standard_name = surface_upward_latent_heat_flux_over_ocean_from_mediator
long_name = sfc latent heat flux input over ocean for coupling
units = W m-2
dimensions = (horizontal_loop_extent)
type = real
kind = kind_phys
intent = in
[dtsfc_med]
standard_name = surface_upward_sensible_heat_flux_over_ocean_from_mediator
long_name = sfc sensible heat flux input over ocean for coupling
units = W m-2
dimensions = (horizontal_loop_extent)
type = real
kind = kind_phys
intent = in
[qsurf]
standard_name = surface_specific_humidity_over_water
long_name = surface air saturation specific humidity over water
Expand Down