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

Add variables we can use to modify the radiative fluxes to the Statein #158

Merged
merged 13 commits into from
Mar 12, 2021
Merged
149 changes: 141 additions & 8 deletions FV3/gfsphysics/GFS_layer/GFS_diagnostics.F90
Original file line number Diff line number Diff line change
Expand Up @@ -200,13 +200,21 @@ subroutine GFS_externaldiag_populate (ExtDiag, Model, Statein, Stateout, Sfcprop
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%cnvfac = cn_one
ExtDiag(idx)%time_avg = .TRUE.
ExtDiag(idx)%time_avg_kind = 'rad_sw'
if (.not. Model%override_surface_radiative_fluxes) then
ExtDiag(idx)%time_avg_kind = 'rad_sw'
endif
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%fluxr(:,4)
enddo
if (Model%override_surface_radiative_fluxes) then
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%dswsfc(:)
enddo
else
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%fluxr(:,4)
enddo
endif

idx = idx + 1
ExtDiag(idx)%axes = 2
Expand All @@ -229,13 +237,21 @@ subroutine GFS_externaldiag_populate (ExtDiag, Model, Statein, Stateout, Sfcprop
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%cnvfac = cn_one
ExtDiag(idx)%time_avg = .TRUE.
ExtDiag(idx)%time_avg_kind = 'rad_sw'
if (.not. Model%override_surface_radiative_fluxes) then
ExtDiag(idx)%time_avg_kind = 'rad_sw'
endif
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%fluxr(:,3)
enddo
if (Model%override_surface_radiative_fluxes) then
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%uswsfc(:)
enddo
else
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%fluxr(:,3)
spencerkclark marked this conversation as resolved.
Show resolved Hide resolved
enddo
endif

idx = idx + 1
ExtDiag(idx)%axes = 2
Expand All @@ -250,6 +266,123 @@ subroutine GFS_externaldiag_populate (ExtDiag, Model, Statein, Stateout, Sfcprop
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%uswsfci(:)
enddo

if (Model%override_surface_radiative_fluxes) then
idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'DLWRF_from_rrtmg'
ExtDiag(idx)%desc = 'surface downward longwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%cnvfac = cn_one
ExtDiag(idx)%time_avg = .TRUE.
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%dlwsfc_rrtmg(:)
enddo

idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'DLWRFI_from_rrtmg'
ExtDiag(idx)%desc = 'instantaneous surface downward longwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%dlwsfci_rrtmg(:)
enddo


idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'ULWRF_from_rrtmg'
ExtDiag(idx)%desc = 'surface upward longwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%cnvfac = cn_one
ExtDiag(idx)%time_avg = .TRUE.
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%ulwsfc_rrtmg(:)
enddo

idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'ULWRFI_from_rrtmg'
ExtDiag(idx)%desc = 'instantaneous surface upward longwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%ulwsfci_rrtmg(:)
enddo
spencerkclark marked this conversation as resolved.
Show resolved Hide resolved

idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'DSWRF_from_rrtmg'
ExtDiag(idx)%desc = 'averaged surface downward shortwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%cnvfac = cn_one
ExtDiag(idx)%time_avg = .TRUE.
ExtDiag(idx)%time_avg_kind = 'rad_sw'
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%fluxr(:,4)
enddo

idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'DSWRFI_from_rrtmg'
ExtDiag(idx)%desc = 'instantaneous surface downward shortwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%dswsfci_rrtmg(:)
enddo

idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'USWRF_from_rrtmg'
ExtDiag(idx)%desc = 'averaged surface upward shortwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%cnvfac = cn_one
ExtDiag(idx)%time_avg = .TRUE.
ExtDiag(idx)%time_avg_kind = 'rad_sw'
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%fluxr(:,3)
enddo

idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'USWRFI_from_rrtmg'
ExtDiag(idx)%desc = 'instantaneous surface upward shortwave flux due to RRTMG'
ExtDiag(idx)%unit = 'W/m**2'
ExtDiag(idx)%mod_name = 'gfs_phys'
ExtDiag(idx)%intpl_method = 'bilinear'
ExtDiag(idx)%coarse_graining_method = 'area_weighted'
allocate (ExtDiag(idx)%data(nblks))
do nb = 1,nblks
ExtDiag(idx)%data(nb)%var2 => IntDiag(nb)%uswsfci_rrtmg(:)
enddo
endif

idx = idx + 1
ExtDiag(idx)%axes = 2
ExtDiag(idx)%name = 'duvb_ave'
Expand Down
56 changes: 43 additions & 13 deletions FV3/gfsphysics/GFS_layer/GFS_physics_driver.F90
Original file line number Diff line number Diff line change
Expand Up @@ -529,7 +529,7 @@ subroutine GFS_physics_driver &
stress, t850, ep1d, gamt, gamq, sigmaf, &
wind, work1, work2, work3, work4, runof, xmu, fm10, fh2, &
tx1, tx2, tx3, tx4, ctei_r, evbs, evcw, trans, sbsno,&
snowc, frland, adjsfcdsw, adjsfcnsw, adjsfcdlw, adjsfculw, &
snowc, frland, adjsfculw, &
adjnirbmu, adjnirdfu, adjvisbmu, adjvisdfu, adjnirbmd, &
adjnirdfd, adjvisbmd, adjvisdfd, xcosz, tseal, &
! adjnirdfd, adjvisbmd, adjvisdfd, gabsbdlw, xcosz, tseal, &
Expand All @@ -542,6 +542,7 @@ subroutine GFS_physics_driver &
! dtsfc_cice, dqsfc_cice, dusfc_cice, dvsfc_cice, ulwsfc_cice, &
!--- for CS-convection
wcbmax
real(kind=kind_phys), target, dimension(size(Grid%xlon,1)) :: adjsfcdlw, adjsfcdsw, adjsfcnsw

! 1 - land, 2 - ice, 3 - ocean
real(kind=kind_phys), dimension(size(Grid%xlon,1),3) :: &
Expand Down Expand Up @@ -672,6 +673,7 @@ subroutine GFS_physics_driver &
real :: pshltr,QCQ,rh02
real(kind=kind_phys), allocatable, dimension(:,:) :: den
real(kind=kind_phys), allocatable, dimension(:,:) :: dqdt_work
real(kind=kind_phys), pointer :: adjsfcdlw_for_lsm(:), adjsfcdsw_for_lsm(:), adjsfcnsw_for_lsm(:)
integer :: nwat

!! Initialize local variables (mainly for debugging purposes, because the
Expand Down Expand Up @@ -716,6 +718,16 @@ subroutine GFS_physics_driver &
dq3dt_initial = Diag%dq3dt
endif

if (Model%override_surface_radiative_fluxes) then
adjsfcdlw_for_lsm => Statein%adjsfcdlw_override
adjsfcdsw_for_lsm => Statein%adjsfcdsw_override
adjsfcnsw_for_lsm => Statein%adjsfcnsw_override
else
adjsfcdlw_for_lsm => adjsfcdlw
adjsfcdsw_for_lsm => adjsfcdsw
adjsfcnsw_for_lsm => adjsfcnsw
endif

!-------
! For COORDE-2019 averaging with fwindow, it was done before
! 3Diag fixes and averaging ingested using "fdaily"-factor
Expand Down Expand Up @@ -1437,11 +1449,10 @@ subroutine GFS_physics_driver &
! net = up - down = sfcemis * (sigma*T**4 - adjsfcdlw)

! --- ... define the downward lw flux absorbed by ground

do i=1,im
if (dry(i)) gabsbdlw3(i,1) = semis3(i,1) * adjsfcdlw(i)
if (icy(i)) gabsbdlw3(i,2) = semis3(i,2) * adjsfcdlw(i)
if (wet(i)) gabsbdlw3(i,3) = semis3(i,3) * adjsfcdlw(i)
if (dry(i)) gabsbdlw3(i,1) = semis3(i,1) * adjsfcdlw_for_lsm(i)
if (icy(i)) gabsbdlw3(i,2) = semis3(i,2) * adjsfcdlw_for_lsm(i)
if (wet(i)) gabsbdlw3(i,3) = semis3(i,3) * adjsfcdlw_for_lsm(i)
enddo

if (Model%lssav) then ! --- ... accumulate/save output variables
Expand Down Expand Up @@ -1502,11 +1513,21 @@ subroutine GFS_physics_driver &
! ' adjsfculw3=',adjsfculw3(ipr,:),' icefr=',Sfcprop%fice(ipr),' tsfc3=',tsfc3(ipr,:)
!
do i=1,im
Diag%dlwsfc(i) = Diag%dlwsfc(i) + adjsfcdlw(i)*dtf
Diag%dlwsfc(i) = Diag%dlwsfc(i) + adjsfcdlw_for_lsm(i)*dtf
Diag%ulwsfc(i) = Diag%ulwsfc(i) + adjsfculw(i)*dtf
Diag%psmean(i) = Diag%psmean(i) + Statein%pgr(i)*dtf ! mean surface pressure
enddo

if (Model%override_surface_radiative_fluxes) then
do i=1,im
Diag%dswsfc(i) = Diag%dswsfc(i) + adjsfcdsw_for_lsm(i)*dtf
Diag%uswsfc(i) = Diag%uswsfc(i) + (adjsfcdsw_for_lsm(i) - adjsfcnsw_for_lsm(i))*dtf

Diag%dlwsfc_rrtmg(i) = Diag%dlwsfc_rrtmg(i) + adjsfcdlw(i)*dtf
Diag%ulwsfc_rrtmg(i) = Diag%ulwsfc_rrtmg(i) + adjsfculw(i)*dtf
spencerkclark marked this conversation as resolved.
Show resolved Hide resolved
enddo
endif

if (Model%ldiag3d) then
if (Model%lsidea) then
do k=1,levs
Expand Down Expand Up @@ -1682,7 +1703,7 @@ subroutine GFS_physics_driver &
Statein%tgrs(:,1), Statein%qgrs(:,1,1), &
Sfcprop%tref, cd3(:,3), cdq3(:,3), Statein%prsl(:,1), &
work3, wet, Grid%xlon, Grid%sinlat, stress3(:,3), &
semis3(:,3), gabsbdlw3(:,3), adjsfcnsw, tprcp3(:,3), &
semis3(:,3), gabsbdlw3(:,3), adjsfcnsw_for_lsm, tprcp3(:,3), &
dtf, kdt, Model%solhr, xcosz, &
wind, flag_iter, &
flag_guess, Model%nstf_name, lprnt, ipr, &
Expand Down Expand Up @@ -1761,7 +1782,7 @@ subroutine GFS_physics_driver &
! --- inputs:
(im, lsoil, Statein%pgr, &
Statein%tgrs(:,1), Statein%qgrs(:,1,1), soiltyp, vegtype, &
sigmaf, semis3(:,1), gabsbdlw3(:,1), adjsfcdsw, adjsfcnsw, dtf,&
sigmaf, semis3(:,1), gabsbdlw3(:,1), adjsfcdsw_for_lsm, adjsfcnsw_for_lsm, dtf,&
! sigmaf, Radtend%semis, gabsbdlw, adjsfcdsw, adjsfcnsw, dtf, &
Sfcprop%tg3, cd3(:,1), cdq3(:,1), Statein%prsl(:,1), work3, &
Diag%zlvl, dry, wind, slopetyp, &
Expand Down Expand Up @@ -1791,7 +1812,7 @@ subroutine GFS_physics_driver &
! --- inputs:
(im, lsoil,kdt, Statein%pgr, Statein%ugrs, Statein%vgrs, &
Statein%tgrs, Statein%qgrs, soiltyp, vegtype, sigmaf, &
semis3(:,1), gabsbdlw3(:,1), adjsfcdsw, adjsfcnsw, dtf, &
semis3(:,1), gabsbdlw3(:,1), adjsfcdsw_for_lsm, adjsfcnsw_for_lsm, dtf, &
! Radtend%semis, gabsbdlw, adjsfcdsw, adjsfcnsw, dtf, &
Sfcprop%tg3, cd3(:,1), cdq3(:,1), Statein%prsl(:,1), work3,&
Diag%zlvl, dry, wind, slopetyp, &
Expand Down Expand Up @@ -1871,7 +1892,7 @@ subroutine GFS_physics_driver &
(im, lsoil, Statein%pgr, &
Statein%tgrs(:,1), Statein%qgrs(:,1,1), dtf, semis3(:,2), &
! Statein%tgrs(:,1), Statein%qgrs(:,1,1), dtf, Radtend%semis, &
gabsbdlw3(:,2), adjsfcnsw, adjsfcdsw, Sfcprop%srflag, &
gabsbdlw3(:,2), adjsfcnsw_for_lsm, adjsfcdsw_for_lsm, Sfcprop%srflag, &
cd3(:,2), cdq3(:,2), &
Statein%prsl(:,1), work3, islmsk, wind, &
flag_iter, lprnt, ipr, Model%min_lakeice, &
Expand Down Expand Up @@ -2074,17 +2095,26 @@ subroutine GFS_physics_driver &

do i=1,im
Diag%epi(i) = ep1d(i)
Diag%dlwsfci(i) = adjsfcdlw(i)
Diag%dlwsfci(i) = adjsfcdlw_for_lsm(i)
Diag%ulwsfci(i) = adjsfculw(i)
Diag%uswsfci(i) = adjsfcdsw(i) - adjsfcnsw(i)
Diag%dswsfci(i) = adjsfcdsw(i)
Diag%uswsfci(i) = adjsfcdsw_for_lsm(i) - adjsfcnsw_for_lsm(i)
Diag%dswsfci(i) = adjsfcdsw_for_lsm(i)
Diag%gfluxi(i) = gflx(i)
Diag%t1(i) = Statein%tgrs(i,1)
Diag%q1(i) = Statein%qgrs(i,1,1)
Diag%u1(i) = Statein%ugrs(i,1)
Diag%v1(i) = Statein%vgrs(i,1)
enddo

if (Model%override_surface_radiative_fluxes) then
do i=1,im
Diag%dlwsfci_rrtmg(i) = adjsfcdlw(i)
Diag%ulwsfci_rrtmg(i) = adjsfculw(i)
Diag%uswsfci_rrtmg(i) = adjsfcdsw(i) - adjsfcnsw(i)
Diag%dswsfci_rrtmg(i) = adjsfcdsw(i)
spencerkclark marked this conversation as resolved.
Show resolved Hide resolved
enddo
endif

! --- ... update near surface fields

call sfc_diag (im, Statein%pgr, Statein%ugrs(:,1), Statein%vgrs(:,1), &
Expand Down
Loading