Loading [MathJax]/extensions/tex2jax.js
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.
 
virtual std::optional< VectorSegmentgetVectorDeformationSegment () const
 
- 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 }
constexpr 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 double p_int_pt = 0.0;
173 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
174
175 pos = {std::nullopt, _element.getID(),
177 NumLib::interpolateCoordinates<ShapeFunction,
179 _element, _ip_data[ip].N))};
180
181 vars.liquid_phase_pressure = p_int_pt;
182 vars.capillary_pressure = -p_int_pt;
183 // setting pG to 1 atm
184 // TODO : rewrite equations s.t. p_L = pG-p_cap
185 vars.gas_phase_pressure = 1.0e5;
186
187 auto const permeability =
189 medium.property(MaterialPropertyLib::permeability)
190 .value(vars, pos, t, dt));
191
192 auto const porosity =
194 .template value<double>(vars, pos, t, dt);
195
196 double const Sw =
198 .template value<double>(vars, pos, t, dt);
199 _saturation[ip] = Sw;
200 vars.liquid_saturation = Sw;
201
202 double const dSw_dpc =
204 .template dValue<double>(
206 pos, t, dt);
207
208 auto const drhow_dp =
209 liquid_phase
211 .template dValue<double>(
212 vars,
214 pos, t, dt);
215 auto const storage =
217 .template value<double>(vars, pos, t, dt);
218 double const mass_mat_coeff =
219 storage * Sw + porosity * Sw * drhow_dp - porosity * dSw_dpc;
220
221 local_M.noalias() += mass_mat_coeff * _ip_data[ip].mass_operator;
222
223 double const k_rel =
224 medium
226 relative_permeability)
227 .template value<double>(vars, pos, t, dt);
228
229 auto const mu =
230 liquid_phase.property(MaterialPropertyLib::viscosity)
231 .template value<double>(vars, pos, t, dt);
232 local_K.noalias() += _ip_data[ip].dNdx.transpose() * permeability *
233 _ip_data[ip].dNdx *
234 _ip_data[ip].integration_weight * (k_rel / mu);
235
237 {
238 auto const rho_w =
239 liquid_phase
241 .template value<double>(vars, pos, t, dt);
242 auto const& body_force = _process_data.specific_body_force;
243 assert(body_force.size() == GlobalDim);
244 NodalVectorType gravity_operator =
245 _ip_data[ip].dNdx.transpose() * permeability * body_force *
246 _ip_data[ip].integration_weight;
247 local_b.noalias() += (k_rel / mu) * rho_w * gravity_operator;
248 }
249 }
251 {
252 for (int idx_ml = 0; idx_ml < local_M.cols(); idx_ml++)
253 {
254 double const mass_lump_val = local_M.col(idx_ml).sum();
255 local_M.col(idx_ml).setZero();
256 local_M(idx_ml, idx_ml) = mass_lump_val;
257 }
258 } // end of mass lumping
259 }
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)
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 &)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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, NumLib::interpolateCoordinates(), 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(), 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 280 of file RichardsFlowFEM.h.

285 {
286 // TODO (tf) Temporary value not used by current material models. Need
287 // extension of secondary variable interface
288 double const dt = std::numeric_limits<double>::quiet_NaN();
289
290 constexpr int process_id = 0; // monolithic scheme.
291 auto const indices =
292 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
293 assert(!indices.empty());
294 auto const local_x = x[process_id]->get(indices);
295
298
299 auto const& medium =
301 auto const& liquid_phase = medium.phase("AqueousLiquid");
302
304 vars.temperature =
305 medium
306 .property(
308 .template value<double>(vars, pos, t, dt);
309
310 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
311 &local_x[0], ShapeFunction::NPOINTS);
312
313 unsigned const n_integration_points =
315
316 cache.clear();
317 auto cache_vec = MathLib::createZeroedMatrix<
318 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
319 cache, GlobalDim, n_integration_points);
320
321 for (unsigned ip = 0; ip < n_integration_points; ++ip)
322 {
323 double p_int_pt = 0.0;
324 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
325
326 pos = {std::nullopt, _element.getID(),
328 NumLib::interpolateCoordinates<ShapeFunction,
330 _element, _ip_data[ip].N))};
331 vars.liquid_phase_pressure = p_int_pt;
332 vars.capillary_pressure = -p_int_pt;
333 // setting pG to 1 atm
334 // TODO : rewrite equations s.t. p_L = pG-p_cap
335 vars.gas_phase_pressure = 1.0e5;
336
337 double const Sw =
339 .template value<double>(vars, pos, t, dt);
340 vars.liquid_saturation = Sw;
341
342 auto const permeability =
344 medium.property(MaterialPropertyLib::permeability)
345 .value(vars, pos, t, dt));
346
347 double const k_rel =
348 medium
350 relative_permeability)
351 .template value<double>(vars, pos, t, dt);
352 auto const mu =
353 liquid_phase.property(MaterialPropertyLib::viscosity)
354 .template value<double>(vars, pos, t, dt);
355 auto const K_mat_coeff = permeability * (k_rel / mu);
356 cache_vec.col(ip).noalias() =
357 -K_mat_coeff * _ip_data[ip].dNdx * p_nodal_values;
359 {
360 auto const rho_w =
361 liquid_phase
363 .template value<double>(vars, pos, t, dt);
364 auto const& body_force = _process_data.specific_body_force;
365 assert(body_force.size() == GlobalDim);
366 // here it is assumed that the vector body_force is directed
367 // 'downwards'
368 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
369 }
370 }
371
372 return cache;
373 }
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, NumLib::interpolateCoordinates(), 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 261 of file RichardsFlowFEM.h.

263 {
264 auto const& N = _ip_data[integration_point].N;
265
266 // assumes N is stored contiguously in memory
267 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
268 }

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: