Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Introduce namelist flag xr_cnvcld to control if suspended grid-mean convective cloud condensate should be included in cloud fraction and optical depth calculation in the GFS suite #184

Merged
merged 2 commits into from
Mar 27, 2024
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
6 changes: 3 additions & 3 deletions physics/Interstitials/UFS_SCM_NEPTUNE/GFS_rrtmg_pre.F90
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,7 @@ module GFS_rrtmg_pre
!>\section rrtmg_pre_gen General Algorithm
subroutine GFS_rrtmg_pre_run (im, levs, lm, lmk, lmp, n_var_lndp, lextop,&
ltp, imfdeepcnv, imfdeepcnv_gf, imfdeepcnv_c3, me, ncnd, ntrac, &
num_p3d, npdf3d, &
num_p3d, npdf3d, xr_cnvcld, &
ncnvcld3d,ntqv, ntcw,ntiw, ntlnc, ntinc, ntrnc, ntsnc, ntccn, top_at_1,&
ntrw, ntsw, ntgl, nthl, ntwa, ntoz, ntsmoke, ntdust, ntcoarsepm, &
ntclamt, nleffr, nieffr, nseffr, lndp_type, kdt, &
Expand Down Expand Up @@ -129,7 +129,7 @@ subroutine GFS_rrtmg_pre_run (im, levs, lm, lmk, lmp, n_var_lndp, lextop,&
uni_cld, effr_in, do_mynnedmf, &
lmfshal, lmfdeep2, pert_clds, lcrick,&
lcnorm, top_at_1, lextop, mraerosol
logical, intent(in) :: rrfs_sd, aero_dir_fdb
logical, intent(in) :: rrfs_sd, aero_dir_fdb, xr_cnvcld

logical, intent(in) :: nssl_ccn_on, nssl_invertccn
integer, intent(in) :: spp_rad
Expand Down Expand Up @@ -981,7 +981,7 @@ subroutine GFS_rrtmg_pre_run (im, levs, lm, lmk, lmp, n_var_lndp, lextop,&
& iovr_dcorr, iovr_exp, iovr_exprand, idcor, idcor_con, &
& idcor_hogan, idcor_oreopoulos, lcrick, lcnorm, &
& imfdeepcnv, imfdeepcnv_gf, imfdeepcnv_c3, do_mynnedmf, &
& lgfdlmprad, &
& lgfdlmprad, xr_cnvcld, &
& uni_cld, lmfshal, lmfdeep2, cldcov, clouds1, &
& effrl, effri, effrr, effrs, effr_in, &
& effrl_inout, effri_inout, effrs_inout, &
Expand Down
7 changes: 7 additions & 0 deletions physics/Interstitials/UFS_SCM_NEPTUNE/GFS_rrtmg_pre.meta
Original file line number Diff line number Diff line change
Expand Up @@ -56,6 +56,13 @@
dimensions = ()
type = logical
intent = in
[xr_cnvcld]
standard_name = flag_for_suspended_convective_clouds_in_Xu_Randall
long_name = flag for using suspended convective clouds in Xu Randall
units = flag
dimensions = ()
type = logical
intent = in
[ltp]
standard_name = extra_top_layer
long_name = extra top layer for radiation
Expand Down
38 changes: 30 additions & 8 deletions physics/Radiation/radiation_clouds.f
Original file line number Diff line number Diff line change
Expand Up @@ -348,7 +348,7 @@ subroutine radiation_clouds_prop &
& iovr_dcorr, iovr_exp, iovr_exprand, idcor, idcor_con, &
& idcor_hogan, idcor_oreopoulos, lcrick, lcnorm, &
& imfdeepcnv, imfdeepcnv_gf, imfdeepcnv_c3, &
& do_mynnedmf, lgfdlmprad, &
& do_mynnedmf, lgfdlmprad, xr_cnvcld, &
& uni_cld, lmfshal, lmfdeep2, cldcov, clouds1, &
& effrl, effri, effrr, effrs, effr_in, &
& effrl_inout, effri_inout, effrs_inout, &
Expand Down Expand Up @@ -538,7 +538,8 @@ subroutine radiation_clouds_prop &


logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, effr_in, &
& do_mynnedmf, lgfdlmprad, top_at_1, lcrick, lcnorm
& do_mynnedmf, lgfdlmprad, top_at_1, lcrick, lcnorm, &
& xr_cnvcld

real (kind=kind_phys), dimension(:,:,:), intent(in) :: ccnd, &
& tracer1
Expand Down Expand Up @@ -727,7 +728,7 @@ subroutine radiation_clouds_prop &
call progcld_thompson_wsm6 (plyr,plvl,tlyr,qlyr,qstl, & ! --- inputs
& rhly,tracer1,xlat,xlon,slmsk,dz,delp, &
& ntrac-1, ntcw-1,ntiw-1,ntrw-1, &
& ntsw-1,ntgl-1,con_ttp, &
& ntsw-1,ntgl-1,con_ttp,xr_cnvcld, &
& IX, NLAY, NLP1, uni_cld, lmfshal, lmfdeep2, &
& cldcov(:,1:NLAY), cnvw, effrl_inout, &
& effri_inout, effrs_inout, &
Expand Down Expand Up @@ -801,7 +802,7 @@ subroutine radiation_clouds_prop &
call progcld_thompson_wsm6 (plyr,plvl,tlyr,qlyr,qstl, & ! --- inputs
& rhly,tracer1,xlat,xlon,slmsk,dz,delp, &
& ntrac-1, ntcw-1,ntiw-1,ntrw-1, &
& ntsw-1,ntgl-1,con_ttp, &
& ntsw-1,ntgl-1,con_ttp,xr_cnvcld, &
& IX, NLAY, NLP1, uni_cld, lmfshal, lmfdeep2, &
& cldcov(:,1:NLAY), cnvw, effrl, effri, effrs, &
& lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
Expand Down Expand Up @@ -1964,7 +1965,7 @@ subroutine progcld_thompson_wsm6 &
& ( plyr,plvl,tlyr,qlyr,qstl,rhly,clw, & ! --- inputs:
& xlat,xlon,slmsk,dz,delp, &
& ntrac,ntcw,ntiw,ntrw,ntsw,ntgl,con_ttp, &
& IX, NLAY, NLP1, &
& xr_cnvcld, IX, NLAY, NLP1, &
& uni_cld, lmfshal, lmfdeep2, cldcov, cnvw, &
& re_cloud,re_ice,re_snow, &
& lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
Expand Down Expand Up @@ -2051,7 +2052,8 @@ subroutine progcld_thompson_wsm6 &
integer, intent(in) :: IX, NLAY, NLP1
integer, intent(in) :: ntrac, ntcw, ntiw, ntrw, ntsw, ntgl

logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, lcnorm
logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, lcnorm, &
& xr_cnvcld

real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
& tlyr, qlyr, qstl, rhly, cldcov, delp, dz, dzlay, &
Expand Down Expand Up @@ -2122,23 +2124,43 @@ subroutine progcld_thompson_wsm6 &
! enddo
! endif

!> - Include grid-mean suspended cloud condensate in Xu-Randall cloud fraction
!> if xr_cnvcld is true:

if(xr_cnvcld)then
do k = 1, NLAY
do i = 1, IX
clwf(i,k) = clw(i,k,ntcw) + clw(i,k,ntiw) + clw(i,k,ntsw)
& + clw(i,k,ntrw) + cnvw(i,k)
enddo
enddo
else
do k = 1, NLAY
do i = 1, IX
clwf(i,k) = clw(i,k,ntcw) + clw(i,k,ntiw) + clw(i,k,ntsw)
& + clw(i,k,ntrw)
enddo
enddo
endif

!> - Compute total-cloud liquid/ice condensate path in \f$ g/m^2 \f$.
!> The total condensate includes convective condensate.
do k = 1, NLAY-1
do i = 1, IX
tem1 = cnvw(i,k)*(1.-tem2d(i,k))
if(xr_cnvcld)then
tem1 = cnvw(i,k)*(1.-tem2d(i,k))
else
tem1 = 0.
endif
cwp(i,k) = max(0.0, (clw(i,k,ntcw)+tem1) *
& gfac * delp(i,k))
if(tem1 > 1.e-12 .and. clw(i,k,ntcw) < 1.e-12)
& rew(i,k)=reliq_def
tem2 = cnvw(i,k)*tem2d(i,k)
if(xr_cnvcld)then
tem2 = cnvw(i,k)*tem2d(i,k)
else
tem2 = 0.
endif
cip(i,k) = max(0.0, (clw(i,k,ntiw) +
& snow2ice*clw(i,k,ntsw) + tem2) *
& gfac * delp(i,k))
Expand Down
Loading