landlab.components.landslides.landslide_probability

Landlab component that simulates landslide probability of failure as well as mean relative wetness and probability of saturation.

Relative wetness and factor-of-safety are based on the infinite slope stability model driven by topographic and soils inputs and recharge provided by user as inputs to the component. For each node, component simulates mean relative wetness as well as the probability of saturation based on Monte Carlo simulation of relative wetness where the probability is the number of iterations with relative wetness >= 1.0 divided by the number of iterations. Probability of failure for each node is also simulated in the Monte Carlo simulation as the number of iterations with factor-of-safety <= 1.0 divided by the number of iterations.

Code author: R.Strauch, E.Istanbulluoglu, & S.S.Nudurupati

University of Washington

Ref 1: Strauch et. al. 2017, ‘A hydro-climatological approach to predicting regional landslide probability using Landlab, Earth Surface Dynamics, In prep.

Ref 2: ‘The Landlab LandslideProbability Component User Manual’ @ https://github.com/RondaStrauch/pub_strauch_etal_esurf/blob/master/LandslideComponentUsersManual.pdf

Created on Thu Aug 20, 2015 Last edit June 7, 2017

class LandslideProbability[source]

Bases: Component

Landslide probability component using the infinite slope stability model.

Landlab component designed to calculate probability of failure at each grid node based on the infinite slope stability model stability index (Factor of Safety).

The driving force for failure is provided by the user in the form of groundwater recharge; four options for providing recharge are supported. The model uses topographic and soil characteristics provided as input by the user.

The main method of the LandslideProbability class is calculate_landslide_probability()`, which calculates the mean soil relative wetness, probability of soil saturation, and probability of failure at each node based on a Monte Carlo simulation.

Usage:

Option 1 - Uniform recharge

LandslideProbability(
    grid,
    number_of_iterations=250,
    groundwater__recharge_distribution="uniform",
    groundwater__recharge_min_value=5.0,
    groundwater__recharge_max_value=121.0,
)

Option 2 - Lognormal recharge

LandslideProbability(
    grid,
    number_of_iterations=250,
    groundwater__recharge_distribution="lognormal",
    groundwater__recharge_mean=30.0,
    groundwater__recharge_standard_deviation=0.25,
)

Option 3 - Lognormal_spatial recharge

LandslideProbability(
    grid,
    number_of_iterations=250,
    groundwater__recharge_distribution="lognormal_spatial",
    groundwater__recharge_mean=np.random.randint(20, 120, grid_size),
    groundwater__recharge_standard_deviation=np.random.rand(grid_size),
)

Option 4 - Data_driven_spatial recharge

LandslideProbability(
    grid,
    number_of_iterations=250,
    groundwater__recharge_distribution="data_driven_spatial",
    groundwater__recharge_HSD_inputs=[HSD_dict, HSD_id_dict, fract_dict],
)

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components.landslides import LandslideProbability
>>> import numpy as np

Create a grid on which to calculate landslide probability.

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

Check the number of core nodes.

>>> grid.number_of_core_nodes
6

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

>>> sorted(LandslideProbability.input_var_names)
['soil__density',
 'soil__internal_friction_angle',
 'soil__maximum_total_cohesion',
 'soil__minimum_total_cohesion',
 'soil__mode_total_cohesion',
 'soil__saturated_hydraulic_conductivity',
 'soil__thickness',
 'soil__transmissivity',
 'topographic__slope',
 'topographic__specific_contributing_area']

Check the units for the fields.

>>> LandslideProbability.var_units("topographic__specific_contributing_area")
'm'

Create an input field.

>>> grid.at_node["topographic__slope"] = np.random.rand(grid.number_of_nodes)

If you are not sure about one of the input or output variables, you can get help for specific variables.

>>> LandslideProbability.var_help("soil__transmissivity")
name: soil__transmissivity
description:
  mode rate of water transmitted through a unit width of saturated
  soil - either provided or calculated with Ksat and soil depth
units: m2/day
unit agnostic: False
at: node
intent: in

Additional required fields for component.

>>> scatter_dat = np.random.randint(1, 10, grid.number_of_nodes)
>>> grid.at_node["topographic__specific_contributing_area"] = np.sort(
...     np.random.randint(30, 900, grid.number_of_nodes).astype(float)
... )
>>> grid.at_node["soil__transmissivity"] = np.sort(
...     np.random.randint(5, 20, grid.number_of_nodes).astype(float), -1
... )
>>> grid.at_node["soil__saturated_hydraulic_conductivity"] = np.sort(
...     np.random.randint(2, 10, grid.number_of_nodes).astype(float), -1
... )
>>> grid.at_node["soil__mode_total_cohesion"] = np.sort(
...     np.random.randint(30, 900, grid.number_of_nodes).astype(float)
... )
>>> grid.at_node["soil__minimum_total_cohesion"] = (
...     grid.at_node["soil__mode_total_cohesion"] - scatter_dat
... )
>>> grid.at_node["soil__maximum_total_cohesion"] = (
...     grid.at_node["soil__mode_total_cohesion"] + scatter_dat
... )
>>> grid.at_node["soil__internal_friction_angle"] = np.sort(
...     np.random.randint(26, 40, grid.number_of_nodes).astype(float)
... )
>>> grid.at_node["soil__thickness"] = np.sort(
...     np.random.randint(1, 10, grid.number_of_nodes).astype(float)
... )
>>> grid.at_node["soil__density"] = 2000.0 * np.ones(grid.number_of_nodes)

Instantiate the ‘LandslideProbability’ component to work on this grid, and run it.

>>> ls_prob = LandslideProbability(grid)
>>> np.allclose(grid.at_node["landslide__probability_of_failure"], 0.0)
True

Run the calculate_landslide_probability method to update output variables with grid

>>> ls_prob.calculate_landslide_probability()

Check the output variable names.

>>> sorted(ls_prob.output_var_names)
['landslide__probability_of_failure',
 'soil__mean_relative_wetness',
 'soil__probability_of_saturation']

Check the output from the component, including array at one node.

>>> np.allclose(grid.at_node["landslide__probability_of_failure"], 0.0)
False
>>> core_nodes = ls_prob.grid.core_nodes

References

Required Software Citation(s) Specific to this Component

Strauch, R., Istanbulluoglu, E., Nudurupati, S., Bandaragoda, C., Gasparini, N., Tucker, G. (2018). A hydroclimatological approach to predicting regional landslide probability using Landlab Earth Surface Dynamics 6(1), 49-75. https://dx.doi.org/10.5194/esurf-6-49-2018

Additional References

None Listed

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

  • number_of_iterations (int, optional) – Number of iterations to run Monte Carlo simulation (default=250).

  • groundwater__recharge_distribution (str, optional) – single word indicating recharge distribution, either ‘uniform’, ‘lognormal’, ‘lognormal_spatial,’ or ‘data_driven_spatial’. (default=’uniform’)

  • groundwater__recharge_min_value (float, optional (mm/d)) – minium groundwater recharge for ‘uniform’ (default=20.)

  • groundwater__recharge_max_value (float, optional (mm/d)) – maximum groundwater recharge for ‘uniform’ (default=120.)

  • groundwater__recharge_mean (float, optional (mm/d)) – mean grounwater recharge for ‘lognormal’ and ‘lognormal_spatial’ (default=None)

  • groundwater__recharge_standard_deviation (float, optional (mm/d)) – standard deviation of grounwater recharge for ‘lognormal’ and ‘lognormal_spatial’ (default=None)

  • groundwater__recharge_HSD_inputs (list, optional) – list of 3 dictionaries in order (default=[]) - HSD_dict {Hydrologic Source Domain (HSD) keys: recharge numpy array values}, {node IDs keys: list of HSD_Id values}, HSD_fractions {node IDS keys: list of HSD fractions values} (none) Note: this input method is a very specific one, and to use this method, one has to refer Ref 1 & Ref 2 mentioned above, as this set of inputs require rigorous pre-processing of data.

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

  • seed (int, optional) – seed for random number generation. if seed is assigned any value other than the default value of zero, it will create different sequence. To create a certain sequence repititively, use the same value as input for seed.

__init__(grid, number_of_iterations=250, g=scipy.constants.g, groundwater__recharge_distribution='uniform', groundwater__recharge_min_value=20.0, groundwater__recharge_max_value=120.0, groundwater__recharge_mean=None, groundwater__recharge_standard_deviation=None, groundwater__recharge_HSD_inputs=(), seed=0)[source]
Parameters:
  • grid (RasterModelGrid) – A raster grid.

  • number_of_iterations (int, optional) – Number of iterations to run Monte Carlo simulation (default=250).

  • groundwater__recharge_distribution (str, optional) – single word indicating recharge distribution, either ‘uniform’, ‘lognormal’, ‘lognormal_spatial,’ or ‘data_driven_spatial’. (default=’uniform’)

  • groundwater__recharge_min_value (float, optional (mm/d)) – minium groundwater recharge for ‘uniform’ (default=20.)

  • groundwater__recharge_max_value (float, optional (mm/d)) – maximum groundwater recharge for ‘uniform’ (default=120.)

  • groundwater__recharge_mean (float, optional (mm/d)) – mean grounwater recharge for ‘lognormal’ and ‘lognormal_spatial’ (default=None)

  • groundwater__recharge_standard_deviation (float, optional (mm/d)) – standard deviation of grounwater recharge for ‘lognormal’ and ‘lognormal_spatial’ (default=None)

  • groundwater__recharge_HSD_inputs (list, optional) – list of 3 dictionaries in order (default=[]) - HSD_dict {Hydrologic Source Domain (HSD) keys: recharge numpy array values}, {node IDs keys: list of HSD_Id values}, HSD_fractions {node IDS keys: list of HSD fractions values} (none) Note: this input method is a very specific one, and to use this method, one has to refer Ref 1 & Ref 2 mentioned above, as this set of inputs require rigorous pre-processing of data.

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

  • seed (int, optional) – seed for random number generation. if seed is assigned any value other than the default value of zero, it will create different sequence. To create a certain sequence repititively, use the same value as input for seed.

static __new__(cls, *args, **kwds)
calculate_factor_of_safety(i)[source]

Method to calculate factor of safety.

Method calculates factor-of-safety stability index by using node specific parameters, creating distributions of these parameters, and calculating the index by sampling these distributions ‘n’ times.

The index is calculated from the ‘infinite slope stabilty factor-of-safety equation’ in the format of Pack RT, Tarboton DG, and Goodwin CN (1998),The SINMAP approach to terrain stability mapping.

Parameters:

i (int) – index of core node ID.

calculate_landslide_probability()[source]

Main method of Landslide Probability class.

Method creates arrays for output variables then loops through all the core nodes to run the method ‘calculate_factor_of_safety.’ Output parameters probability of failure, mean relative wetness, and probability of saturation are assigned as fields to nodes.

cite_as = '\n    @article{strauch2018hydroclimatological,\n      author = {Strauch, Ronda and Istanbulluoglu, Erkan and Nudurupati,\n      Sai Siddhartha and Bandaragoda, Christina and Gasparini, Nicole M and\n      Tucker, Gregory E},\n      title = {{A hydroclimatological approach to predicting regional landslide\n      probability using Landlab}},\n      issn = {2196-6311},\n      doi = {10.5194/esurf-6-49-2018},\n      pages = {49--75},\n      number = {1},\n      volume = {6},\n      journal = {Earth Surface Dynamics},\n      year = {2018}\n    }\n    '
property coords

Return the coordinates of nodes on grid attached to the component.

property current_time

Current time.

Some components may keep track of the current time. In this case, the current_time attribute is incremented. Otherwise it is set to None.

Return type:

current_time

definitions = (('landslide__probability_of_failure', 'number of times FS is <=1 out of number of iterations user selected'), ('soil__density', 'wet bulk density of soil'), ('soil__internal_friction_angle', 'critical angle just before failure due to friction between particles'), ('soil__maximum_total_cohesion', 'maximum of combined root and soil cohesion at node'), ('soil__mean_relative_wetness', 'Indicator of soil wetness; relative depth perched water table within the soil layer'), ('soil__minimum_total_cohesion', 'minimum of combined root and soil cohesion at node'), ('soil__mode_total_cohesion', 'mode of combined root and soil cohesion at node'), ('soil__probability_of_saturation', 'number of times relative wetness is >=1 out of number of iterations user selected'), ('soil__saturated_hydraulic_conductivity', 'mode rate of water transmitted through soil - provided if transmissivity is NOT provided to calculate tranmissivity with soil depth'), ('soil__thickness', 'soil depth to restrictive layer'), ('soil__transmissivity', 'mode rate of water transmitted through a unit width of saturated soil - either provided or calculated with Ksat and soil depth'), ('topographic__slope', 'gradient of the ground surface'), ('topographic__specific_contributing_area', 'specific contributing (upslope area/cell face ) that drains to node'))
classmethod from_path(grid, path)

Create a component from an input file.

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

  • path (str or file_like) – Path to a parameter file, contents of a parameter file, or a file-like object.

Returns:

A newly-created component.

Return type:

Component

property grid

Return the grid attached to the component.

initialize_optional_output_fields()

Create fields for a component based on its optional field outputs, if declared in _optional_var_names.

This method will create new fields (without overwrite) for any fields output by the component as optional. New fields are initialized to zero. New fields are created as arrays of floats, unless the component also contains the specifying property _var_type.

initialize_output_fields(values_per_element=None)

Create fields for a component based on its input and output var names.

This method will create new fields (without overwrite) for any fields output by, but not supplied to, the component. New fields are initialized to zero. Ignores optional fields. New fields are created as arrays of floats, unless the component specifies the variable type.

Parameters:

values_per_element (int (optional)) – On occasion, it is necessary to create a field that is of size (n_grid_elements, values_per_element) instead of the default size (n_grid_elements,). Use this keyword argument to acomplish this task.

input_var_names = ('soil__density', 'soil__internal_friction_angle', 'soil__maximum_total_cohesion', 'soil__minimum_total_cohesion', 'soil__mode_total_cohesion', 'soil__saturated_hydraulic_conductivity', 'soil__thickness', 'soil__transmissivity', 'topographic__slope', 'topographic__specific_contributing_area')
name = 'Landslide Probability'
optional_var_names = ()
output_var_names = ('landslide__probability_of_failure', 'soil__mean_relative_wetness', 'soil__probability_of_saturation')
property shape

Return the grid shape attached to the component, if defined.

unit_agnostic = False
units = (('landslide__probability_of_failure', 'None'), ('soil__density', 'kg/m3'), ('soil__internal_friction_angle', 'degrees'), ('soil__maximum_total_cohesion', 'Pa or kg/m-s2'), ('soil__mean_relative_wetness', 'None'), ('soil__minimum_total_cohesion', 'Pa or kg/m-s2'), ('soil__mode_total_cohesion', 'Pa or kg/m-s2'), ('soil__probability_of_saturation', 'None'), ('soil__saturated_hydraulic_conductivity', 'm/day'), ('soil__thickness', 'm'), ('soil__transmissivity', 'm2/day'), ('topographic__slope', 'tan theta'), ('topographic__specific_contributing_area', 'm'))
classmethod var_definition(name)

Get a description of a particular field.

Parameters:

name (str) – A field name.

Returns:

A description of each field.

Return type:

tuple of (name, *description*)

classmethod var_help(name)

Print a help message for a particular field.

Parameters:

name (str) – A field name.

classmethod var_loc(name)

Location where a particular variable is defined.

Parameters:

name (str) – A field name.

Returns:

The location (‘node’, ‘link’, etc.) where a variable is defined.

Return type:

str

var_mapping = (('landslide__probability_of_failure', 'node'), ('soil__density', 'node'), ('soil__internal_friction_angle', 'node'), ('soil__maximum_total_cohesion', 'node'), ('soil__mean_relative_wetness', 'node'), ('soil__minimum_total_cohesion', 'node'), ('soil__mode_total_cohesion', 'node'), ('soil__probability_of_saturation', 'node'), ('soil__saturated_hydraulic_conductivity', 'node'), ('soil__thickness', 'node'), ('soil__transmissivity', 'node'), ('topographic__slope', 'node'), ('topographic__specific_contributing_area', 'node'))
classmethod var_type(name)

Returns the dtype of a field (float, int, bool, str…).

Parameters:

name (str) – A field name.

Returns:

The dtype of the field.

Return type:

dtype

classmethod var_units(name)

Get the units of a particular field.

Parameters:

name (str) – A field name.

Returns:

Units for the given field.

Return type:

str