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

Update calculation of time-averaged radiation variables #874

Closed
Closed
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
5 changes: 3 additions & 2 deletions atmos_model.F90
Original file line number Diff line number Diff line change
Expand Up @@ -975,6 +975,7 @@ subroutine update_atmos_model_state (Atmos, rc)
!--- local variables
integer :: i, localrc, sec_lastfhzerofh
integer :: isec, seconds, isec_fhzero
integer :: dtatm_temp
logical :: tmpflag_fhzero
real(kind=GFS_kind_phys) :: time_int, time_intfull
!
Expand Down Expand Up @@ -1021,11 +1022,11 @@ subroutine update_atmos_model_state (Atmos, rc)
endif
if (mpp_pe() == mpp_root_pe()) write(6,*) ' gfs diags time since last bucket empty: ',time_int/3600.,'hrs'
call atmosphere_nggps_diag(Atmos%Time)
call get_time ( Atmos%Time_step, dtatm_temp)
call fv3atm_diag_output(Atmos%Time, GFS_Diag, Atm_block, GFS_control%nx, GFS_control%ny, &
GFS_control%levs, 1, 1, 1.0_GFS_kind_phys, time_int, time_intfull, &
GFS_control%fhswr, GFS_control%fhlwr)
GFS_control%fhswr, GFS_control%fhlwr, dtatm_temp)
endif

!--- find current fhzero
if( GFS_Control%fhzero_array(1) > 0. ) then
fhzero_loop: do i=1,size(GFS_Control%fhzero_array)
Expand Down
64 changes: 55 additions & 9 deletions io/fv3atm_history_io.F90
Original file line number Diff line number Diff line change
Expand Up @@ -67,6 +67,7 @@ module fv3atm_history_io_mod
real(kind=kind_phys),dimension(:,:,:),pointer:: uwork3d => null()
logical :: uwork_set = .false.
character(128) :: uwindname = "(noname)"
real(4), dimension(:,:,:),allocatable :: rad_save, rad_dt_save, rad_accum

!--- miscellaneous other variables
logical :: use_wrtgridcomp_output = .FALSE.
Expand Down Expand Up @@ -115,20 +116,20 @@ end subroutine fv3atm_diag_register
!! This routine transfers diagnostic data to the FMS diagnostic
!! manager for eventual output to the history files.
subroutine fv3atm_diag_output(time, diag, atm_block, nx, ny, levs, ntcw, ntoz, &
dt, time_int, time_intfull, time_radsw, time_radlw)
dt, time_int, time_intfull, time_radsw, time_radlw,dt_atmos)
!--- subroutine interface variable definitions
type(time_type), intent(in) :: time
type(GFS_externaldiag_type), intent(in) :: diag(:)
type (block_control_type), intent(in) :: atm_block
integer, intent(in) :: nx, ny, levs, ntcw, ntoz
integer, intent(in) :: nx, ny, levs, ntcw, ntoz, dt_atmos
real(kind=kind_phys), intent(in) :: dt
real(kind=kind_phys), intent(in) :: time_int
real(kind=kind_phys), intent(in) :: time_intfull
real(kind=kind_phys), intent(in) :: time_radsw
real(kind=kind_phys), intent(in) :: time_radlw

call shared_history_data%output(time, diag, atm_block, nx, ny, levs, ntcw, ntoz, &
dt, time_int, time_intfull, time_radsw, time_radlw)
dt, time_int, time_intfull, time_radsw, time_radlw, dt_atmos)

end subroutine fv3atm_diag_output

Expand Down Expand Up @@ -268,9 +269,15 @@ subroutine history_type_register(hist, Diag, Time, Atm_block, Model, xlon, xlat,
allocate(hist%buffer_phys_bl(hist%isco:hist%ieco,hist%jsco:hist%jeco,nrgst_bl))
allocate(hist%buffer_phys_nb(hist%isco:hist%ieco,hist%jsco:hist%jeco,nrgst_nb))
allocate(hist%buffer_phys_windvect(3,hist%isco:hist%ieco,hist%jsco:hist%jeco,nrgst_vctbl))
allocate(hist%rad_save(hist%tot_diag_idx,hist%isco:hist%ieco,hist%jsco:hist%jeco))
allocate(hist%rad_dt_save(hist%tot_diag_idx,hist%isco:hist%ieco,hist%jsco:hist%jeco))
allocate(hist%rad_accum(hist%tot_diag_idx,hist%isco:hist%ieco,hist%jsco:hist%jeco))
hist%buffer_phys_bl = zero
hist%buffer_phys_nb = zero
hist%buffer_phys_windvect = zero
hist%rad_save = zero
hist%rad_dt_save = zero
hist%rad_accum = zero
if(mpp_pe() == mpp_root_pe()) print *,'in fv3atm_diag_register, nrgst_bl=',nrgst_bl,' nrgst_nb=',nrgst_nb, &
' nrgst_vctbl=',nrgst_vctbl, 'hist%isco=',hist%isco,hist%ieco,'hist%jsco=',hist%jsco,hist%jeco,' hist%num_axes_phys=', hist%num_axes_phys

Expand All @@ -282,7 +289,7 @@ end subroutine history_type_register
!! implementation of the public fv3atm_diag_output routine. Never
!! call this directly.
subroutine history_type_output(hist, time, diag, atm_block, nx, ny, levs, ntcw, ntoz, &
dt, time_int, time_intfull, time_radsw, time_radlw)
dt, time_int, time_intfull, time_radsw, time_radlw,dt_atmos)
!--- subroutine interface variable definitions
class(history_type) :: hist
type(time_type), intent(in) :: time
Expand All @@ -294,6 +301,7 @@ subroutine history_type_output(hist, time, diag, atm_block, nx, ny, levs, ntcw,
real(kind=kind_phys), intent(in) :: time_intfull
real(kind=kind_phys), intent(in) :: time_radsw
real(kind=kind_phys), intent(in) :: time_radlw
integer, intent(in) :: dt_atmos
!--- local variables
integer :: i, j, k, idx, nb, ix, ii, jj, levo_3d
character(len=2) :: xtra
Expand All @@ -305,13 +313,15 @@ subroutine history_type_output(hist, time, diag, atm_block, nx, ny, levs, ntcw,
real(kind=kind_phys), dimension(:,:,:), allocatable :: var3
#endif
real(kind=kind_phys) :: rtime_int, rtime_intfull, lcnvfac
real(kind=kind_phys) :: rtime_radsw, rtime_radlw
real(kind=kind_phys) :: rtime_radsw, rtime_radlw, this_radfh
real(kind=kind_phys) :: rad_old, rad_dt, rad_save

rtime_int = one/time_int
rtime_intfull = one/time_intfull
rtime_radsw = one/time_radsw
rtime_radlw = one/time_radlw

!intradrem = mod(time_int,time_radsw)
! if(mpp_pe()==mpp_root_pe())print *,'in,fv3atm_io. time avg, time_int=',time_int
history_loop: do idx = 1,hist%tot_diag_idx
has_id: if (diag(idx)%id > 0) then
Expand All @@ -321,13 +331,13 @@ subroutine history_type_output(hist, time, diag, atm_block, nx, ny, levs, ntcw,
lcnvfac = lcnvfac*rtime_intfull
! if(mpp_pe()==mpp_root_pe())print *,'in,fv3atm_io. full time avg, field=',trim(Diag(idx)%name),' time=',time_intfull
else if ( trim(diag(idx)%time_avg_kind) == 'rad_lw' ) then
lcnvfac = lcnvfac*min(rtime_radlw,rtime_int)
lcnvfac = lcnvfac*rtime_radlw/int(time_int/dt_atmos)
! if(mpp_pe()==mpp_root_pe())print *,'in,fv3atm_io. rad longwave avg, field=',trim(Diag(idx)%name),' time=',time_radlw
else if ( trim(diag(idx)%time_avg_kind) == 'rad_sw' ) then
lcnvfac = lcnvfac*min(rtime_radsw,rtime_int)
lcnvfac = lcnvfac*rtime_radsw/int(time_int/dt_atmos)
! if(mpp_pe()==mpp_root_pe())print *,'in,fv3atm_io. rad shortwave avg, field=',trim(Diag(idx)%name),' time=',time_radsw
else if ( trim(diag(idx)%time_avg_kind) == 'rad_swlw_min' ) then
lcnvfac = lcnvfac*min(max(rtime_radsw,rtime_radlw),rtime_int)
lcnvfac = lcnvfac*max(rtime_radsw,rtime_radlw)/int(time_int/dt_atmos)
! if(mpp_pe()==mpp_root_pe())print *,'in,fv3atm_io. rad swlw min avg, field=',trim(Diag(idx)%name),' time=',time_radlw,time_radsw,time_int
else
lcnvfac = lcnvfac*rtime_int
Expand Down Expand Up @@ -445,7 +455,43 @@ subroutine history_type_output(hist, time, diag, atm_block, nx, ny, levs, ntcw,
ii = i + Atm_block%isc -1
nb = Atm_block%blkno(ii,jj)
ix = Atm_block%ixp(ii,jj)
var2(i,j) = Diag(idx)%data(nb)%var2(ix)*lcnvfac
if ( trim(diag(idx)%time_avg_kind) == 'rad_sw' .or. trim(diag(idx)%time_avg_kind) == 'rad_lw' .or. trim(diag(idx)%time_avg_kind) == 'rad_swlw_min') then
select case (trim(diag(idx)%time_avg_kind))
case ('rad_sw')
this_radfh = time_radsw
case ('rad_lw')
this_radfh = time_radlw
case ('rad_swlw_min')
this_radfh = min(time_radsw,time_radlw)
case default

end select

if (time_int == dt_atmos ) then
var2(i,j) = Diag(idx)%data(nb)%var2(ix)*lcnvfac
hist%rad_accum(idx,ii,jj) = Diag(idx)%data(nb)%var2(ix)
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The computation of accumulated radiation or physics stays in CCPP physics (fluxr), I suggest making the accumulation related change in ccpp, not in the IO file (fv3atm_history_io.F90). We only apply time scale or mask update in the fv3atm_history_io.F90.

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I looked the CCPP code, and see that the fluxr is accumulated every radiation hour (3600s) (multiplied by coszdg(i) / coszen(i)) in the code:
physics/physics/Interstitials/UFS_SCM_NEPTUNE/GFS_rrtmg_post.F90
so we need modify the bucket weighting in order to calculate the dswrf_ave and ulwrf correctly


if (mod(time_int,this_radfh) == dt_atmos) then
hist%rad_dt_save(idx,ii,jj) = Diag(idx)%data(nb)%var2(ix)
hist%rad_save(idx,ii,jj) = Diag(idx)%data(nb)%var2(ix)
endif
elseif (mod(time_int,this_radfh) == dt_atmos) then
rad_save = hist%rad_save(idx,ii,jj)
hist%rad_save(idx,ii,jj) = Diag(idx)%data(nb)%var2(ix)
rad_dt = Diag(idx)%data(nb)%var2(ix) - rad_save
hist%rad_dt_save(idx,ii,jj) = rad_dt
rad_old = hist%rad_accum(idx,ii,jj)
var2(i,j) = (rad_old+rad_dt)*lcnvfac
hist%rad_accum(idx,ii,jj) = rad_old+rad_dt
else
rad_old = hist%rad_accum(idx,ii,jj)
rad_dt = hist%rad_dt_save(idx,ii,jj)
var2(i,j) = (rad_old+rad_dt) *lcnvfac
hist%rad_accum(idx,ii,jj) = rad_old+rad_dt
endif
else
var2(i,j) = Diag(idx)%data(nb)%var2(ix)*lcnvfac
endif
enddo
enddo
endif if_mask
Expand Down