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

Detailed Description

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

Definition at line 79 of file RichardsFlowFEM.h.

#include <RichardsFlowFEM.h>

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

Public Member Functions

 LocalAssemblerData (MeshLib::Element const &element, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, 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. More...
 
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::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
void setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, 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_xdot, 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_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, 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_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, 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_dot, 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, double const t, double const dt)
 
void postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
 
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &)
 
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. More...
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType< ShapeFunction, GlobalDim >
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using LocalAssemblerTraits = ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim >
 
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
 
using GlobalDimNodalMatrixType = typename ShapeMatricesType::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
 
IntegrationMethod 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 , typename IntegrationMethod , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 92 of file RichardsFlowFEM.h.

◆ GlobalDimNodalMatrixType

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

Definition at line 88 of file RichardsFlowFEM.h.

◆ GlobalDimVectorType

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

Definition at line 93 of file RichardsFlowFEM.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalAssemblerTraits = ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>
private

Definition at line 84 of file RichardsFlowFEM.h.

◆ NodalMatrixType

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

Definition at line 90 of file RichardsFlowFEM.h.

◆ NodalRowVectorType

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

Definition at line 87 of file RichardsFlowFEM.h.

◆ NodalVectorType

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

Definition at line 91 of file RichardsFlowFEM.h.

◆ ShapeMatrices

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

Definition at line 82 of file RichardsFlowFEM.h.

◆ ShapeMatricesType

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

Definition at line 81 of file RichardsFlowFEM.h.

Constructor & Destructor Documentation

◆ LocalAssemblerData()

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

Definition at line 96 of file RichardsFlowFEM.h.

101  : _element(element),
102  _process_data(process_data),
103  _integration_method(integration_order),
104  _saturation(
105  std::vector<double>(_integration_method.getNumberOfPoints()))
106  {
107  // This assertion is valid only if all nodal d.o.f. use the same shape
108  // matrices.
109  assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
110  (void)local_matrix_size;
111 
112  unsigned const n_integration_points =
113  _integration_method.getNumberOfPoints();
114  _ip_data.reserve(n_integration_points);
115 
116  auto const shape_matrices =
118  GlobalDim>(element, is_axially_symmetric,
120 
121  for (unsigned ip = 0; ip < n_integration_points; ip++)
122  {
123  auto const& sm = shape_matrices[ip];
124  const double integration_factor = sm.integralMeasure * sm.detJ;
125  _ip_data.emplace_back(
126  sm.N, sm.dNdx,
127  _integration_method.getWeightedPoint(ip).getWeight() *
128  integration_factor,
129  sm.N.transpose() * sm.N * integration_factor *
130  _integration_method.getWeightedPoint(ip).getWeight());
131  }
132  }
RichardsFlowProcessData const & _process_data
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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)
const unsigned NUM_NODAL_DOF

References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, NumLib::initShapeMatrices(), and ProcessLib::RichardsFlow::NUM_NODAL_DOF.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, 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 134 of file RichardsFlowFEM.h.

140  {
141  auto const local_matrix_size = local_x.size();
142  // This assertion is valid only if all nodal d.o.f. use the same shape
143  // matrices.
144  assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
145 
146  auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
147  local_M_data, local_matrix_size, local_matrix_size);
148  auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
149  local_K_data, local_matrix_size, local_matrix_size);
150  auto local_b = MathLib::createZeroedVector<NodalVectorType>(
151  local_b_data, local_matrix_size);
152 
153  unsigned const n_integration_points =
154  _integration_method.getNumberOfPoints();
156  pos.setElementID(_element.getID());
157 
158  auto const& medium =
159  *_process_data.media_map->getMedium(_element.getID());
160  auto const& liquid_phase = medium.phase("AqueousLiquid");
162  vars[static_cast<int>(MaterialPropertyLib::Variable::temperature)] =
163  medium
164  .property(
166  .template value<double>(vars, pos, t, dt);
167 
168  for (unsigned ip = 0; ip < n_integration_points; ip++)
169  {
170  pos.setIntegrationPoint(ip);
171  double p_int_pt = 0.0;
172  NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
173 
174  vars[static_cast<int>(
176  vars[static_cast<int>(
178 
179  auto const permeability =
180  MaterialPropertyLib::formEigenTensor<GlobalDim>(
181  medium.property(MaterialPropertyLib::permeability)
182  .value(vars, pos, t, dt));
183 
184  auto const porosity =
186  .template value<double>(vars, pos, t, dt);
187 
188  double const Sw =
189  medium
191  .template value<double>(vars, pos, t, dt);
192  _saturation[ip] = Sw;
193  vars[static_cast<int>(
195 
196  double const dSw_dpc =
197  medium
199  .template dValue<double>(
201  pos, t, dt);
202 
203  auto const drhow_dp =
204  liquid_phase
206  .template dValue<double>(
208  pos, t, dt);
209  auto const storage =
211  .template value<double>(vars, pos, t, dt);
212  double const mass_mat_coeff =
213  storage * Sw + porosity * Sw * drhow_dp - porosity * dSw_dpc;
214 
215  local_M.noalias() += mass_mat_coeff * _ip_data[ip].mass_operator;
216 
217  double const k_rel =
218  medium
221  .template value<double>(vars, pos, t, dt);
222 
223  auto const mu =
224  liquid_phase.property(MaterialPropertyLib::viscosity)
225  .template value<double>(vars, pos, t, dt);
226  local_K.noalias() += _ip_data[ip].dNdx.transpose() * permeability *
227  _ip_data[ip].dNdx *
228  _ip_data[ip].integration_weight * (k_rel / mu);
229 
231  {
232  auto const rho_w =
233  liquid_phase
235  .template value<double>(vars, pos, t, dt);
236  auto const& body_force = _process_data.specific_body_force;
237  assert(body_force.size() == GlobalDim);
238  NodalVectorType gravity_operator =
239  _ip_data[ip].dNdx.transpose() * permeability * body_force *
240  _ip_data[ip].integration_weight;
241  local_b.noalias() += (k_rel / mu) * rho_w * gravity_operator;
242  }
243  }
245  {
246  for (int idx_ml = 0; idx_ml < local_M.cols(); idx_ml++)
247  {
248  double const mass_lump_val = local_M.col(idx_ml).sum();
249  local_M.col(idx_ml).setZero();
250  local_M(idx_ml, idx_ml) = mass_lump_val;
251  }
252  } // end of mass lumping
253  }
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
typename ShapeMatricesType::NodalVectorType NodalVectorType
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
Definition: VariableType.h:108
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
Definition: Interpolation.h:79
std::unique_ptr< MaterialPropertyLib::MaterialSpatialDistributionMap > media_map

References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_saturation, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::density, MeshLib::Element::getID(), ProcessLib::RichardsFlow::RichardsFlowProcessData::has_gravity, ProcessLib::RichardsFlow::RichardsFlowProcessData::has_mass_lumping, MaterialPropertyLib::liquid_saturation, ProcessLib::RichardsFlow::RichardsFlowProcessData::media_map, ProcessLib::RichardsFlow::NUM_NODAL_DOF, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::shapeFunctionInterpolate(), ProcessLib::RichardsFlow::RichardsFlowProcessData::specific_body_force, MaterialPropertyLib::storage, MaterialPropertyLib::temperature, and MaterialPropertyLib::viscosity.

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
std::vector<double> const& ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, 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 274 of file RichardsFlowFEM.h.

279  {
280  // TODO (tf) Temporary value not used by current material models. Need
281  // extension of secondary variable interface
282  double const dt = std::numeric_limits<double>::quiet_NaN();
283 
284  constexpr int process_id = 0; // monolithic scheme.
285  auto const indices =
286  NumLib::getIndices(_element.getID(), *dof_table[process_id]);
287  assert(!indices.empty());
288  auto const local_x = x[process_id]->get(indices);
289 
291  pos.setElementID(_element.getID());
292 
293  auto const& medium =
294  *_process_data.media_map->getMedium(_element.getID());
295  auto const& liquid_phase = medium.phase("AqueousLiquid");
296 
298  vars[static_cast<int>(MaterialPropertyLib::Variable::temperature)] =
299  medium
300  .property(
302  .template value<double>(vars, pos, t, dt);
303 
304  auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
305  &local_x[0], ShapeFunction::NPOINTS);
306 
307  unsigned const n_integration_points =
308  _integration_method.getNumberOfPoints();
309 
310  cache.clear();
311  auto cache_vec = MathLib::createZeroedMatrix<
312  Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
313  cache, GlobalDim, n_integration_points);
314 
315  for (unsigned ip = 0; ip < n_integration_points; ++ip)
316  {
317  double p_int_pt = 0.0;
318  NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
319  vars[static_cast<int>(
321  vars[static_cast<int>(
323 
324  double const Sw =
325  medium
327  .template value<double>(vars, pos, t, dt);
328  vars[static_cast<int>(
330 
331  auto const permeability =
332  MaterialPropertyLib::formEigenTensor<GlobalDim>(
333  medium.property(MaterialPropertyLib::permeability)
334  .value(vars, pos, t, dt));
335 
336  double const k_rel =
337  medium
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  }
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, MaterialPropertyLib::capillary_pressure, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MeshLib::Element::getID(), NumLib::getIndices(), ProcessLib::RichardsFlow::RichardsFlowProcessData::has_gravity, MaterialPropertyLib::liquid_saturation, ProcessLib::RichardsFlow::RichardsFlowProcessData::media_map, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), ProcessLib::RichardsFlow::RichardsFlowProcessData::specific_body_force, MaterialPropertyLib::temperature, and MaterialPropertyLib::viscosity.

◆ getIntPtSaturation()

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

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 255 of file RichardsFlowFEM.h.

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

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

Member Data Documentation

◆ _element

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
MeshLib::Element const& ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element
private

◆ _integration_method

◆ _ip_data

◆ _process_data

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
RichardsFlowProcessData const& ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data
private

◆ _saturation

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

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