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::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:
ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits

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),
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.
113
114 unsigned const n_integration_points =
115 _integration_method.getNumberOfPoints();
117
118 auto const shape_matrices =
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,
129 _integration_method.getWeightedPoint(ip).getWeight() *
131 sm.N.transpose() * sm.N * integration_factor *
132 _integration_method.getWeightedPoint(ip).getWeight());
133 }
134 }
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
RichardsFlowProcessData const & _process_data

References _element, _integration_method, _ip_data, _process_data, _saturation, 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.
147
154
155 unsigned const n_integration_points =
156 _integration_method.getNumberOfPoints();
158 pos.setElementID(_element.getID());
159
160 auto const& medium =
161 *_process_data.media_map.getMedium(_element.getID());
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;
174
175 pos = {std::nullopt, _element.getID(),
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 =
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 rho_w =
211 .template value<double>(vars, pos, t, dt);
212
213 auto const drhow_dp =
216 .template dValue<double>(
217 vars,
219 pos, t, dt);
220 auto const storage =
222 .template value<double>(vars, pos, t, dt);
223 double const mass_mat_coeff = storage * Sw +
224 porosity * Sw / rho_w * drhow_dp -
226
227 local_M.noalias() += mass_mat_coeff * _ip_data[ip].mass_operator;
228
229 double const k_rel =
230 medium
233 .template value<double>(vars, pos, t, dt);
234
235 auto const mu =
237 .template value<double>(vars, pos, t, dt);
238 local_K.noalias() += _ip_data[ip].dNdx.transpose() * permeability *
239 _ip_data[ip].dNdx *
240 _ip_data[ip].integration_weight * (k_rel / mu);
241
242 if (_process_data.has_gravity)
243 {
244 auto const& body_force = _process_data.specific_body_force;
245 assert(body_force.size() == GlobalDim);
247 _ip_data[ip].dNdx.transpose() * permeability * body_force *
248 _ip_data[ip].integration_weight;
249 local_b.noalias() += (k_rel / mu) * rho_w * gravity_operator;
250 }
251 }
252 if (_process_data.has_mass_lumping)
253 {
254 for (int idx_ml = 0; idx_ml < local_M.cols(); idx_ml++)
255 {
256 double const mass_lump_val = local_M.col(idx_ml).sum();
257 local_M.col(idx_ml).setZero();
259 }
260 } // end of mass lumping
261 }
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 &)

References _element, _integration_method, _ip_data, _process_data, _saturation, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::VariableArray::capillary_pressure, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::gas_phase_pressure, NumLib::interpolateCoordinates(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, ProcessLib::RichardsFlow::NUM_NODAL_DOF, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), 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 282 of file RichardsFlowFEM.h.

287 {
288 // TODO (tf) Temporary value not used by current material models. Need
289 // extension of secondary variable interface
291
292 constexpr int process_id = 0; // monolithic scheme.
293 auto const indices =
295 assert(!indices.empty());
296 auto const local_x = x[process_id]->get(indices);
297
299 pos.setElementID(_element.getID());
300
301 auto const& medium =
302 *_process_data.media_map.getMedium(_element.getID());
303 auto const& liquid_phase = medium.phase("AqueousLiquid");
304
306 vars.temperature =
307 medium
308 .property(
310 .template value<double>(vars, pos, t, dt);
311
314
315 unsigned const n_integration_points =
316 _integration_method.getNumberOfPoints();
317
318 cache.clear();
322
323 for (unsigned ip = 0; ip < n_integration_points; ++ip)
324 {
325 double p_int_pt = 0.0;
327
328 pos = {std::nullopt, _element.getID(),
332 _element, _ip_data[ip].N))};
333 vars.liquid_phase_pressure = p_int_pt;
334 vars.capillary_pressure = -p_int_pt;
335 // setting pG to 1 atm
336 // TODO : rewrite equations s.t. p_L = pG-p_cap
337 vars.gas_phase_pressure = 1.0e5;
338
339 double const Sw =
341 .template value<double>(vars, pos, t, dt);
342 vars.liquid_saturation = Sw;
343
344 auto const permeability =
347 .value(vars, pos, t, dt));
348
349 double const k_rel =
350 medium
353 .template value<double>(vars, pos, t, dt);
354 auto const mu =
356 .template value<double>(vars, pos, t, dt);
357 auto const K_mat_coeff = permeability * (k_rel / mu);
358 cache_vec.col(ip).noalias() =
360 if (_process_data.has_gravity)
361 {
362 auto const rho_w =
365 .template value<double>(vars, pos, t, dt);
366 auto const& body_force = _process_data.specific_body_force;
367 assert(body_force.size() == GlobalDim);
368 // here it is assumed that the vector body_force is directed
369 // 'downwards'
370 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
371 }
372 }
373
374 return cache;
375 }
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References _element, _integration_method, _ip_data, _process_data, MaterialPropertyLib::VariableArray::capillary_pressure, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::gas_phase_pressure, NumLib::getIndices(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::permeability, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), 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

Implements ProcessLib::RichardsFlow::RichardsFlowLocalAssemblerInterface.

Definition at line 272 of file RichardsFlowFEM.h.

277 {
278 assert(!_saturation.empty());
279 return _saturation;
280 }

References _saturation.

◆ 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 263 of file RichardsFlowFEM.h.

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

References _ip_data.

Member Data Documentation

◆ _element

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

Definition at line 378 of file RichardsFlowFEM.h.

Referenced by LocalAssemblerData(), assemble(), and getIntPtDarcyVelocity().

◆ _integration_method

template<typename ShapeFunction, int GlobalDim>
NumLib::GenericIntegrationMethod const& ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method
private

Definition at line 381 of file RichardsFlowFEM.h.

Referenced by LocalAssemblerData(), assemble(), and getIntPtDarcyVelocity().

◆ _ip_data

◆ _process_data

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

Definition at line 379 of file RichardsFlowFEM.h.

Referenced by LocalAssemblerData(), assemble(), and getIntPtDarcyVelocity().

◆ _saturation

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

Definition at line 388 of file RichardsFlowFEM.h.

Referenced by LocalAssemblerData(), assemble(), and getIntPtSaturation().


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