sequence.shoreline module#
Find the shoreline of a SequenceModelGrid.
This module contains methods for calculating a grid’s shoreline and shelf edge.
- class sequence.shoreline.ShorelineFinder(*args, **kwds)[source]#
Bases:
Component
Find a grid’s shoreline.
- __init__(grid: SequenceModelGrid, alpha: float = 0.0005)[source]#
Create a shoreline finder.
- Parameters:
grid (SequenceModelGrid) – A Sequence grid.
alpha (float, optional) – Parameter used to find the shoreline when interpolating between nodes.
- sequence.shoreline.find_shelf_edge(x: ndarray[Any, dtype[floating]], dz: ndarray[Any, dtype[floating]], x_of_shore: ndarray[Any, dtype[floating]] | float = 0.0, alpha: float = 0.0005) float [source]#
Find the shelf edge based on deposit thickness.
- Parameters:
- Returns:
The x-coordinate of the shelf edge.
- Return type:
- Raises:
ShelfEdgeError – If the data do not contain the shelf edge.
Examples
>>> import numpy as np >>> from sequence.shoreline import find_shelf_edge
Create a depositional profile where the maximum deposition is at x=1.0.
>>> x = np.arange(50.) / 5.0 >>> deposit = x * np.exp(-x) >>> find_shelf_edge(x, deposit, x_of_shore=0.0, alpha=1e3) 1.0
- sequence.shoreline.find_shelf_edge_by_curvature(x: ndarray[Any, dtype[floating]], z: ndarray[Any, dtype[floating]], sea_level: float = 0.0) float [source]#
Find the x-coordinate of the shelf edge.
The shelf edge is the location where the curvature of sea-floor elevations is a minimum.
- Parameters:
- Returns:
The x-coordinate of the shelf edge.
- Return type:
Examples
>>> import numpy as np >>> from sequence.shoreline import find_shelf_edge_by_curvature
>>> x = np.arange(50.) >>> z = - np.arctan(x / 5.0 - 5.0) >>> find_shelf_edge_by_curvature(x, z, sea_level=z.max()) 22.0
- sequence.shoreline.find_shoreline(x: ndarray[Any, dtype[floating]], z: ndarray[Any, dtype[floating]], sea_level: float = 0.0, kind: str = 'cubic') ndarray[Any, dtype[floating]] | float [source]#
Find the shoreline of a profile.
- Parameters:
x (array of float) – X-positions of profile.
z (array of float) – Elevations along the profile.
sea_level (float, optional) – Elevation of sea level.
kind (str, optional) – Interpolation method used to find shoreline. Values are the same as those used for scipy.interpolate.interp1d. Default is ‘cubic’.
- Returns:
X-position of the shoreline.
- Return type:
Examples
>>> from sequence.shoreline import find_shoreline >>> import numpy as np
Create a linearly-dipping profile.
>>> x = np.arange(10.) >>> z = - x + 5. >>> z array([ 5., 4., 3., 2., 1., 0., -1., -2., -3., -4.])
Find the shoreline.
>>> find_shoreline(x, z, kind='linear') 5.0 >>> find_shoreline(x, z, sea_level=.25, kind='linear') 4.75
If sea level is higher/lower than the max/min elevation, return the first/last x value.
>>> find_shoreline(x, z, sea_level=100., kind='linear') 0.0 >>> find_shoreline(x, z, sea_level=-100., kind='linear') 9.0
- sequence.shoreline.find_shoreline_index(x: ndarray[Any, dtype[floating]], z: ndarray[Any, dtype[floating]], sea_level: float = 0.0) ndarray[Any, dtype[integer]] | int [source]#
Find the landward-index of the shoreline.
- Parameters:
- Returns:
Index into z landward of the shoreline.
- Return type:
- Raises:
ValueError – If the profile does not contain a shoreline.
Examples
>>> from sequence.shoreline import find_shoreline_index >>> import numpy as np
Create a linearly-dipping profile.
>>> x = np.arange(10.) >>> z = - x + 5. >>> z array([ 5., 4., 3., 2., 1., 0., -1., -2., -3., -4.])
Find the shoreline.
>>> find_shoreline_index(x, z) 6 >>> find_shoreline_index(x, z, sea_level=.25) 5
If sea level is higher/lower than the max/min elevation, raise a ShorelineError.
>>> find_shoreline_index(x, z, sea_level=100.) Traceback (most recent call last): sequence.errors.ShorelineError: No shoreline found. The profile is all below sea level >>> find_shoreline_index(x, z, sea_level=-100.) Traceback (most recent call last): sequence.errors.ShorelineError: No shoreline found. The profile is all above sea level