landlab.components.tidal_flow.tidal_flow_calculator

Calculate cycle-averaged tidal flow field using approach of Mariotti (2018)

class TidalFlowCalculator[source]

Bases: Component

Component that calculates average flow over a tidal cycle.

The TidalFlowCalculator component calculates a tidal flow velocity field on a grid using the method of Mariotti et al. (2018). The grid can be raster or hex. In essence, the method uses a numerical solution to the steady diffusion equation. The resulting velocity fields—one for flood tide, and on for ebb tide—represent the average velocity over a tidal half-cycle. Velocity is obtained by starting with the total depth of water that is added or removed during flood tide or ebb tide, respectively, at each grid point. The idea is to solve for the water-surface elevation field that produces this total inflow or outflow of water at each point, using a linear approximation for shallow-water flow. The math is given in Mariotti (2018), and is summarized below.

Tidal-averaged water depth: with z as bed surface elevation, and r as tidal range, the water depth h is:

\[h = [\max(0, r/2 − z) + \max(0, −r/2 − z)]/2\]

Horizontal flow velocity (2d vector): with \(\eta\) as water-surface elevation, n as Manning roughness, and \(\chi\) as a scale velocity (here unity), depth-averaged flow velocity U is:

\[U = \frac{h^{4/3}}{n^2\chi} \nabla \eta\]

Tidal inundation / drainage rate, I: at any given point, this is the depth of water added (flood tide) or drained (ebb tide) divided by the tidal half period, T/2. It is defined as:

\[I = [r/2 − \max(−r/2, \min(z, r/2))] / (T/2)\]

Mass conservation:

\[\nabla \cdot (h U) = I\]

Because h is assumed constant, this becomes a steady diffusion equation:

\[\nabla^2 \eta = \frac{I n^{4/3} \chi}{h^{7/3}}\]

This is a Poisson (steady diffusion) equation, which is solved numerically at grid nodes using a finite-volume method. The water-surface gradient is then used to calculate the velocity field, using the above equation for U.

Parameters:
  • grid (RasterModelGrid or HexModelGrid) – A Landlab grid object.

  • tidal_range (float, optional) – Tidal range (2x tidal amplitude) (m) (default 1)

  • tidal_period (float, optional) – Tidal perioid (s) (default M2 tidal period = 12 h 25 min)

  • roughness (float, array, or field name; optional) – Manning roughness coefficient (“n”) (s/m^1/3) (default 0.01)

  • mean_sea_level (float, optional) – Mean sea level (m) (default 0)

  • scale_velocity (float, optional) – Scale velocity (see Mariotti, 2018) (m/s) (default 1)

  • min_water_depth (float, optional) – Minimum depth for calculating diffusion coefficient (m) (default 0.01)

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.components import TidalFlowCalculator
>>> grid = RasterModelGrid((3, 5), xy_spacing=2.0)  # 1 row core nodes
>>> grid.set_closed_boundaries_at_grid_edges(False, True, True, True)
>>> z = grid.add_zeros("topographic__elevation", at="node")
>>> z[:] = -50.0  # mean water depth is 50 m below MSL
>>> tfc = TidalFlowCalculator(
...     grid, tidal_range=2.0, tidal_period=4.0e4, roughness=0.01
... )
>>> tfc.run_one_step()
>>> int(round(grid.at_link["ebb_tide_flow__velocity"][10] * 1.0e6))
4

References

Mariotti, G. (2018) Marsh channel morphological response to sea level rise and sediment supply. Estuarine, Coastal and Shelf Science, 209, 89–101, https://doi.org/10.1016/j.ecss.2018.05.016.

Initialize TidalFlowCalculator.

__init__(grid, tidal_range=1.0, tidal_period=44712.0, roughness=0.01, mean_sea_level=0.0, scale_velocity=1.0, min_water_depth=0.01)[source]

Initialize TidalFlowCalculator.

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

Calculate and store the rate of inundation/draining at each node, averaged over a tidal half-cycle.

Examples

>>> grid = RasterModelGrid((3, 5))
>>> z = grid.add_zeros("topographic__elevation", at="node")
>>> z[5:10] = [10.0, 0.25, 0.0, -0.25, -10.0]
>>> period = 4.0e4  # tidal period in s, for convenient calculation
>>> tfc = TidalFlowCalculator(grid, tidal_period=period)
>>> rate = tfc.calc_tidal_inundation_rate()
>>> 0.5 * rate[5:10] * period  # depth in m
array([0.  , 0.25, 0.5 , 0.75, 1.  ])
>>> rate[5:10]  # rate in m/s
array([0.00e+00, 1.25e-05, 2.50e-05, 3.75e-05, 5.00e-05])

Notes

This calculates I in Mariotti (2018) using his equation (1).

cite_as = ''
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 = (('ebb_tide_flow__velocity', 'Horizontal flow velocity along links during ebb tide'), ('flood_tide_flow__velocity', 'Horizontal flow velocity along links during flood tide'), ('mean_water__depth', 'Tidal mean water depth'), ('topographic__elevation', 'Land surface topographic elevation'))
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 = ('topographic__elevation',)
property mean_sea_level

Mean sea level.

name = 'TidalFlowCalculator'
optional_var_names = ()
output_var_names = ('ebb_tide_flow__velocity', 'flood_tide_flow__velocity', 'mean_water__depth')
property roughness

Roughness coefficient (Manning’s n).

run_one_step()[source]

Calculate the tidal flow field and water-surface elevation.

property shape

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

property tidal_period

Tidal period.

property tidal_range

Tidal range.

unit_agnostic = False
units = (('ebb_tide_flow__velocity', 'm/s'), ('flood_tide_flow__velocity', 'm/s'), ('mean_water__depth', 'm'), ('topographic__elevation', '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 = (('ebb_tide_flow__velocity', 'link'), ('flood_tide_flow__velocity', 'link'), ('mean_water__depth', 'node'), ('topographic__elevation', '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