2.0.0b10
catchment modelling framework
Loading...
Searching...
No Matches
Reach Class Reference

A reach represents the section of a riover and is a specialization of an open water storage. More...

+ Inheritance diagram for Reach:
+ Collaboration diagram for Reach:

Detailed Description

A reach represents the section of a riover and is a specialization of an open water storage.

The OpenWaterStorage attributes and methods are extended by topological features, for the creation of a network of reaches.

Public Member Functions

real conc (cmf::math::Time t, const cmf::water::solute &solute) const
 Returns the water quality of the flux_node, if it is not overridden this is the mix of the incoming fluxes.
 
real conc (const cmf::water::solute &_Solute) const
 Returns the concentration of the given solute.
 
void conc (const cmf::water::solute &_Solute, real NewConcetration)
 Sets a new concentration.
 
void connect_to_surfacewater (cmf::upslope::Cell *cell, real width, bool diffusive)
 Connects the surfacewater of a cell with this reach.
 
cmf::water::flux_connectionconnection_to (const cmf::water::flux_node &target)
 Returns the connection between this and target.
 
double distance_to_cell (cmf::upslope::Cell *cell) const
 Returns the distance (d) for connections between this reach and a cell.
 
virtual real dxdt (const cmf::math::Time &time)
 Returns the derivate of the state variable at time time.
 
real flux_to (const cmf::water::flux_node &target, cmf::math::Time t)
 Returns the actual flux between this and target (positive sign means "from this into target")
 
cmf::geometry::point get_3d_flux (cmf::math::Time t)
 Returns the sum of all flux vectors.
 
real get_depth () const
 Returns the water table depth.
 
bool get_diffusive () const
 Returns if flow is calculated using a diffusive wave model.
 
cmf::water::flux_node::ptr get_downstream () const
 Returns the reach downstream of this (or null if there is no reach downstream)
 
virtual const IVolumeHeightFunctionget_height_function () const
 The functional relation between volume, depth and exposed area.
 
real get_length () const
 Returns the length of the reach.
 
real get_potential (cmf::math::Time t=cmf::math::never) const override
 Returns the water potential of the node in m waterhead.
 
cmf::projectget_project () const
 Returns the project, this node is part of.
 
Channel get_reachtype () const
 Returns the channel shape.
 
ptr get_root ()
 Returns the reach most downstream from this reach

 
real get_state () const
 Returns the current state of the variable.
 
char get_state_variable_content () const
 A character indicating the integrated variable (either 'V' for Volume or 'h' for head)
 
ptr get_upstream (int index) const
 Returns a reach upstream of this.
 
virtual real get_volume () const
 Returns the volume of water in this storage in m3
 
real get_width () const
 Returns the average width of the reach.
 
virtual bool is_connected (const cmf::math::StateVariable &other) const
 Returns True if this waterstorage is effected by another state.
 
bool is_storage () const override
 Returns true, since this is a storage.
 
real operator() (cmf::math::Time t) const
 returns the waterblance
 
bool remove_connection (cmf::water::flux_node::ptr To)
 Remove the connection.
 
void set_dead_end ()
 Deletes any downstream connection.
 
void set_diffusive (bool use_diffusive_wave)
 Sets all Manning kind connections to either diffusive or kinematic.
 
void set_downstream (ptr new_downstream, bool use_meanchannel=false)
 Connects the reach to another one downstream.
 
void set_height_function (const IChannel &val)
 Sets the channel shape.
 
void set_outlet (cmf::water::flux_node::ptr outlet)
 Connects the reach to an outlet, e.g. a boundary condition.
 
void set_potential (real newpotential) override
 Sets the potential of this flux node.
 
void set_state (real newState)
 Gives access to the state variable.
 
void set_state_variable_content (char content)
 A character indicating the integrated variable (either 'V' for Volume or 'h' for head)
 
virtual void set_volume (real newwatercontent)
 Sets the volume of water in this storage in m3
 
SoluteStorage & Solute (const cmf::water::solute _Solute)
 Returns the water quality of the water storage.
 
int upstream_count () const
 Returns the number of reaches upstream of this.
 
real waterbalance (cmf::math::Time t, const flux_connection *Without=0) const
 Returns the sum of all fluxes (positive and negative) at time t.
 
real wet_area () const
 Returns the exposed surface area in m2.
 

Static Public Member Functions

static ptr cast (cmf::water::flux_node::ptr node)
 Casts a flux node to an open water storage.
 
static ptr create (cmf::project &_project, const cmf::river::IVolumeHeightFunction &base_geo)
 Creates an open water storage with any type of a volume.
 
static ptr create (cmf::project &_project, real Area)
 Creates an open water storage with a prismatic volume

 

Public Attributes

std::string Name
 The Name of this node.
 
const int node_id
 The Id of the node.
 
cmf::geometry::point position
 The spatial position of the node.
 

Protected Member Functions

void MarkStateChangeHandled ()
 Sets the updated flag (m_StateIsNew) to false.
 
bool StateIsChanged ()
 Returns if the state was currently updated.
 

Overrides of flux_node

virtual bool RecalcFluxes (cmf::math::Time t)
 Pure flux_nodes do not influence fluxes, therefore no recalculation of fluxes is required by flux_node.
 
virtual double is_empty () const
 Returns true if the node has no water.
 

Member Function Documentation

◆ connect_to_surfacewater()

void connect_to_surfacewater ( cmf::upslope::Cell * cell,
real width,
bool diffusive )

Connects the surfacewater of a cell with this reach.

Parameters
cellThe cell with the surface water to be connected with this reach
widthThe flow width from cell to this reach
diffusiveIf diffusive is false, a ManningKinematic connection is used, else a ManningDiffusive connection

◆ distance_to_cell()

double distance_to_cell ( cmf::upslope::Cell * cell) const

Returns the distance (d) for connections between this reach and a cell.

If the effective inner cell distance (defined as \( R_{Cell} = 0.5\frac{\sqrt{A}}{\pi}\)) is smaller than the distance between the center points, the cell radius is returned

◆ get_potential()

real get_potential ( cmf::math::Time = cmf::math::never) const
overridevirtualinherited

Returns the water potential of the node in m waterhead.

The base class water storage always returns the height of the location

Reimplemented from flux_node.

◆ RecalcFluxes()

virtual bool RecalcFluxes ( cmf::math::Time t)
virtualinherited

Pure flux_nodes do not influence fluxes, therefore no recalculation of fluxes is required by flux_node.

WaterStorage overrides this, since state changes require an update of the fluxes

Reimplemented from flux_node.

◆ waterbalance()

real waterbalance ( cmf::math::Time t,
const flux_connection * Without = 0 ) const
inherited

Returns the sum of all fluxes (positive and negative) at time t.


Single fluxes can be excluded from the calculation

Parameters
tTime of the query
WithoutA flux_connection that is excluded from the waterbalance (e.g. to prevent closed circuits)