OGS
ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >

Definition at line 80 of file RichardsFlowFEM.h.

#include <RichardsFlowFEM.h>

Inheritance diagram for ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 LocalAssemblerData (MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, RichardsFlowProcessData const &process_data)
 
void assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
std::vector< double > const & getIntPtSaturation (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
 
std::vector< double > const & getIntPtDarcyVelocity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
- Public Member Functions inherited from ProcessLib::RichardsFlow::RichardsFlowLocalAssemblerInterface
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
 
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
 
virtual void preAssemble (double const, double const, std::vector< double > const &)
 
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
 
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using LocalAssemblerTraits
 
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
 
using GlobalDimNodalMatrixType
 
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
 
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
 

Private Attributes

MeshLib::Element const & _element
 
RichardsFlowProcessData const & _process_data
 
NumLib::GenericIntegrationMethod const & _integration_method
 
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
 
std::vector< double > _saturation
 

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 93 of file RichardsFlowFEM.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::GlobalDimNodalMatrixType
private
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 89 of file RichardsFlowFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 94 of file RichardsFlowFEM.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalAssemblerTraits
private
Initial value:

Definition at line 85 of file RichardsFlowFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
private

Definition at line 91 of file RichardsFlowFEM.h.

◆ NodalRowVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 88 of file RichardsFlowFEM.h.

◆ NodalVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 92 of file RichardsFlowFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 83 of file RichardsFlowFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 82 of file RichardsFlowFEM.h.

Constructor & Destructor Documentation

◆ LocalAssemblerData()

template<typename ShapeFunction , int GlobalDim>
ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalAssemblerData ( MeshLib::Element const & element,
std::size_t const local_matrix_size,
NumLib::GenericIntegrationMethod const & integration_method,
bool is_axially_symmetric,
RichardsFlowProcessData const & process_data )
inline

Definition at line 97 of file RichardsFlowFEM.h.

103 : _element(element),
104 _process_data(process_data),
105 _integration_method(integration_method),
107 std::vector<double>(_integration_method.getNumberOfPoints()))
108 {
109 // This assertion is valid only if all nodal d.o.f. use the same shape
110 // matrices.
111 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
112 (void)local_matrix_size;
113
114 unsigned const n_integration_points =
116 _ip_data.reserve(n_integration_points);
117
118 auto const shape_matrices =
120 GlobalDim>(element, is_axially_symmetric,
122
123 for (unsigned ip = 0; ip < n_integration_points; ip++)
124 {
125 auto const& sm = shape_matrices[ip];
126 const double integration_factor = sm.integralMeasure * sm.detJ;
127 _ip_data.emplace_back(
128 sm.N, sm.dNdx,
130 integration_factor,
131 sm.N.transpose() * sm.N * integration_factor *
133 }
134 }
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
RichardsFlowProcessData const & _process_data
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)

References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_ip_data, NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), and ProcessLib::RichardsFlow::NUM_NODAL_DOF.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int GlobalDim>
void ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & ,
std::vector< double > & local_M_data,
std::vector< double > & local_K_data,
std::vector< double > & local_b_data )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 136 of file RichardsFlowFEM.h.

142 {
143 auto const local_matrix_size = local_x.size();
144 // This assertion is valid only if all nodal d.o.f. use the same shape
145 // matrices.
146 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
147
149 local_M_data, local_matrix_size, local_matrix_size);
151 local_K_data, local_matrix_size, local_matrix_size);
153 local_b_data, local_matrix_size);
154
155 unsigned const n_integration_points =
159
160 auto const& medium =
162 auto const& liquid_phase = medium.phase("AqueousLiquid");
164 vars.temperature =
165 medium
166 .property(
168 .template value<double>(vars, pos, t, dt);
169
170 for (unsigned ip = 0; ip < n_integration_points; ip++)
171 {
172 pos.setIntegrationPoint(ip);
173 double p_int_pt = 0.0;
174 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
175
176 vars.liquid_phase_pressure = p_int_pt;
177 vars.capillary_pressure = -p_int_pt;
178 // setting pG to 1 atm
179 // TODO : rewrite equations s.t. p_L = pG-p_cap
180 vars.gas_phase_pressure = 1.0e5;
181
182 auto const permeability =
184 medium.property(MaterialPropertyLib::permeability)
185 .value(vars, pos, t, dt));
186
187 auto const porosity =
189 .template value<double>(vars, pos, t, dt);
190
191 double const Sw =
193 .template value<double>(vars, pos, t, dt);
194 _saturation[ip] = Sw;
195 vars.liquid_saturation = Sw;
196
197 double const dSw_dpc =
199 .template dValue<double>(
201 pos, t, dt);
202
203 auto const drhow_dp =
204 liquid_phase
206 .template dValue<double>(
207 vars,
209 pos, t, dt);
210 auto const storage =
212 .template value<double>(vars, pos, t, dt);
213 double const mass_mat_coeff =
214 storage * Sw + porosity * Sw * drhow_dp - porosity * dSw_dpc;
215
216 local_M.noalias() += mass_mat_coeff * _ip_data[ip].mass_operator;
217
218 double const k_rel =
219 medium
221 relative_permeability)
222 .template value<double>(vars, pos, t, dt);
223
224 auto const mu =
225 liquid_phase.property(MaterialPropertyLib::viscosity)
226 .template value<double>(vars, pos, t, dt);
227 local_K.noalias() += _ip_data[ip].dNdx.transpose() * permeability *
228 _ip_data[ip].dNdx *
229 _ip_data[ip].integration_weight * (k_rel / mu);
230
232 {
233 auto const rho_w =
234 liquid_phase
236 .template value<double>(vars, pos, t, dt);
237 auto const& body_force = _process_data.specific_body_force;
238 assert(body_force.size() == GlobalDim);
239 NodalVectorType gravity_operator =
240 _ip_data[ip].dNdx.transpose() * permeability * body_force *
241 _ip_data[ip].integration_weight;
242 local_b.noalias() += (k_rel / mu) * rho_w * gravity_operator;
243 }
244 }
246 {
247 for (int idx_ml = 0; idx_ml < local_M.cols(); idx_ml++)
248 {
249 double const mass_lump_val = local_M.col(idx_ml).sum();
250 local_M.col(idx_ml).setZero();
251 local_M(idx_ml, idx_ml) = mass_lump_val;
252 }
253 } // end of mass lumping
254 }
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
typename ShapeMatricesType::NodalVectorType NodalVectorType
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map

References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_saturation, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::VariableArray::capillary_pressure, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::gas_phase_pressure, MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::RichardsFlow::RichardsFlowProcessData::has_gravity, ProcessLib::RichardsFlow::RichardsFlowProcessData::has_mass_lumping, MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, ProcessLib::RichardsFlow::RichardsFlowProcessData::media_map, ProcessLib::RichardsFlow::NUM_NODAL_DOF, MaterialPropertyLib::permeability, MaterialPropertyLib::Medium::phase(), MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), ProcessLib::RichardsFlow::RichardsFlowProcessData::specific_body_force, MaterialPropertyLib::storage, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction , int GlobalDim>
std::vector< double > const & ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::RichardsFlow::RichardsFlowLocalAssemblerInterface.

Definition at line 275 of file RichardsFlowFEM.h.

280 {
281 // TODO (tf) Temporary value not used by current material models. Need
282 // extension of secondary variable interface
283 double const dt = std::numeric_limits<double>::quiet_NaN();
284
285 constexpr int process_id = 0; // monolithic scheme.
286 auto const indices =
287 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
288 assert(!indices.empty());
289 auto const local_x = x[process_id]->get(indices);
290
293
294 auto const& medium =
296 auto const& liquid_phase = medium.phase("AqueousLiquid");
297
299 vars.temperature =
300 medium
301 .property(
303 .template value<double>(vars, pos, t, dt);
304
305 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
306 &local_x[0], ShapeFunction::NPOINTS);
307
308 unsigned const n_integration_points =
310
311 cache.clear();
312 auto cache_vec = MathLib::createZeroedMatrix<
313 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
314 cache, GlobalDim, n_integration_points);
315
316 for (unsigned ip = 0; ip < n_integration_points; ++ip)
317 {
318 double p_int_pt = 0.0;
319 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
320 vars.liquid_phase_pressure = p_int_pt;
321 vars.capillary_pressure = -p_int_pt;
322 // setting pG to 1 atm
323 // TODO : rewrite equations s.t. p_L = pG-p_cap
324 vars.gas_phase_pressure = 1.0e5;
325
326 double const Sw =
328 .template value<double>(vars, pos, t, dt);
329 vars.liquid_saturation = Sw;
330
331 auto const permeability =
333 medium.property(MaterialPropertyLib::permeability)
334 .value(vars, pos, t, dt));
335
336 double const k_rel =
337 medium
339 relative_permeability)
340 .template value<double>(vars, pos, t, dt);
341 auto const mu =
342 liquid_phase.property(MaterialPropertyLib::viscosity)
343 .template value<double>(vars, pos, t, dt);
344 auto const K_mat_coeff = permeability * (k_rel / mu);
345 cache_vec.col(ip).noalias() =
346 -K_mat_coeff * _ip_data[ip].dNdx * p_nodal_values;
348 {
349 auto const rho_w =
350 liquid_phase
352 .template value<double>(vars, pos, t, dt);
353 auto const& body_force = _process_data.specific_body_force;
354 assert(body_force.size() == GlobalDim);
355 // here it is assumed that the vector body_force is directed
356 // 'downwards'
357 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
358 }
359 }
360
361 return cache;
362 }
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, MaterialPropertyLib::VariableArray::capillary_pressure, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::gas_phase_pressure, MeshLib::Element::getID(), NumLib::getIndices(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::RichardsFlow::RichardsFlowProcessData::has_gravity, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, ProcessLib::RichardsFlow::RichardsFlowProcessData::media_map, MaterialPropertyLib::permeability, MaterialPropertyLib::Medium::phase(), MaterialPropertyLib::reference_temperature, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), ProcessLib::RichardsFlow::RichardsFlowProcessData::specific_body_force, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

◆ getIntPtSaturation()

template<typename ShapeFunction , int GlobalDim>
std::vector< double > const & ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtSaturation ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > &  ) const
inlineoverridevirtual

◆ getShapeMatrix()

template<typename ShapeFunction , int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 256 of file RichardsFlowFEM.h.

258 {
259 auto const& N = _ip_data[integration_point].N;
260
261 // assumes N is stored contiguously in memory
262 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
263 }

References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_ip_data.

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _process_data

◆ _saturation

template<typename ShapeFunction , int GlobalDim>
std::vector<double> ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_saturation
private

The documentation for this class was generated from the following file: