2.0.0b10
catchment modelling framework
Loading...
Searching...
No Matches

An additional water storage for a soil layer to model matrix water and macro pore water seperately. More...

+ Inheritance diagram for MacroPore:
+ Collaboration diagram for MacroPore:

Detailed Description

An additional water storage for a soil layer to model matrix water and macro pore water seperately.

Deprecated
The MacroPore model is still very experimental and not stable. Only for tryouts!

If present, the soil layer water storage holds the matrix water and the MacroPore holds the water in the macro pore. Use cmf::upslope::Macropore::create to create a macropore storage.

Use cmf::upslope::connections::GradientMacroFlow or cmf::upslope::connections::KinematicMacroFlow to model water flow between macro pores and a lateral connection (lateral subsurface fluxes) like cmf::upslope::connections::Richards_lateral to connect the macro pore with the matrix.

import cmf
c=p.NewCell(0,0,0,1000,True)
# Add two layers
c.add_layer(0.1, cmf.VanGenuchtenMualem())
c.add_layer(0.2, cmf.VanGenuchtenMualem())
# Create Macropore storages
macropores = [cmf.MacroPore.create(l,porefraction=0.05,Ksat=10.,density=0.05) for l in c.layers]
# Macro pore infiltration
cmf.KinematicMacroFlow(c.surfacewater, macropores[0])
# Macro pore percolation
cmf.KinematicMacroFlow(macropores[0], macropores[1])
# Macro / Micro exchange
for mp in macropores:
cmf.Richards_lateral(mp.layer, mp, FlowWidth = mp.get_flowwidth(), Distance = mp.density/2)
The study area, holding all cells, outlets and streams.
Definition project.h:53

Public Member Functions

real conc (cmf::math::Time t, const cmf::water::solute &_Solute) const override
 Returns the current WaterQuality (concentration of all solutes)
 
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.
 
cmf::water::flux_connectionconnection_to (const cmf::water::flux_node &target)
 Returns the connection between this and target.
 
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_capacity () const
 Returns the capacity of the macropores in m3.
 
cmf::upslope::Cellget_cell () const
 The cell of this macropore.
 
virtual real get_crackwidth () const
 Returns the crack width for a prismatic crackstructure.
 
real get_filled_fraction () const
 Get the relative water content in the macro pore \(\theta_{macro} = V_{macro}/V_{max}\).
 
real get_flowwidth () const
 The approximate length of the aggregate boundaries.
 
virtual real get_K () const
 Returns the actual conductivity.
 
virtual real get_K (cmf::geometry::point direction) const
 Returns the actual anisotropic conductivity along a direction \(K = (k_f \cdot d) K\).
 
SoilLayer::ptr get_layer () const
 Gets the soil layer (matrix water storage) for this macropore storage.
 
real get_porefraction () const
 The fraction of the macro pores in m3/m3. This adds to the porosity of the layer.
 
real get_potential (cmf::math::Time t=cmf::math::never) const
 Returns the actual water level in the macropore in m above reference.
 
cmf::projectget_project () const
 Returns the project, this node is part of.
 
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)
 
real get_volume () const
 Returns the actual stored volume in this macropore in m3.
 
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_potential (real waterhead)
 Sets the water level in the macropore. Be aware of not setting it below the lower boundary.
 
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)
 
void set_volume (real volume)
 Sets the volume of stored water in m3.
 
SoluteStorage & Solute (const cmf::water::solute _Solute)
 Returns the water quality of the water storage.
 
real waterbalance (cmf::math::Time t, const flux_connection *Without=0) const
 Returns the sum of all fluxes (positive and negative) at time t.
 

Static Public Member Functions

static MacroPore::ptr create (SoilLayer::ptr layer, real porefraction=0.05, real Ksat=10, real density=0.05, real porefraction_wilt=-1., real K_shape=0.0)
 Creates a macropore water storage for a soil layer.
 

Public Attributes

real crack_wetness
 The layer's cracking wetness, crack dynamic starts below this value. The value is ignored if porefraction_min = porefraction_max.
 
real density
 The mean distance between the macro pores in m.
 
real K_shape
 The shape of the conductivity \(f_K\) in relation to the matric potential of the micropore system.
 
real Ksat
 The saturated conductivity of the macropores in m/day.
 
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

◆ create()

static MacroPore::ptr create ( SoilLayer::ptr layer,
real porefraction = 0.05,
real Ksat = 10,
real density = 0.05,
real porefraction_wilt = -1.,
real K_shape = 0.0 )
static

Creates a macropore water storage for a soil layer.

Returns
a MacroPore water storage
Parameters
layerThe soil layer holding the macro pore
porefractionThe relative macro pore content of the soil layer
KsatThe saturated conductivity of the macro pores
densityThe mean distance between macropores or mean aggregate size
porefraction_wiltThe fraction of macropore space at wilting point. Use this parameter for swelling soils
K_shapeA parameter to change the conductivity of the macropore with the matrix potential

◆ get_capacity()

real get_capacity ( ) const

Returns the capacity of the macropores in m3.

\[V_{max} = \Phi_{macro} A \Delta z\]

where:

  • \(V_{max}\) is the water capacity of the macropore
  • \(\Phi_{macro}\) is the fraction of macro pores in the soil in \(m^3 macro pores/m^3 soil\)
  • \(A \Delta z\) is the total volume of the soil layer (including all pores) in \(m^3\)

◆ get_crackwidth()

virtual real get_crackwidth ( ) const
virtual

Returns the crack width for a prismatic crackstructure.

For a prismatic crack structure, the porefraction in m3/m3 equals the vertical crack area in m2/m2. The length of equally spaced cracks is in one direction the inverse of the density and twice the length for two directions.

\[ l_{crack} [m/m^2]= 2 \frac {1}{d[m]}\]

If we again ignore the fact that the spacing of the cracking crossings is counted double, the crack width is:

\[ w_{crack}[m] = \frac{A_{crack}[m^2/m^2]}{l_{crack}[m/m^2]} \]

Combining both eq. above:

\[ w_{crack}[m] = A_{crack}[m^2/m^2]\frac{d[m]}{2} \]

◆ get_flowwidth()

real get_flowwidth ( ) const

The approximate length of the aggregate boundaries.

\[l = \frac{2}{d_{macro}} A\]

where:

  • \(l\) is the length of the aggregate boundaries (in m)
  • \(2\) is the number of directions
  • \(d_{macro}\) is the mean distance between macropores (density) in m
  • \(A\) is the area of the cell

◆ 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)

Member Data Documentation

◆ K_shape

real K_shape

The shape of the conductivity \(f_K\) in relation to the matric potential of the micropore system.

\[K = K_{macro} \exp(f_K \Psi_M)\]