Skip to content

Commit

Permalink
Fix for updating stochastic physics on separate time-step. (#199)
Browse files Browse the repository at this point in the history
This bug fix allows the random patterns in the stochastic physics persist the for a period of time (defined as SKEBINT,SPPTINT, etc.) before calculating new patterns.
The fix is to move the allocation of the saved variables into the init section of stochastic_physics_wrapper, and remove the deallocates in the run section.
  • Loading branch information
pjpegion committed Nov 30, 2020
1 parent fd668c4 commit a898cca
Show file tree
Hide file tree
Showing 2 changed files with 68 additions and 48 deletions.
5 changes: 4 additions & 1 deletion atmos_model.F90
Original file line number Diff line number Diff line change
Expand Up @@ -100,7 +100,7 @@ module atmos_model_mod
use IPD_driver, only: IPD_initialize, IPD_initialize_rst
use CCPP_driver, only: CCPP_step, non_uniform_blocks

use stochastic_physics_wrapper_mod, only: stochastic_physics_wrapper
use stochastic_physics_wrapper_mod, only: stochastic_physics_wrapper,stochastic_physics_wrapper_end
#else
use IPD_driver, only: IPD_initialize, IPD_initialize_rst, IPD_step
use physics_abstraction_layer, only: time_vary_step, radiation_step1, physics_step1, physics_step2
Expand Down Expand Up @@ -962,6 +962,9 @@ subroutine atmos_model_end (Atmos)
!---- termination routine for atmospheric model ----

call atmosphere_end (Atmos % Time, Atmos%grid, restart_endfcst)

call stochastic_physics_wrapper_end(IPD_Control)

if(restart_endfcst) then
call FV3GFS_restart_write (IPD_Data, IPD_Restart, Atm_block, &
IPD_Control, Atmos%domain)
Expand Down
111 changes: 64 additions & 47 deletions stochastic_physics/stochastic_physics_wrapper.F90
Original file line number Diff line number Diff line change
Expand Up @@ -92,25 +92,43 @@ subroutine stochastic_physics_wrapper (GFS_Control, GFS_Data, Atm_block, ierr)
return
endif
end if
allocate(xlat(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
allocate(xlon(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
if (GFS_Control%do_sppt) then
allocate(sppt_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
end if
if (GFS_Control%do_shum) then
allocate(shum_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
end if
if (GFS_Control%do_skeb) then
allocate(skebu_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
allocate(skebv_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
end if
if ( GFS_Control%lndp_type .EQ. 2 ) then ! this scheme updates through forecast
allocate(sfc_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%n_var_lndp))
end if
if (GFS_Control%lndp_type .EQ. 2) then ! save wts, and apply lndp scheme
allocate(smc(1:Atm_block%nblks,maxval(GFS_Control%blksz),GFS_Control%lsoil))
allocate(slc(1:Atm_block%nblks,maxval(GFS_Control%blksz),GFS_Control%lsoil))
allocate(stc(1:Atm_block%nblks,maxval(GFS_Control%blksz),GFS_Control%lsoil))
allocate(stype(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
allocate(vfrac(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
endif

do nb=1,Atm_block%nblks
xlat(nb,1:GFS_Control%blksz(nb)) = GFS_Data(nb)%Grid%xlat(:)
xlon(nb,1:GFS_Control%blksz(nb)) = GFS_Data(nb)%Grid%xlon(:)
end do

if ( GFS_Control%lndp_type .EQ. 1 ) then ! this scheme sets perts once
! Copy blocked data into contiguous arrays; no need to copy sfc_wts in (intent out)
allocate(xlat(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
allocate(xlon(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
allocate(sfc_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),GFS_Control%n_var_lndp))
do nb=1,Atm_block%nblks
xlat(nb,1:GFS_Control%blksz(nb)) = GFS_Data(nb)%Grid%xlat(:)
xlon(nb,1:GFS_Control%blksz(nb)) = GFS_Data(nb)%Grid%xlon(:)
end do
call run_stochastic_physics(GFS_Control%levs, GFS_Control%kdt, GFS_Control%phour, GFS_Control%blksz, xlat=xlat, xlon=xlon, &
sppt_wts=sppt_wts, shum_wts=shum_wts, skebu_wts=skebu_wts, skebv_wts=skebv_wts, sfc_wts=sfc_wts, &
nthreads=nthreads)
! Copy contiguous data back; no need to copy xlat/xlon, these are intent(in) - just deallocate
do nb=1,Atm_block%nblks
GFS_Data(nb)%Coupling%sfc_wts(:,:) = sfc_wts(nb,1:GFS_Control%blksz(nb),:)
end do
deallocate(xlat)
deallocate(xlon)
deallocate(sfc_wts)
end if
! Consistency check for cellular automata
Expand All @@ -126,27 +144,6 @@ subroutine stochastic_physics_wrapper (GFS_Control, GFS_Data, Atm_block, ierr)
else initalize_stochastic_physics

if (GFS_Control%do_sppt .OR. GFS_Control%do_shum .OR. GFS_Control%do_skeb .OR. (GFS_Control%lndp_type .EQ. 2) ) then
! Copy blocked data into contiguous arrays; no need to copy weights in (intent(out))
allocate(xlat(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
allocate(xlon(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
do nb=1,Atm_block%nblks
xlat(nb,1:GFS_Control%blksz(nb)) = GFS_Data(nb)%Grid%xlat(:)
xlon(nb,1:GFS_Control%blksz(nb)) = GFS_Data(nb)%Grid%xlon(:)
end do
if (GFS_Control%do_sppt) then
allocate(sppt_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
end if
if (GFS_Control%do_shum) then
allocate(shum_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
end if
if (GFS_Control%do_skeb) then
allocate(skebu_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
allocate(skebv_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%levs))
end if
if ( GFS_Control%lndp_type .EQ. 2 ) then ! this scheme updates through forecast
allocate(sfc_wts(1:Atm_block%nblks,maxval(GFS_Control%blksz),1:GFS_Control%n_var_lndp))
end if

call run_stochastic_physics(GFS_Control%levs, GFS_Control%kdt, GFS_Control%phour, GFS_Control%blksz, xlat=xlat, xlon=xlon, &
sppt_wts=sppt_wts, shum_wts=shum_wts, skebu_wts=skebu_wts, skebv_wts=skebv_wts, sfc_wts=sfc_wts, &
nthreads=nthreads)
Expand All @@ -155,32 +152,23 @@ subroutine stochastic_physics_wrapper (GFS_Control, GFS_Data, Atm_block, ierr)
do nb=1,Atm_block%nblks
GFS_Data(nb)%Coupling%sppt_wts(:,:) = sppt_wts(nb,1:GFS_Control%blksz(nb),:)
end do
deallocate(sppt_wts)
end if
if (GFS_Control%do_shum) then
do nb=1,Atm_block%nblks
GFS_Data(nb)%Coupling%shum_wts(:,:) = shum_wts(nb,1:GFS_Control%blksz(nb),:)
end do
deallocate(shum_wts)
end if
if (GFS_Control%do_skeb) then
do nb=1,Atm_block%nblks
GFS_Data(nb)%Coupling%skebu_wts(:,:) = skebu_wts(nb,1:GFS_Control%blksz(nb),:)
GFS_Data(nb)%Coupling%skebv_wts(:,:) = skebv_wts(nb,1:GFS_Control%blksz(nb),:)
end do
deallocate(skebu_wts)
deallocate(skebv_wts)
end if
if (GFS_Control%lndp_type .EQ. 2) then ! save wts, and apply lndp scheme
do nb=1,Atm_block%nblks
GFS_Data(nb)%Coupling%sfc_wts(:,:) = sfc_wts(nb,1:GFS_Control%blksz(nb),:)
end do

allocate(smc(1:Atm_block%nblks,maxval(GFS_Control%blksz),GFS_Control%lsoil))
allocate(slc(1:Atm_block%nblks,maxval(GFS_Control%blksz),GFS_Control%lsoil))
allocate(stc(1:Atm_block%nblks,maxval(GFS_Control%blksz),GFS_Control%lsoil))
allocate(stype(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
allocate(vfrac(1:Atm_block%nblks,maxval(GFS_Control%blksz)))
do nb=1,Atm_block%nblks
stype(nb,1:GFS_Control%blksz(nb)) = GFS_Data(nb)%Sfcprop%stype(:)
smc(nb,1:GFS_Control%blksz(nb),:) = GFS_Data(nb)%Sfcprop%smc(:,:)
Expand All @@ -202,21 +190,13 @@ subroutine stochastic_physics_wrapper (GFS_Control, GFS_Data, Atm_block, ierr)
write(6,*) 'call to GFS_apply_lndp failed'
return
endif
deallocate(stype)
deallocate(sfc_wts)
do nb=1,Atm_block%nblks
GFS_Data(nb)%Sfcprop%smc(:,:) = smc(nb,1:GFS_Control%blksz(nb),:)
GFS_Data(nb)%Sfcprop%slc(:,:) = slc(nb,1:GFS_Control%blksz(nb),:)
GFS_Data(nb)%Sfcprop%stc(:,:) = stc(nb,1:GFS_Control%blksz(nb),:)
GFS_Data(nb)%Sfcprop%vfrac(:) = vfrac(nb,1:GFS_Control%blksz(nb))
enddo
deallocate(smc)
deallocate(slc)
deallocate(stc)
deallocate(vfrac)
endif ! lndp block
deallocate(xlat)
deallocate(xlon)
end if

endif initalize_stochastic_physics
Expand Down Expand Up @@ -309,4 +289,41 @@ subroutine stochastic_physics_wrapper (GFS_Control, GFS_Data, Atm_block, ierr)

end subroutine stochastic_physics_wrapper


subroutine stochastic_physics_wrapper_end (GFS_Control)

use GFS_typedefs, only: GFS_control_type, GFS_data_type
use stochastic_physics, only: finalize_stochastic_physics

implicit none

type(GFS_control_type), intent(inout) :: GFS_Control

if (GFS_Control%do_sppt .OR. GFS_Control%do_shum .OR. GFS_Control%do_skeb .OR. (GFS_Control%lndp_type .GT. 0) ) then
if (allocated(xlat)) deallocate(xlat)
if (allocated(xlon)) deallocate(xlon)
if (GFS_Control%do_sppt) then
if (allocated(sppt_wts)) deallocate(sppt_wts)
end if
if (GFS_Control%do_shum) then
if (allocated(shum_wts)) deallocate(shum_wts)
end if
if (GFS_Control%do_skeb) then
if (allocated(skebu_wts)) deallocate(skebu_wts)
if (allocated(skebv_wts)) deallocate(skebv_wts)
end if
if ( GFS_Control%lndp_type .EQ. 2 ) then ! this scheme updates through forecast
if (allocated(sfc_wts)) deallocate(sfc_wts)
end if
if (GFS_Control%lndp_type .EQ. 2) then ! save wts, and apply lndp scheme
if (allocated(smc)) deallocate(smc)
if (allocated(slc)) deallocate(slc)
if (allocated(stc)) deallocate(stc)
if (allocated(stype)) deallocate(stype)
if (allocated(vfrac)) deallocate(vfrac)
endif
call finalize_stochastic_physics()
endif
end subroutine stochastic_physics_wrapper_end

end module stochastic_physics_wrapper_mod

0 comments on commit a898cca

Please sign in to comment.