landlab.components.lithology.litholayers

Create a LithoLayers component with different properties.

class LithoLayers[source]

Bases: Lithology

Create LithoLayers component.

A LithoLayers is a three dimentional representation of material operated on by landlab components. Material can be removed through erosion or added to through deposition. Rock types can have multiple attributes (e.g. age, erodability or other parameter values, etc).

If the tracked properties are model grid fields, they will be updated to the surface values of the Lithology. If the properties are not grid fields then at-node grid fields will be created with their names.

It is constructed by specifying a series of depths below the surface, an anchor point, a series of rock type ids, and the functional form of a surface. Depths and IDs are both specified in order of closest to the surface to furthest from the surface.

Additionally, an attribute dictionary specifies the properties of each rock type. This dictionary is expected to have the form of:

attrs = {"K_sp": {1: 0.001, 2: 0.0001}, "D": {1: 0.01, 2: 0.001}}

Where 'K_sp' and 'D' are properties to track, and 1 and 2 are rock type IDs. The rock type IDs can be any type that is valid as a python dictionary key.

References

Required Software Citation(s) Specific to this Component

Barnhart, K., Hutton, E., Gasparini, N., Tucker, G. (2018). Lithology: A Landlab submodule for spatially variable rock properties. Journal of Open Source Software 3(30), 979 - 2. https://dx.doi.org/10.21105/joss.00979

Additional References

None Listed

Create a new instance of a LithoLayers.

Parameters:
  • grid (Landlab ModelGrid)

  • z0s (ndarray of shape (n_layers, )) – Values of layer depth from surface at horizontal location (x0, y0).

  • ids (ndarray of shape (n_layers, )) – Values of rock type IDs corresponding to each layer specified in z0s.

  • attrs (dict) – Rock type property dictionary. See class docstring for example of required format.

  • x0 (float, optional) – x value of anchor point for all layers.

  • y0 (float, optional) – y value of anchor point for all layers.

  • function (function, optional) – Functional form of layers as a function of two variables, x and y. Default value is lambda x, y: 0*x + 0*y for flatlying layers.

  • layer_type (str, optional) – Type of Landlab layers object used to store the layers. If MaterialLayers (default) is specified, then erosion removes material and does not create a layer of thickness zero. If EventLayers is used, then erosion removes material and creates layers of thickness zero. Thus, EventLayers may be appropriate if the user is interested in chronostratigraphy.

  • dz_advection (float, (n_nodes, ) shape array, or at-node field array optional) – Change in rock elevation due to advection by some external process. This can be changed using the property setter.

  • rock_id (value or (n_nodes, ) shape array, optional) – Rock type id for new material if deposited. This can be changed using the property setter.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import LithoLayers
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("node", "topographic__elevation")

Create a LithoLayers with flatlying layers that altrnate between layers of type 1 and type 2 rock.

>>> z0s = [-4, -3, -2, -1, 0, 1, 2, 3, 4]
>>> ids = [1, 2, 1, 2, 1, 2, 1, 2, 1]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = LithoLayers(mg, z0s, ids, attrs)
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [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., 0., 0., 0., 0., 0.],
       [0., 0., 0., 0., 0., 0., 0., 0., 0.]])

Now create a set of layers that dip. Our anchor point will be the default value of (x0, y0) = (0, 0)

>>> lith = LithoLayers(mg, z0s, ids, attrs, function=lambda x, y: x + y)
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [0., 1., 1., 1., 1., 1., 1., 1., 1.],
       [0., 0., 1., 0., 1., 1., 1., 1., 1.],
       [0., 0., 0., 0., 0., 1., 0., 1., 1.],
       [0., 0., 0., 0., 0., 0., 0., 0., 1.],
       [0., 0., 0., 0., 0., 0., 0., 0., 0.]])

We can get the surface values, and as we’d expect, they alternate as the dipping layers are exposed at the surface.

>>> lith["K_sp"]
array([0.0001, 0.001 , 0.0001, 0.001 , 0.0001, 0.001 , 0.0001, 0.001 , 0.0001])
__init__(grid, z0s, ids, attrs, x0=0, y0=0, function=<function LithoLayers.<lambda>>, layer_type='EventLayers', dz_advection=0, rock_id=None)[source]

Create a new instance of a LithoLayers.

Parameters:
  • grid (Landlab ModelGrid)

  • z0s (ndarray of shape (n_layers, )) – Values of layer depth from surface at horizontal location (x0, y0).

  • ids (ndarray of shape (n_layers, )) – Values of rock type IDs corresponding to each layer specified in z0s.

  • attrs (dict) – Rock type property dictionary. See class docstring for example of required format.

  • x0 (float, optional) – x value of anchor point for all layers.

  • y0 (float, optional) – y value of anchor point for all layers.

  • function (function, optional) – Functional form of layers as a function of two variables, x and y. Default value is lambda x, y: 0*x + 0*y for flatlying layers.

  • layer_type (str, optional) – Type of Landlab layers object used to store the layers. If MaterialLayers (default) is specified, then erosion removes material and does not create a layer of thickness zero. If EventLayers is used, then erosion removes material and creates layers of thickness zero. Thus, EventLayers may be appropriate if the user is interested in chronostratigraphy.

  • dz_advection (float, (n_nodes, ) shape array, or at-node field array optional) – Change in rock elevation due to advection by some external process. This can be changed using the property setter.

  • rock_id (value or (n_nodes, ) shape array, optional) – Rock type id for new material if deposited. This can be changed using the property setter.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import LithoLayers
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("node", "topographic__elevation")

Create a LithoLayers with flatlying layers that altrnate between layers of type 1 and type 2 rock.

>>> z0s = [-4, -3, -2, -1, 0, 1, 2, 3, 4]
>>> ids = [1, 2, 1, 2, 1, 2, 1, 2, 1]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = LithoLayers(mg, z0s, ids, attrs)
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [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., 0., 0., 0., 0., 0.],
       [0., 0., 0., 0., 0., 0., 0., 0., 0.]])

Now create a set of layers that dip. Our anchor point will be the default value of (x0, y0) = (0, 0)

>>> lith = LithoLayers(mg, z0s, ids, attrs, function=lambda x, y: x + y)
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [0., 1., 1., 1., 1., 1., 1., 1., 1.],
       [0., 0., 1., 0., 1., 1., 1., 1., 1.],
       [0., 0., 0., 0., 0., 1., 0., 1., 1.],
       [0., 0., 0., 0., 0., 0., 0., 0., 1.],
       [0., 0., 0., 0., 0., 0., 0., 0., 0.]])

We can get the surface values, and as we’d expect, they alternate as the dipping layers are exposed at the surface.

>>> lith["K_sp"]
array([0.0001, 0.001 , 0.0001, 0.001 , 0.0001, 0.001 , 0.0001, 0.001 , 0.0001])
static __new__(cls, *args, **kwds)
add_layer(thickness, rock_id=None)

Add a new layer to Lithology.

Parameters:
  • thickness (float or (n_nodes,) array) – Positive values deposit material on to Lithology while negative values erode Lithology.

  • rock_id (single value or n_nodes long itterable, optional if only erosion occurs) – Rock type ID for new deposits. Can be single value or an number- of-nodes array.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]

We can instantiate Lithology with rock type properties we know we will use in the future.

>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001, 3: 0.01}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)

Add a layer of thickness 3 and rock type 3.

>>> lith.add_layer(3, rock_id=3)

The value of K_sp at node is now updated to the value of rock type 3

>>> mg.at_node["K_sp"]
array([0.01, 0.01, 0.01, 0.01, 0.01, 0.01, 0.01, 0.01, 0.01])

A negative value will erode. We can also pass a (n_nodes,) long array to erode unevenly. If all parts of the layer erode, then no `rock_id needs to be passed.

>>> erosion_amount = [-2.0, -2.0, -2.0, -4.0, -4.0, -4.0, -6.0, -6.0, -6.0]
>>> lith.add_layer(erosion_amount)
>>> mg.at_node["K_sp"]
array([0.01  , 0.01  , 0.01  , 0.0001, 0.0001, 0.0001, 0.001 ,
       0.001 , 0.001 ])

Now different layers are exposed at the surface and the value of K_sp is spatially variable.

add_property(attrs)

Add new property to Lithology.

Parameters:

attrs (dict) – Rock attribute dictionary for the new property(s).

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.add_property({"D": {1: 0.03, 2: 0.004}})
>>> lith.tracked_properties
['D', 'K_sp']
>>> mg.at_node["D"]
array([0.03, 0.03, 0.03, 0.03, 0.03, 0.03, 0.03, 0.03, 0.03])
add_rock_type(attrs)

Add rock type to Lithology.

Parameters:

attrs (dict) – Rock attribute dictionary for the new rock type(s).

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.add_rock_type({"K_sp": {4: 0.03, 6: 0.004}})
>>> lith.ids
[1, 2, 4, 6]
>>> lith.properties
{'K_sp': {1: 0.001, 2: 0.0001, 4: 0.03, 6: 0.004}}
cite_as = '\n    @article{barnhart2018lithology,\n        title = "Lithology: A Landlab submodule for spatially variable rock properties",\n        journal = "Journal of Open Source Software",\n        volume = "",\n        pages = "",\n        year = "2018",\n        doi = "10.21105/joss.00979",\n        author = {Katherine R. Barnhart and Eric Hutton and Nicole M. Gasparini\n                  and Gregory E. Tucker},\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 = ()
property dz

Thickness of each layer in the Lithology at each node.

The thickness of each layer in the Lithology as an array of shape (number_of_layers, number_of_nodes).

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [4., 4., 4., 4., 4., 4., 4., 4., 4.],
       [2., 2., 2., 2., 2., 2., 2., 2., 2.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.]])
property dz_advection

Rate of vertical advection.

Parameters:

dz_advection (float, (n_nodes, ) shape array, or at-node field array optional) – Change in rock elevation due to advection by some external process. This can be changed using the property setter. Dimensions are in length, not length per time.

Return type:

current rate of vertical advection

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.

property ids

Rock type IDs used by Lithology.

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 = ()
name = 'LithoLayers'
optional_var_names = ()
output_var_names = ()
property properties

Properties dictionary used by Lithology.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.properties
{'K_sp': {1: 0.001, 2: 0.0001}}
rock_cube_to_xarray(depths)

Construct a 3D rock cube of rock type ID as an xarray dataset.

Create an xarray dataset in (x, y, z) that shows the rock type with depth relative to the current topographic surface.

Here the z dimension is depth relative to the current topographic surface, NOT depth relative to an absolute datum.

Note also that when this method is called, it will construct the current values of lithology with depth, NOT the initial values.

Parameters:

depths (array)

Returns:

ds

Return type:

xarray dataset

property rock_id

Rock type for deposition.

Parameters:

rock_id (value or (n_nodes, ) shape array, optional) – Rock type id for new material if deposited. This can be changed using the property setter.

Return type:

current type of rock being deposited (if deposition occurs)

run_one_step()

Update Lithology.

The run_one_step method calculates elevation change of the Lithology surface (taking into account any advection due to external processes) and then either deposits or erodes based on elevation change.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_ones("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [4., 4., 4., 4., 4., 4., 4., 4., 4.],
       [2., 2., 2., 2., 2., 2., 2., 2., 2.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.]])
>>> lith.thickness
array([8., 8., 8., 8., 8., 8., 8., 8., 8.])

If we erode the surface, and then update Lithology, the thickness will change.

>>> z -= 0.5
>>> lith.run_one_step()
>>> lith.thickness
array([7.5, 7.5, 7.5, 7.5, 7.5, 7.5, 7.5, 7.5, 7.5])

The default of Lithology is to use MaterialLayers from the Landlab layers submodule. This means that when we erode, we will remove a layer from the layers datastructure if it has no material anywhere.

>>> lith.dz
array([[1. , 1. , 1. , 1. , 1. , 1. , 1. , 1. , 1. ],
       [4. , 4. , 4. , 4. , 4. , 4. , 4. , 4. , 4. ],
       [2. , 2. , 2. , 2. , 2. , 2. , 2. , 2. , 2. ],
       [0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5]])

We can see the value of the rock type at the surface.

>>> mg.at_node["rock_type__id"]
array([1., 1., 1., 1., 1., 1., 1., 1., 1.])

If you deposit, a valid rock_id must be provided. If the rock type is the same as the current surface value everywhere, then the layers will be combined. This rock_id can be provided as part of the init of Lithology or by setting a property (as shown below).

>>> z += 1.5
>>> lith.rock_id = 1
>>> lith.run_one_step()
>>> lith.thickness
array([9., 9., 9., 9., 9., 9., 9., 9., 9.])
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [4., 4., 4., 4., 4., 4., 4., 4., 4.],
       [2., 2., 2., 2., 2., 2., 2., 2., 2.],
       [2., 2., 2., 2., 2., 2., 2., 2., 2.]])

This contrasts with the behavior of Lithology if we use EventLayers. Next we repeat this example with EventLayers. Note that no matter which method you use, the values of the model grid fields will be the same. These two methods differ only in the details of the data structure they use to store the layer information.

>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_ones("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs, layer_type="EventLayers")
>>> lith.dz
array([[1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [4., 4., 4., 4., 4., 4., 4., 4., 4.],
       [2., 2., 2., 2., 2., 2., 2., 2., 2.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.]])
>>> lith.thickness
array([8., 8., 8., 8., 8., 8., 8., 8., 8.])

If we erode the surface, and then update Lithology, the thickness will change. However, with EventLayers, the lith.dz structure will be different. It will have a layer with thickness zero that represents the event of erosion.

>>> z -= 0.5
>>> lith.run_one_step()
>>> lith.thickness
array([7.5, 7.5, 7.5, 7.5, 7.5, 7.5, 7.5, 7.5, 7.5])
>>> lith.dz
array([[1. , 1. , 1. , 1. , 1. , 1. , 1. , 1. , 1. ],
       [4. , 4. , 4. , 4. , 4. , 4. , 4. , 4. , 4. ],
       [2. , 2. , 2. , 2. , 2. , 2. , 2. , 2. , 2. ],
       [0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5],
       [0. , 0. , 0. , 0. , 0. , 0. , 0. , 0. , 0. ]])

We can see the value of the rock type at the surface. As expected, it is just the same as if we used MaterialLayers.

>>> mg.at_node["rock_type__id"]
array([1., 1., 1., 1., 1., 1., 1., 1., 1.])

If you deposit, a valid rock_id must be provided. Unlike MaterialLayers, these two layers will not be combined, even if they have the same properties.

>>> z += 1.5
>>> lith.rock_id = 1
>>> lith.run_one_step()
>>> lith.thickness
array([9., 9., 9., 9., 9., 9., 9., 9., 9.])
>>> lith.dz
array([[1. , 1. , 1. , 1. , 1. , 1. , 1. , 1. , 1. ],
       [4. , 4. , 4. , 4. , 4. , 4. , 4. , 4. , 4. ],
       [2. , 2. , 2. , 2. , 2. , 2. , 2. , 2. , 2. ],
       [0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5],
       [0. , 0. , 0. , 0. , 0. , 0. , 0. , 0. , 0. ],
       [1.5, 1.5, 1.5, 1.5, 1.5, 1.5, 1.5, 1.5, 1.5]])
property shape

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

property thickness

Total thickness of the Lithology at each node.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.thickness
array([8., 8., 8., 8., 8., 8., 8., 8., 8.])
property tracked_properties

Properties tracked by Lithology.

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.tracked_properties
['K_sp']
unit_agnostic = True
units = ()
update_rock_properties(at, rock_id, value)

Update rock type attribute.

Parameters:
  • at (str) – Attribute name

  • rock_id (value) – Rock type ID

  • value (value) – New value for rock type attribute

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> mg.at_node["K_sp"]
array([0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001])
>>> lith.update_rock_properties("K_sp", 1, 0.03)
>>> mg.at_node["K_sp"]
array([0.03, 0.03, 0.03, 0.03, 0.03, 0.03, 0.03, 0.03, 0.03])
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 = ()
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

property z_bottom

Thickness from the surface to the bottom of each layer in Lithology.

Thickness from the topographic surface to the bottom of each layer as an array of shape (number_of_layers, number_of_nodes).

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.z_bottom
array([[8., 8., 8., 8., 8., 8., 8., 8., 8.],
       [7., 7., 7., 7., 7., 7., 7., 7., 7.],
       [3., 3., 3., 3., 3., 3., 3., 3., 3.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.]])
property z_top

Thickness from the surface to the top of each layer in Lithology.

Thickness from the topographic surface to the top of each layer as an array of shape (number_of_layers, number_of_nodes).

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import Lithology
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.add_zeros("topographic__elevation", at="node")
>>> thicknesses = [1, 2, 4, 1]
>>> ids = [1, 2, 1, 2]
>>> attrs = {"K_sp": {1: 0.001, 2: 0.0001}}
>>> lith = Lithology(mg, thicknesses, ids, attrs)
>>> lith.z_top
array([[7., 7., 7., 7., 7., 7., 7., 7., 7.],
       [3., 3., 3., 3., 3., 3., 3., 3., 3.],
       [1., 1., 1., 1., 1., 1., 1., 1., 1.],
       [0., 0., 0., 0., 0., 0., 0., 0., 0.]])