Skip to content

Commit

Permalink
Move setting of qflx_infl to its own subroutine
Browse files Browse the repository at this point in the history
This is cleaner and safer than having it set initially in one routine
and then updated in a later routine.
  • Loading branch information
billsacks committed Dec 28, 2017
1 parent 6c979ad commit 5c75419
Show file tree
Hide file tree
Showing 3 changed files with 46 additions and 13 deletions.
6 changes: 5 additions & 1 deletion src/biogeophys/HydrologyNoDrainageMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -62,7 +62,8 @@ subroutine HydrologyNoDrainage(bounds, &
use SnowHydrologyMod , only : SnowCompaction, CombineSnowLayers, DivideSnowLayers, SnowCapping
use SnowHydrologyMod , only : SnowWater, BuildSnowFilter
use SoilHydrologyMod , only : CLMVICMap, SetSoilWaterFractions
use SoilHydrologyMod , only : SetQflxInputs, RouteInfiltrationExcess, UpdateH2osfc, TotalSurfaceRunoff
use SoilHydrologyMod , only : SetQflxInputs, RouteInfiltrationExcess, UpdateH2osfc
use SoilHydrologyMod , only : Infiltration, TotalSurfaceRunoff
use SoilHydrologyMod , only : UpdateUrbanPonding
use SoilHydrologyMod , only : WaterTable, PerchedWaterTable
use SoilHydrologyMod , only : ThetaBasedWaterTable, RenewCondensation
Expand Down Expand Up @@ -202,6 +203,9 @@ subroutine HydrologyNoDrainage(bounds, &
energyflux_inst, soilhydrology_inst, &
waterflux_inst, waterstate_inst)

call Infiltration(bounds, num_hydrologyc, filter_hydrologyc, &
waterflux_inst)

call TotalSurfaceRunoff(bounds, num_hydrologyc, filter_hydrologyc, &
num_urbanc, filter_urbanc, &
waterflux_inst, soilhydrology_inst, saturated_excess_runoff_inst, waterstate_inst)
Expand Down
49 changes: 37 additions & 12 deletions src/biogeophys/SoilHydrologyMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -38,6 +38,7 @@ module SoilHydrologyMod
public :: SetSoilWaterFractions ! Set diagnostic variables related to the fraction of water and ice in each layer
public :: SetQflxInputs ! Set the flux of water into the soil from the top
public :: UpdateH2osfc ! Calculate fluxes out of h2osfc and update the h2osfc state
public :: Infiltration ! Calculate total infiltration
public :: TotalSurfaceRunoff ! Calculate total surface runoff
public :: UpdateUrbanPonding ! Update the state variable representing ponding on urban surfaces
public :: WaterTable ! Calculate water table before imposing drainage
Expand Down Expand Up @@ -281,7 +282,7 @@ subroutine RouteInfiltrationExcess(bounds, num_hydrologyc, filter_hydrologyc, &
!-----------------------------------------------------------------------

associate( &
qflx_infl => waterflux_inst%qflx_infl_col , & ! Output: [real(r8) (:) ] infiltration (mm H2O /s)
qflx_in_soil_limited => waterflux_inst%qflx_in_soil_limited_col , & ! Output: [real(r8) (:) ] surface input to soil, limited by max infiltration rate (mm H2O /s)
qflx_in_h2osfc => waterflux_inst%qflx_in_h2osfc_col , & ! Output: [real(r8) (:) ] total surface input to h2osfc
qflx_infl_excess_surf => waterflux_inst%qflx_infl_excess_surf_col , & ! Output: [real(r8) (:) ] surface runoff due to infiltration excess (mm H2O /s)
qflx_in_soil => waterflux_inst%qflx_in_soil_col , & ! Input: [real(r8) (:) ] surface input to soil (mm/s)
Expand All @@ -295,7 +296,7 @@ subroutine RouteInfiltrationExcess(bounds, num_hydrologyc, filter_hydrologyc, &
do fc = 1, num_hydrologyc
c = filter_hydrologyc(fc)
if (lun%itype(col%landunit(c)) == istsoil .or. lun%itype(col%landunit(c))==istcrop) then
qflx_infl(c) = qflx_in_soil(c) - qflx_infl_excess(c)
qflx_in_soil_limited(c) = qflx_in_soil(c) - qflx_infl_excess(c)
if (h2osfcflag /= 0) then
qflx_in_h2osfc(c) = qflx_top_soil_to_h2osfc(c) + qflx_infl_excess(c)
qflx_infl_excess_surf(c) = 0._r8
Expand All @@ -306,7 +307,7 @@ subroutine RouteInfiltrationExcess(bounds, num_hydrologyc, filter_hydrologyc, &
end if
else
! non-vegetated landunits (i.e. urban) use original CLM4 code
qflx_infl(c) = qflx_in_soil(c)
qflx_in_soil_limited(c) = qflx_in_soil(c)
qflx_in_h2osfc(c) = 0._r8
qflx_infl_excess_surf(c) = 0._r8
end if
Expand Down Expand Up @@ -344,9 +345,6 @@ subroutine UpdateH2osfc(bounds, num_hydrologyc, filter_hydrologyc, &
integer :: c,l,fc ! indices
real(r8) :: dtime ! land model time step (sec)
real(r8) :: h2osfc_partial(bounds%begc:bounds%endc) ! partially-updated h2osfc
! FIXME(wjs, 2017-08-22) put qflx_h2osfc_drain in waterflux_type (this will be
! needed when we move the update of qflx_infl to a separate routine)
real(r8) :: qflx_h2osfc_drain(bounds%begc:bounds%endc) ! bottom drainage from h2osfc
logical :: truncate_h2osfc_to_zero(bounds%begc:bounds%endc)
!-----------------------------------------------------------------------

Expand All @@ -359,7 +357,7 @@ subroutine UpdateH2osfc(bounds, num_hydrologyc, filter_hydrologyc, &

qflx_in_h2osfc => waterflux_inst%qflx_in_h2osfc_col , & ! Input: [real(r8) (:) ] total surface input to h2osfc
qflx_h2osfc_surf => waterflux_inst%qflx_h2osfc_surf_col , & ! Output: [real(r8) (:) ] surface water runoff (mm H2O /s)
qflx_infl => waterflux_inst%qflx_infl_col , & ! Output: [real(r8) (:) ] infiltration (mm H2O /s)
qflx_h2osfc_drain => waterflux_inst%qflx_h2osfc_drain_col , & ! Output: [real(r8) (:) ] bottom drainage from h2osfc (mm H2O /s)

h2osfc_thresh => soilhydrology_inst%h2osfc_thresh_col, & ! Input: [real(r8) (:) ] level at which h2osfc "percolates"
h2osfcflag => soilhydrology_inst%h2osfcflag & ! Input: integer
Expand Down Expand Up @@ -415,16 +413,43 @@ subroutine UpdateH2osfc(bounds, num_hydrologyc, filter_hydrologyc, &

end do

! FIXME(wjs, 2017-08-22) Move this update to a separate routine that calculations
! qflx_infl from its components.
end associate

end subroutine UpdateH2osfc

!-----------------------------------------------------------------------
subroutine Infiltration(bounds, num_hydrologyc, filter_hydrologyc, &
waterflux_inst)
!
! !DESCRIPTION:
! Calculate total infiltration
!
! !ARGUMENTS:
type(bounds_type) , intent(in) :: bounds
integer , intent(in) :: num_hydrologyc ! number of column soil points in column filter
integer , intent(in) :: filter_hydrologyc(:) ! column filter for soil points
type(waterflux_type) , intent(inout) :: waterflux_inst
!
! !LOCAL VARIABLES:
integer :: fc, c

character(len=*), parameter :: subname = 'Infiltration'
!-----------------------------------------------------------------------

associate( &
qflx_infl => waterflux_inst%qflx_infl_col , & ! Output: [real(r8) (:) ] infiltration (mm H2O /s)
qflx_in_soil_limited => waterflux_inst%qflx_in_soil_limited_col , & ! Input: [real(r8) (:) ] surface input to soil, limited by max infiltration rate (mm H2O /s)
qflx_h2osfc_drain => waterflux_inst%qflx_h2osfc_drain_col & ! Input: [real(r8) (:) ] bottom drainage from h2osfc (mm H2O /s)
)

do fc = 1, num_hydrologyc
c = filter_hydrologyc(fc)
qflx_infl(c) = qflx_infl(c) + qflx_h2osfc_drain(c)
qflx_infl(c) = qflx_in_soil_limited(c) + qflx_h2osfc_drain(c)
end do

end associate
end associate

end subroutine UpdateH2osfc
end subroutine Infiltration

!-----------------------------------------------------------------------
subroutine QflxH2osfcSurf(bounds, num_hydrologyc, filter_hydrologyc, &
Expand Down
4 changes: 4 additions & 0 deletions src/biogeophys/WaterfluxType.F90
Original file line number Diff line number Diff line change
Expand Up @@ -81,6 +81,8 @@ module WaterfluxType
real(r8), pointer :: qflx_rain_plus_snomelt_col(:) ! col rain plus snow melt falling on the soil (mm/s)
real(r8), pointer :: qflx_top_soil_col (:) ! col net water input into soil from top (mm/s)
real(r8), pointer :: qflx_in_soil_col (:) ! col surface input to soil (mm/s)
real(r8), pointer :: qflx_in_soil_limited_col (:) ! col surface input to soil, limited by max infiltration rate (mm/s)
real(r8), pointer :: qflx_h2osfc_drain_col (:) ! col bottom drainage from h2osfc (mm/s)
real(r8), pointer :: qflx_top_soil_to_h2osfc_col(:) ! col portion of qflx_top_soil going to h2osfc, minus evaporation (mm/s)
real(r8), pointer :: qflx_in_h2osfc_col(:) ! col total surface input to h2osfc
real(r8), pointer :: qflx_h2osfc_to_ice_col (:) ! col conversion of h2osfc to ice
Expand Down Expand Up @@ -226,6 +228,8 @@ subroutine InitAllocate(this, bounds)
allocate(this%qflx_rain_plus_snomelt_col(begc:endc)) ; this%qflx_rain_plus_snomelt_col(:) = nan
allocate(this%qflx_top_soil_col (begc:endc)) ; this%qflx_top_soil_col (:) = nan
allocate(this%qflx_in_soil_col (begc:endc)) ; this%qflx_in_soil_col (:) = nan
allocate(this%qflx_in_soil_limited_col (begc:endc)) ; this%qflx_in_soil_limited_col (:) = nan
allocate(this%qflx_h2osfc_drain_col (begc:endc)) ; this%qflx_h2osfc_drain_col (:) = nan
allocate(this%qflx_top_soil_to_h2osfc_col(begc:endc)) ; this%qflx_top_soil_to_h2osfc_col(:) = nan
allocate(this%qflx_in_h2osfc_col (begc:endc)) ; this%qflx_in_h2osfc_col(:) = nan
allocate(this%qflx_h2osfc_to_ice_col (begc:endc)) ; this%qflx_h2osfc_to_ice_col (:) = nan
Expand Down

0 comments on commit 5c75419

Please sign in to comment.