OverlandFlow: Model shallow water flow over topography using the numerical approximation of de Almeida#

Landlab component that simulates overland flow.

This component simulates overland flow using the 2-D numerical model of shallow-water flow over topography using the de Almeida et al., 2012 algorithm for storage-cell inundation modeling.

Code author: Jordan Adams

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> from landlab.components.overland_flow import OverlandFlow

Create a grid on which to calculate overland flow.

>>> grid = RasterModelGrid((4, 5))

The grid will need some data to provide the overland flow component. To check the names of the fields that provide input to the overland flow component use the input_var_names class property.

>>> OverlandFlow.input_var_names
('surface_water__depth', 'topographic__elevation')

Create fields of data for each of these input variables.

>>> grid.at_node["topographic__elevation"] = [
...     [0.0, 0.0, 0.0, 0.0, 0.0],
...     [1.0, 1.0, 1.0, 1.0, 1.0],
...     [2.0, 2.0, 2.0, 2.0, 2.0],
...     [3.0, 3.0, 3.0, 3.0, 3.0],
... ]
>>> grid.at_node["surface_water__depth"] = [
...     [0.0, 0.0, 0.0, 0.0, 0.0],
...     [0.0, 0.0, 0.0, 0.0, 0.0],
...     [0.0, 0.0, 0.0, 0.0, 0.0],
...     [0.1, 0.1, 0.1, 0.1, 0.1],
... ]

Instantiate the OverlandFlow component to work on this grid, and run it.

>>> of = OverlandFlow(grid, steep_slopes=True)
>>> of.run_one_step()

After calculating the overland flow, new fields have been added to the grid. Use the output_var_names property to see the names of the fields that have been changed.

>>> of.output_var_names
('surface_water__depth', 'surface_water__discharge', 'water_surface__gradient')

The surface_water__depth field is defined at nodes.

>>> of.var_loc("surface_water__depth")
'node'
>>> grid.at_node["surface_water__depth"]
array([  1.00000000e-05,   1.00000000e-05,   1.00000000e-05,
         1.00000000e-05,   1.00000000e-05,   1.00000000e-05,
         1.00000000e-05,   1.00000000e-05,   1.00000000e-05,
         1.00000000e-05,   1.00000000e-05,   2.00100000e-02,
         2.00100000e-02,   2.00100000e-02,   1.00000000e-05,
         1.00010000e-01,   1.00010000e-01,   1.00010000e-01,
         1.00010000e-01,   1.00010000e-01])

The surface_water__discharge field is defined at links. Because our initial topography was a dipping plane, there is no water discharge in the horizontal direction, only toward the bottom of the grid.

>>> of.var_loc("surface_water__discharge")
'link'
>>> q = grid.at_link["surface_water__discharge"]
>>> np.all(q[grid.horizontal_links] == 0.0)
True
>>> np.all(q[grid.vertical_links] <= 0.0)
True

The water_surface__gradient is also defined at links.

>>> of.var_loc("water_surface__gradient")
'link'
>>> grid.at_link["water_surface__gradient"]
array([ 0. ,  0. ,  0. ,  0. ,
        0. ,  1. ,  1. ,  1. ,  0. ,
        0. ,  0. ,  0. ,  0. ,
        0. ,  1. ,  1. ,  1. ,  0. ,
        0. ,  0. ,  0. ,  0. ,
        0. ,  1.1,  1.1,  1.1,  0. ,
        0. ,  0. ,  0. ,  0. ])
class OverlandFlow(*args, **kwds)[source]#

Bases: Component

Simulate overland flow using de Almeida approximations.

Landlab component that simulates overland flow using the de Almeida et al., 2012 approximations of the 1D shallow water equations to be used for 2D flood inundation modeling.

This component calculates discharge, depth and shear stress after some precipitation event across any raster grid. Default input file is named “overland_flow_input.txt’ and is contained in the landlab.components.overland_flow folder.

The primary method of this class is run_one_step.

References

Required Software Citation(s) Specific to this Component

Adams, J., Gasparini, N., Hobley, D., Tucker, G., Hutton, E., Nudurupati, S., Istanbulluoglu, E. (2017). The Landlab v1. 0 OverlandFlow component: a Python tool for computing shallow-water flow across watersheds. Geoscientific Model Development 10(4), 1645. https://dx.doi.org/10.5194/gmd-10-1645-2017

Additional References

de Almeida, G., Bates, P., Freer, J., Souvignet, M. (2012). Improving the stability of a simple formulation of the shallow water equations for 2-D flood modeling. Water Resources Research 48(5) https://dx.doi.org/10.1029/2011wr011570

Create an overland flow component.

Parameters:
  • grid (RasterModelGrid) – A landlab grid.

  • h_init (float, optional) – Thicknes of initial thin layer of water to prevent divide by zero errors (m).

  • alpha (float, optional) – Time step coeffcient, described in Bates et al., 2010 and de Almeida et al., 2012.

  • mannings_n (float, optional) – Manning’s roughness coefficient.

  • g (float, optional) – Acceleration due to gravity (m/s^2).

  • theta (float, optional) – Weighting factor from de Almeida et al., 2012.

  • rainfall_intensity (float or array of float, optional) – Rainfall intensity. Default is zero.

  • steep_slopes (bool, optional) – Modify the algorithm to handle steeper slopes at the expense of speed. If model runs become unstable, consider setting to True.

__init__(grid, default_fixed_links=False, h_init=1e-05, alpha=0.7, mannings_n=0.03, g=9.80665, theta=0.8, rainfall_intensity=0.0, steep_slopes=False)[source]#

Create an overland flow component.

Parameters:
  • grid (RasterModelGrid) – A landlab grid.

  • h_init (float, optional) – Thicknes of initial thin layer of water to prevent divide by zero errors (m).

  • alpha (float, optional) – Time step coeffcient, described in Bates et al., 2010 and de Almeida et al., 2012.

  • mannings_n (float, optional) – Manning’s roughness coefficient.

  • g (float, optional) – Acceleration due to gravity (m/s^2).

  • theta (float, optional) – Weighting factor from de Almeida et al., 2012.

  • rainfall_intensity (float or array of float, optional) – Rainfall intensity. Default is zero.

  • steep_slopes (bool, optional) – Modify the algorithm to handle steeper slopes at the expense of speed. If model runs become unstable, consider setting to True.

calc_time_step()[source]#

Calculate time step.

Adaptive time stepper from Bates et al., 2010 and de Almeida et al., 2012

discharge_mapper(input_discharge, convert_to_volume=False)[source]#

Maps discharge value from links onto nodes.

This method takes the discharge values on links and determines the links that are flowing INTO a given node. The fluxes moving INTO a given node are summed.

This method ignores all flow moving OUT of a given node.

This takes values from the OverlandFlow component (by default) in units of [L^2/T]. If the convert_to_cms flag is raised as True, this method converts discharge to units [L^3/T] - as of Aug 2016, only operates for square RasterModelGrid instances.

The output array is of length grid.number_of_nodes and can be used with the Landlab imshow_grid plotter.

Returns a numpy array (discharge_vals)

property dt#

Component timestep.

Type:

dt

property h#

The depth of water at each node.

overland_flow(dt=None)[source]#

Generate overland flow across a grid.

For one time step, this generates ‘overland flow’ across a given grid by calculating discharge at each node.

Using the depth slope product, shear stress is calculated at every node.

Outputs water depth, discharge and shear stress values through time at every point in the input grid.

property rainfall_intensity#

the rainfall rate [m/s]

Must be positive.

Type:

rainfall_intensity

run_one_step(dt=None)[source]#

Generate overland flow across a grid.

For one time step, this generates ‘overland flow’ across a given grid by calculating discharge at each node.

Using the depth slope product, shear stress is calculated at every node.

Outputs water depth, discharge and shear stress values through time at every point in the input grid.

set_up_neighbor_arrays()[source]#

Create and initialize link neighbor arrays.

Set up arrays of neighboring horizontal and vertical links that are needed for the de Almeida solution.

Find active link neighbors for every fixed link.

Specialized link ID function used to ID the active links that neighbor fixed links in the vertical and horizontal directions.

If the user wants to assign fixed gradients or values to the fixed links dynamically, this function identifies the nearest active_link neighbor.

Each fixed link can either have 0 or 1 active neighbor. This function finds if and where that active neighbor is and stores those IDs in an array.

Parameters:

grid (RasterModelGrid) – A landlab grid.

Returns:

Flat array of links.

Return type:

ndarray of int, shape (*, )

Examples

>>> from landlab import NodeStatus, RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.status_at_node[:5] = NodeStatus.FIXED_GRADIENT
>>> grid.status_at_node[::5] = NodeStatus.FIXED_GRADIENT
>>> grid.status_at_node.reshape(grid.shape)
array([[2, 2, 2, 2, 2],
       [2, 0, 0, 0, 1],
       [2, 0, 0, 0, 1],
       [2, 1, 1, 1, 1]], dtype=uint8)
>>> grid.fixed_links
array([ 5,  6,  7,  9, 18])
>>> grid.active_links
array([10, 11, 12, 14, 15, 16, 19, 20, 21, 23, 24, 25])
>>> find_active_neighbors_for_fixed_links(grid)
array([14, 15, 16, 10, 19])
>>> rmg = RasterModelGrid((4, 7))
>>> rmg.at_node["topographic__elevation"] = rmg.zeros(at="node")
>>> rmg.at_link["topographic__slope"] = rmg.zeros(at="link")
>>> rmg.status_at_node[rmg.perimeter_nodes] = rmg.BC_NODE_IS_FIXED_GRADIENT
>>> find_active_neighbors_for_fixed_links(rmg)
array([20, 21, 22, 23, 24, 14, 17, 27, 30, 20, 21, 22, 23, 24])

OverlandFlowBates: Model shallow water flow over topography using the numerical approximation of Bates#

generate_overland_flow.py.

This component simulates overland flow using the 2-D numerical model of shallow-water flow over topography using the Bates et al. (2010) algorithm for storage-cell inundation modeling.

Written by Jordan Adams, based on code written by Greg Tucker.

Last updated: April 21, 2016

class OverlandFlowBates(*args, **kwds)[source]#

Bases: Component

Simulate overland flow using Bates et al. (2010).

Landlab component that simulates overland flow using the Bates et al., (2010) approximations of the 1D shallow water equations to be used for 2D flood inundation modeling.

This component calculates discharge, depth and shear stress after some precipitation event across any raster grid. Default input file is named “overland_flow_input.txt’ and is contained in the landlab.components.overland_flow folder.

Parameters:
  • grid (RasterGridModel) – A grid.

  • input_file (str) –

    Contains necessary and optional inputs. If not given, default input file is used.

    • Manning’s n is required.

    • Storm duration is needed if rainfall_duration is not passed in the initialization

    • Rainfall intensity is needed if rainfall_intensity is not passed in the initialization

    • Model run time can be provided in initialization. If not it is set to the storm duration

  • h_init (float, optional) – Some initial depth in the channels. Default = 0.001 m

  • g (float, optional) – Gravitational acceleration, \(m / s^2\)

  • alpha (float, optional) – Non-dimensional time step factor from Bates et al., (2010)

  • rho (integer, optional) – Density of water, \(kg / m^3\)

  • ten_thirds (float, optional) – Precalculated value of \(10 / 3\) which is used in the implicit shallow water equation.

Examples

>>> DEM_name = "DEM_name.asc"
>>> (rg, z) = read_esri_ascii(DEM_name)  
>>> of = OverlandFlowBates(rg)  

References

Required Software Citation(s) Specific to this Component

None Listed

Additional References

Bates, P., Horritt, M., Fewtrell, T. (2010). A simple inertial formulation of the shallow water equations for efficient two-dimensional flood inundation modelling Journal of Hydrology 387(1-2), 33-45. https://dx.doi.org/10.1016/j.jhydrol.2010.03.027

calc_time_step()[source]#
property dt#

component timestep.

Type:

dt

property h#

The depth of water at each node.

overland_flow(dt=None)[source]#

For one time step, this generates ‘overland flow’ across a given grid by calculating discharge at each node.

Using the depth slope product, shear stress is calculated at every node.

Outputs water depth, discharge and shear stress values through time at every point in the input grid.

Parameters:
  • grid (RasterModelGrid) – A grid.

  • dt (float, optional) – Time step. Either set when called or the component will do it for you.

property surface_water__discharge#

The discharge of water on active links.

KinematicWaveRengers#

class KinematicWaveRengers(*args, **kwds)[source]#

Bases: Component

This code is based on an overland flow model by Francis Rengers and colleagues, after Julien et al., 1995. It uses an explicit face-centered solution to a depth-varying Manning’s equation, broadly following, e.g., Mugler et al., 2011. It was implemented in Landlab by DEJH, March ‘16. Please cite Rengers et al., 2016, Model Predictions of Water Runoff in Steep Catchments after Wildfire, WRR.

Note: You will not have a good day if you have pits present in your topo before routing flow with this component. Fill pits before instantiating the component (or call update_topographic_params once you have filled after instantiation).

Note this module assumes that the topography DOES NOT change during the run. If it does, call update_topographic_params to update the component to the new topo.

Boundary condition control can be… interesting with this component. Be sure to close boundaries you do not wish water to leave - or enter! - through. To allow free water discharge from the grid edge it is recommended to use fixed gradient boundary conditions at the open edges. The component will then set the fixed gradient as equal to the underlying topographic gradient throughout the run.

It is also possible to fix the water depth at the open edge, but this is not really recommended.

Construction:

KinematicWaveRengers(grid, mannings_n=0.03, critical_flow_depth=0.003,
                     mannings_epsilon=0.33333333, dt_max=0.3,
                     max_courant=0.2, min_surface_water_depth=1.e-8)
Parameters:
  • grid (RasterModelGrid) – A grid.

  • mannings_n (float) – A value to use for Manning’s n in the Manning discharge equation.

  • critical_flow_depth (float (m)) – An index flow depth for the depth-varying Manning’s equation, controlling the depth at which the effective Manning’s n begins to increase.

  • mannings_epsilon (float) – An exponent for the depth-varying Manning’s equation, controlling the rate of increase of effective Manning’s n at small flow depths.

  • dt_max (float or None (s)) – The largest permitted internal timestep for the component. If the Courant criterion produces a more restrictive condition, that will be used instead.

  • max_courant (float) – The maximum permitted Courant number for the courant stability criterion.

  • min_surface_water_depth (float (m)) – A water depth below which surface water thickness may never fall, to ensure model stabilty.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import KinematicWaveRengers
>>> mg = RasterModelGrid((5, 10), 10.0)
>>> mg.status_at_node[mg.nodes_at_left_edge] = mg.BC_NODE_IS_FIXED_GRADIENT
>>> mg.status_at_node[mg.nodes_at_top_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.status_at_node[mg.nodes_at_bottom_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.status_at_node[mg.nodes_at_right_edge] = mg.BC_NODE_IS_CLOSED
>>> _ = mg.add_field("node", "topographic__elevation", 0.05 * mg.node_x)
>>> _ = mg.add_empty("node", "surface_water__depth")
>>> mg.at_node["surface_water__depth"].fill(1.0e-8)
>>> dt = 60.0  # 1 min intervals
>>> rain_intensities = (1.0e-5, 1.0e-5, 1.0e-5, 1.0e-5, 1.0e-5)
>>> kw = KinematicWaveRengers(mg)
>>> for i in rain_intensities:
...     kw.run_one_step(dt, rainfall_intensity=i)
...
>>> mg.at_node["surface_water__depth"]
array([  1.00000000e-08,   1.00000000e-08,   1.00000000e-08,
         1.00000000e-08,   1.00000000e-08,   1.00000000e-08,
         1.00000000e-08,   1.00000000e-08,   1.00000000e-08,
         1.00000000e-08,   2.95578314e-03,   2.95578314e-03,
         2.90945761e-03,   2.82912876e-03,   2.70127141e-03,
         2.51202011e-03,   2.24617193e-03,   1.88032853e-03,
         1.35451064e-03,   1.00000000e-08,   2.95578314e-03,
         2.95578314e-03,   2.90945761e-03,   2.82912876e-03,
         2.70127141e-03,   2.51202011e-03,   2.24617193e-03,
         1.88032853e-03,   1.35451064e-03,   1.00000000e-08,
         2.95578314e-03,   2.95578314e-03,   2.90945761e-03,
         2.82912876e-03,   2.70127141e-03,   2.51202011e-03,
         2.24617193e-03,   1.88032853e-03,   1.35451064e-03,
         1.00000000e-08,   1.00000000e-08,   1.00000000e-08,
         1.00000000e-08,   1.00000000e-08,   1.00000000e-08,
         1.00000000e-08,   1.00000000e-08,   1.00000000e-08,
         1.00000000e-08,   1.00000000e-08])

Initialize the kinematic wave approximation overland flow component.

__init__(grid, mannings_n=0.03, critical_flow_depth=0.003, mannings_epsilon=0.33333333, dt_max=0.3, max_courant=0.2, min_surface_water_depth=1e-08)[source]#

Initialize the kinematic wave approximation overland flow component.

calc_grads_and_timesteps(update_topography, track_min_depth)[source]#

Perform the first part of the calculation for the main run, mainly velocities and fluxes. The main objective of this part of the calculation is to derive the stable timestep for the run.

Parameters:
  • update_topography (bool) – If False, the underlying surface topography is assumed unchanged since the last run.

  • track_min_depth (bool) – If True, the internal list _water_balance will be appended with the volumetric fractional change in mass balance during the run.

Returns:

internal_dt – The internally stable timestep that will be used on this loop.

Return type:

float

property internal_timestep#

Return the internal timestep last used by the kinematic wave component.

run_one_step(dt, rainfall_intensity=1e-05, update_topography=False, track_min_depth=False)[source]#

Update fields with current hydrologic conditions.

Parameters:
  • rain_intensity (float or array (m/s)) – The rainfall intensity across the grid (water input rate at each node).

  • update_topography (bool) – Set to true if the topography of the grid evolves during the run.

  • track_min_depth (bool) – At very low rainfall inputs, there is a possibility this component could allow creation of small amounts of water mass. Set to true to track this mass, and use the water_balance property to investigate its evolution through time.

update_topographic_params()[source]#

If the topo changes during the run, change the held params used by run_one_step.

property water_balance#

Return a list of the fractional gain/loss of water mass during the run, if it was tracked using the track_min_depth flag.

KinwaveImplicitOverlandFlow#

Landlab component for overland flow using a local implicit solution to the kinematic-wave approximation.

Created on Fri May 27 14:26:13 2016

@author: gtucker

class KinwaveImplicitOverlandFlow(*args, **kwds)[source]#

Bases: Component

Calculate shallow water flow over topography.

Landlab component that implements a two-dimensional kinematic wave model. This is a form of the 2D shallow-water equations in which energy slope is assumed to equal bed slope. The solution method is locally implicit, and works as follows. At each time step, we iterate from upstream to downstream over the topography. Because we are working downstream, we can assume that we know the total water inflow to a given cell. We solve the following mass conservation equation at each cell:

\[(H^{t+1} - H^t)/\Delta t = Q_{in}/A - Q_{out}/A + R\]

where \(H\) is water depth, \(t\) indicates time step number, \(\Delta t\) is time step duration, \(Q_{in}\) is total inflow discharge, \(Q_{out}\) is total outflow discharge, \(A\) is cell area, and \(R\) is local runoff rate (precipitation minus infiltration; could be negative if runon infiltration is occurring).

The specific outflow discharge leaving a cell along one of its faces is:

\[q = (1/C_r) H^\alpha S^{1/2}\]

where \(C_r\) is a roughness coefficient (such as Manning’s n), \(\alpha\) is an exponent equal to \(5/3\) for the Manning equation and \(3/2\) for the Chezy family, and \(S\) is the downhill-positive gradient of the link that crosses this particular face. Outflow discharge is zero for links that are flat or “uphill” from the given node. Total discharge out of a cell is then the sum of (specific discharge x face width) over all outflow faces

\[Q_{out} = \sum_{i=1}^N (1/C_r) H^\alpha S_i^{1/2} W_i\]

where \(N\) is the number of outflow faces (i.e., faces where the ground slopes downhill away from the cell’s node), and \(W_i\) is the width of face \(i\).

We use the depth at the cell’s node, so this simplifies to:

\[Q_{out} = (1/C_r) H'^\alpha \sum_{i=1}^N S_i^{1/2} W_i\]

We define \(H\) in the above as a weighted sum of the “old” (time step \(t\)) and “new” (time step \(t+1\)) depth values:

\[H' = w H^{t+1} + (1-w) H^t\]

If \(w=1\), the method is fully implicit. If \(w=0\), it is a simple forward explicit method.

When we combine these equations, we have an equation that includes the unknown \(H^{t+1}\) and a bunch of terms that are known. If \(w\ne 0\), it is a nonlinear equation in \(H^{t+1}\), and must be solved iteratively. We do this using a root-finding method in the scipy.optimize library.

Examples

>>> from landlab import RasterModelGrid
>>> rg = RasterModelGrid((4, 5), xy_spacing=10.0)
>>> z = rg.add_zeros("topographic__elevation", at="node")
>>> kw = KinwaveImplicitOverlandFlow(rg)
>>> round(kw.runoff_rate * 1.0e7, 2)
2.78
>>> kw.vel_coef  # default value
100.0
>>> rg.at_node["surface_water__depth"][6:9]
array([ 0.,  0.,  0.])

References

Required Software Citation(s) Specific to this Component

None Listed

Additional References

None Listed

Initialize the KinwaveImplicitOverlandFlow.

Parameters:
  • grid (ModelGrid) – Landlab ModelGrid object

  • runoff_rate (float, optional (defaults to 1 mm/hr)) – Precipitation rate, mm/hr. The value provided is divided by 3600000.0.

  • roughness (float, defaults to 0.01) – Manning roughness coefficient; units depend on depth_exp.

  • changing_topo (boolean, optional (defaults to False)) – Flag indicating whether topography changes between time steps

  • depth_exp (float (defaults to 1.5)) – Exponent on water depth in velocity equation (3/2 for Darcy/Chezy, 5/3 for Manning)

  • weight (float (defaults to 1.0)) – Weighting on depth at new time step versus old time step (1 = all implicit; 0 = explicit)

__init__(grid, runoff_rate=1.0, roughness=0.01, changing_topo=False, depth_exp=1.5, weight=1.0)[source]#

Initialize the KinwaveImplicitOverlandFlow.

Parameters:
  • grid (ModelGrid) – Landlab ModelGrid object

  • runoff_rate (float, optional (defaults to 1 mm/hr)) – Precipitation rate, mm/hr. The value provided is divided by 3600000.0.

  • roughness (float, defaults to 0.01) – Manning roughness coefficient; units depend on depth_exp.

  • changing_topo (boolean, optional (defaults to False)) – Flag indicating whether topography changes between time steps

  • depth_exp (float (defaults to 1.5)) – Exponent on water depth in velocity equation (3/2 for Darcy/Chezy, 5/3 for Manning)

  • weight (float (defaults to 1.0)) – Weighting on depth at new time step versus old time step (1 = all implicit; 0 = explicit)

property depth#

The depth of water at each node.

run_one_step(dt)[source]#

Calculate water flow for a time period dt.

property runoff_rate#

Runoff rate.

Parameters:

runoff_rate (float, optional (defaults to 1 mm/hr)) – Precipitation rate, mm/hr. The value provide is divided by 3600000.0.

Return type:

The current value of the runoff rate.

property vel_coef#

Velocity coefficient.

water_fn(x, a, b, c, d, e)[source]#

Evaluates the solution to the water-depth equation.

Called by scipy.newton() to find solution for \(x\) using Newton’s method.

Parameters:
  • x (float) – Water depth at new time step.

  • a (float) – “alpha” parameter (see below)

  • b (float) – Weighting factor on new versus old time step. \(b=1\) means purely implicit solution with all weight on \(H\) at new time step. \(b=0\) (not recommended) would mean purely explicit.

  • c (float) – Water depth at old time step (time step \(t\) instead of \(t+1\))

  • d (float) – Depth-discharge exponent; normally either 5/3 (Manning) or 3/2 (Chezy)

  • e (float) – Water inflow volume per unit cell area in one time step.

This equation represents the implicit solution for water depth \(H\) at the next time step. In the code below, it is formulated in a generic way. Written using more familiar terminology, the equation is:

\[H - H_0 + \alpha ( w H + (w-1) H_0)^d - \Delta t (R + Q_{in} / A)\]
\[\alpha = \frac{\Delta t \sum S^{1/2}}{C_f A}\]

where \(H\) is water depth at the given node at the new time step, \(H_0\) is water depth at the prior time step, \(w\) is a weighting factor, \(d\) is the depth-discharge exponent (2/3 or 1/2), \(\Delta t\) is time-step duration, \(R\) is local runoff rate, \(Q_{in}\) is inflow discharge, \(A\) is cell area, \(C_f\) is a dimensional roughness coefficient, and \(\sum S^{1/2}\) represents the sum of square-root-of-downhill-gradient over all outgoing (downhill) links.

KinwaveOverlandFlowModel#

Landlab component for overland flow using the kinematic-wave approximation.

Created on Fri May 27 14:26:13 2016

@author: gtucker

class KinwaveOverlandFlowModel(*args, **kwds)[source]#

Bases: Component

Calculate water flow over topography.

Landlab component that implements a two-dimensional kinematic wave model. This is an extremely simple, unsophisticated model, originally built simply to demonstrate the component creation process. Limitations to the present version include: infiltration is handled very crudely, the called is responsible for picking a stable time step size (no adaptive time stepping is used in the run_one_step method), precipitation rate is constant for a given duration (then zero), and all parameters are uniform in space. Also, the terrain is assumed to be stable over time. Caveat emptor!

Examples

>>> from landlab import RasterModelGrid
>>> rg = RasterModelGrid((4, 5), xy_spacing=10.0)
>>> z = rg.add_zeros("topographic__elevation", at="node")
>>> s = rg.add_zeros("topographic__gradient", at="link")
>>> kw = KinwaveOverlandFlowModel(rg)
>>> kw.vel_coef
100.0
>>> rg.at_node["surface_water__depth"]
array([ 0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,
        0.,  0.,  0.,  0.,  0.,  0.,  0.])

References

Required Software Citation(s) Specific to this Component

None Listed

Additional References

None Listed

Initialize the KinwaveOverlandFlowModel.

Parameters:
  • grid (ModelGrid) – Landlab ModelGrid object

  • precip_rate (float, optional (defaults to 1 mm/hr)) – Precipitation rate, mm/hr

  • precip_duration (float, optional (defaults to 1 hour)) – Duration of precipitation, hours

  • infilt_rate (float, optional (defaults to 0)) – Maximum rate of infiltration, mm/hr

  • roughness (float, defaults to 0.01) – Manning roughness coefficient, s/m^1/3

__init__(grid, precip_rate=1.0, precip_duration=1.0, infilt_rate=0.0, roughness=0.01)[source]#

Initialize the KinwaveOverlandFlowModel.

Parameters:
  • grid (ModelGrid) – Landlab ModelGrid object

  • precip_rate (float, optional (defaults to 1 mm/hr)) – Precipitation rate, mm/hr

  • precip_duration (float, optional (defaults to 1 hour)) – Duration of precipitation, hours

  • infilt_rate (float, optional (defaults to 0)) – Maximum rate of infiltration, mm/hr

  • roughness (float, defaults to 0.01) – Manning roughness coefficient, s/m^1/3

run_one_step(dt)[source]#

Calculate water flow for a time period dt.

Default units for dt are seconds.

property vel_coef#

Velocity coefficient.

(1/roughness)

LinearDiffusionOverlandFlowRouter#

Landlab component for overland flow using the linearized diffusion-wave approximation.

Created on Fri May 27 14:26:13 2016

@author: gtucker

class LinearDiffusionOverlandFlowRouter(*args, **kwds)[source]#

Bases: Component

Calculate water flow over topography.

Landlab component that implements a two-dimensional, linearized diffusion-wave model. The diffusion-wave approximation is a simplification of the shallow-water equations that omits the momentum terms. The flow velocity is calculated using the local water-surface slope as an approximation of the energy slope. With this linearized form, flow velocity is calculated using a linearized Manning equation, with the water-surface slope being used as the slope factor. There are two governing equations, one that represents conservation of water mass:

..math:

\frac{\partial H}{\partial t} = (P - I) - \nabla\cdot\mathbf{q}

where \(H(x,y,t)\) is local water depth, \(t\) is time, \(P\) is precipitation rate, \(I\) is infiltration rate, and \(\mathbf{q}\) is specific water discharge, which equals depth times depth-averaged velocity. The other governing equation represents specific discharge as a function of gravity, pressure, and friction:

..math:

\mathbf{q} = \frac{H^{4/3}}{n^2 U_c} \nabla w

where \(n\) is the friction factor (“Manning’s n”), \(U_c\) is a characteristic flow velocity, and \(w\) is water-surface height, which is the sum of topographic elevation plus water depth.

Infiltration rate should decline smoothly to zero as surface water depth approaches zero. To ensure this, infiltration rate is calculated as

..math:

I = I_c \left( 1 - e^{-H / H_i} ) \right)

where \(H_i\) is a characteristic water depth. The concept here is that when \(H \le H_i\), small spatial variations in water depth will leave parts of the ground un-ponded and therefore not subject to any infiltration. Mathematically, \(I \approx 0.95 I_c\) when \(H = 3H_i\).

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 3))
>>> _ = grid.add_zeros("topographic__elevation", at="node")
>>> olf = LinearDiffusionOverlandFlowRouter(grid, roughness=0.1)
>>> round(olf.vel_coef)
100
>>> olf.rain_rate
1e-05
>>> olf.rain_rate = 1.0e-4
>>> olf.run_one_step(dt=10.0)
>>> grid.at_node["surface_water__depth"][4]
0.001

References

Required Software Citation(s) Specific to this Component

None Listed

Additional References

None Listed

Initialize the LinearDiffusionOverlandFlowRouter.

Parameters:
  • grid (ModelGrid) – Landlab ModelGrid object

  • roughness (float, defaults to 0.01) – Manning roughness coefficient, s/m^1/3

  • rain_rate (float, optional (defaults to 36 mm/hr)) – Rainfall rate, m/s

  • infilt_depth_scale (float, defaults to 0.001) – Depth scale for water infiltration, m

  • infilt_rate (float, optional (defaults to 0)) – Maximum rate of infiltration, m/s

  • velocity_scale (float, defaults to 1) – Characteristic flow velocity, m/s

  • cfl_factor (float, optional (defaults to 0.2)) – Time-step control factor: fraction of maximum estimated time-step that is actually used (must be <=1)

__init__(grid, roughness=0.01, rain_rate=1e-05, infilt_rate=0.0, infilt_depth_scale=0.001, velocity_scale=1.0, cfl_factor=0.2)[source]#

Initialize the LinearDiffusionOverlandFlowRouter.

Parameters:
  • grid (ModelGrid) – Landlab ModelGrid object

  • roughness (float, defaults to 0.01) – Manning roughness coefficient, s/m^1/3

  • rain_rate (float, optional (defaults to 36 mm/hr)) – Rainfall rate, m/s

  • infilt_depth_scale (float, defaults to 0.001) – Depth scale for water infiltration, m

  • infilt_rate (float, optional (defaults to 0)) – Maximum rate of infiltration, m/s

  • velocity_scale (float, defaults to 1) – Characteristic flow velocity, m/s

  • cfl_factor (float, optional (defaults to 0.2)) – Time-step control factor: fraction of maximum estimated time-step that is actually used (must be <=1)

property rain_rate#

Rainfall rate

run_one_step(dt)[source]#

Calculate water flow for a time period dt.

Default units for dt are seconds. We use a time-step subdivision algorithm that ensures step size is always below CFL limit.

update_for_one_iteration(iter_dt)[source]#

Update state variables for one iteration of duration iter_dt.

property vel_coef#

Velocity coefficient.

(1/(roughness^2 x velocity_scale)