Note

This page was generated from a jupyter notebook.

Green-Ampt infiltration and kinematic wave overland flow

This tutorial shows how to create a simple model of rainfall, infiltration, runoff, and overland flow, using two hydrologic components: SoilInfiltrationGreenAmpt and KinwaveImplicitOverlandFlow.

(Greg Tucker, September 2021)

[ ]:
from tqdm.notebook import trange

from landlab.components import KinwaveImplicitOverlandFlow, SoilInfiltrationGreenAmpt
from landlab.io import esri_ascii

Theory

The Green-Ampt method was introduced by Green and Ampt (1911) as a means of approximating the rate of water infiltration into soil from a layer of surface water. The method represents infiltration in terms of a wetting front that descends into the soil as infiltration progresses. A description of the method can be found in many hydrology textbooks, and in various online resources. The following is a brief summary, using the notation of Julien et al. (1995). The dimensions of each variable are indicated in square brackets, using the common convention that [L] means length, [M] is mass, and [T] is time.

The Green-Ampt method approximates the rate of water infiltration into the soil, \(f\) (dimensions of [L/T], representing water volume per unit surface area). Infiltration is driven by two effects: gravitational force, and downward suction (the “paper towel effect”) due to a gradient in moisture at the wetting front. The method treats the infiltration rate as a function of the following parameters:

  • \(K\) - saturated hydraulic conductivity [L/T]

  • \(H_f\) - capillary pressure head at the wetting front [L]

  • \(\phi\) - total soil porosity [-]

  • \(\theta_r\) - residual saturation [-]

  • \(\theta_e\) - effective porosity \(= \phi - \theta_r\) [-]

  • \(\theta_i\) - initial soil moisture content [-]

  • \(M_d\) - moisture deficit \(=\theta_e - \theta_i\) [-]

  • \(F\) - total infiltrated water depth [L]

The equation for infiltration rate is:

\[f = K \left( 1 + \frac{H_fM_d}{F} \right)\]

The first term in parentheses represents gravity and the second represents pore suction. If there were no pore suction effect, water would simply infiltrate downward at a rate equal to the hydraulic conductivity, \(K\). The suction effect increases this, but it becomes weaker as the cumulative infiltration depth \(F\) grows. Effectively, the second term approximates the pore-pressure gradient, which declines as the wetting front descends.

The version used in this component adds a term for the weight of the surface water with depth \(H\):

\[f = K \left( 1 + \frac{H_fM_d}{F} + \frac{H}{F} \right)\]

The component uses a simple forward-difference numerical scheme, with time step duration \(\Delta t\), in which the infiltration depth during one step is the lesser of the rate calculated above times \(\Delta t\), or the available surface water, \(H\):

\[\Delta F = \min( f\Delta t, H)\]

Note that the cumulative infitration \(F\) must be greater than zero in order to avoid division by zero; therefore, one should initialize the soil_water_infiltration__depth to a small positive value.

Example

Read in topography from a sample DEM

This is a lidar digital elevation model (DEM) from the West Bijou Creek escarpment on the Colorado High Plains, coarsened to 5 m grid resolution.

Note: it is convenient to use local grid coordinates rather than UTM coordinates, which are what the DEM provides. Therefore, after reading topography data into a grid called demgrid, which uses UTM coordinates, we copy over the elevation data into a second grid (grid) of the same dimensions that uses local coordinates (i.e., the lower left corner is (0, 0)).

[ ]:
# Create Landlab model grid and assign the DEM elevations to it,
# then display the terrain.
with open("bijou_gully_subset_5m_edit_dx_filled.asc") as fp:
    grid = esri_ascii.load(fp, name="topographic__elevation", at="node")
grid.imshow(grid.at_node["topographic__elevation"], colorbar_label="Elevation (m)")

Simulate a heavy 5-minute storm

The next bits of code use the SoilInfiltrationGreenAmpt and KinwaveImplicitOverlandFlow components to model infiltration and runoff during a 5-minute, 90 mm/hr storm.

[ ]:
# Create and initialize required input fields for infiltration
# component: depth of surface water, and depth (water volume per
# area) of infiltrated water.
depth = grid.add_zeros("surface_water__depth", at="node")
infilt = grid.add_zeros("soil_water_infiltration__depth", at="node")
infilt[:] = 1.0e-4  # small amount infiltrated (0.1 mm)

# Instantiate an infiltration component
ga = SoilInfiltrationGreenAmpt(grid)

# Instantiate an overland flow component
kw = KinwaveImplicitOverlandFlow(
    grid, runoff_rate=90.0, roughness=0.1, depth_exp=5.0 / 3.0
)
[ ]:
# Set time step and storm duration
dt = 10.0  # time step, sec
storm_duration = 300.0  # storm duration, sec

nsteps = int(storm_duration / dt)
[ ]:
# Run it for 10 minutes of heavy rain
for i in trange(nsteps):
    kw.run_one_step(dt)
    ga.run_one_step(dt)

Plot the cumulative infiltration

The plot below illustrates how the convergence of water in the branches of the gully network leads to greater infiltration, with less infiltration on steeper slopes and higher points in the landscape.

[ ]:
grid.imshow(1000.0 * infilt, colorbar_label="Infiltration depth (mm)", cmap="GnBu")

Optional parameters

The SoilInfiltrationGreenAmpt component provides a variety parameters that can be set by the user. A list and description of these can be found in the component’s __init__ docstring, which is printed below:

[ ]:
print(SoilInfiltrationGreenAmpt.__init__.__doc__)

References

Green, W. H., & Ampt, G. A. (1911). Studies on Soil Phyics. The Journal of Agricultural Science, 4(1), 1-24.

Julien, P. Y., Saghafian, B., and Ogden, F. L. (1995) Raster-based hydrologic modeling of spatially-varied surface runoff, J. Am. Water Resour. As., 31, 523–536, doi:10.1111/j.17521688.1995.tb04039.x.

Rengers, F. K., McGuire, L. A., Kean, J. W., Staley, D. M., and Hobley, D. (2016) Model simulations of flood and debris flow timing in steep catchments after wildfire, Water Resour. Res., 52, 6041–6061, doi:10.1002/2015WR018176.


Generated by nbsphinx from a Jupyter notebook.