30template <
typename NodalRowVectorType,
typename GlobalDimNodalMatrixType,
31 typename NodalMatrixType>
35 GlobalDimNodalMatrixType
const& dNdx_,
36 double const& integration_weight_,
37 NodalMatrixType
const mass_operator_)
45 NodalRowVectorType
const N;
46 GlobalDimNodalMatrixType
const dNdx;
61 std::vector<GlobalVector*>
const& x,
62 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
63 std::vector<double>& cache)
const = 0;
67 std::vector<GlobalVector*>
const& x,
68 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
69 std::vector<double>& cache)
const = 0;
72template <
typename ShapeFunction,
int GlobalDim>
92 std::size_t
const local_matrix_size,
94 bool is_axially_symmetric,
104 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
105 (void)local_matrix_size;
107 unsigned const n_integration_points =
109 _ip_data.reserve(n_integration_points);
111 auto const shape_matrices =
113 GlobalDim>(element, is_axially_symmetric,
116 for (
unsigned ip = 0; ip < n_integration_points; ip++)
118 auto const& sm = shape_matrices[ip];
119 const double integration_factor = sm.integralMeasure * sm.detJ;
124 sm.N.transpose() * sm.N * integration_factor *
130 std::vector<double>
const& local_x,
131 std::vector<double>
const& ,
132 std::vector<double>& local_M_data,
133 std::vector<double>& local_K_data,
134 std::vector<double>& local_b_data)
override
136 auto const local_matrix_size = local_x.size();
139 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
142 local_M_data, local_matrix_size, local_matrix_size);
144 local_K_data, local_matrix_size, local_matrix_size);
146 local_b_data, local_matrix_size);
148 unsigned const n_integration_points =
155 auto const& liquid_phase =
162 .template value<double>(vars, pos, t, dt);
164 for (
unsigned ip = 0; ip < n_integration_points; ip++)
166 double p_int_pt = 0.0;
169 pos = {std::nullopt,
_element.getID(),
181 auto const permeability =
184 .value(vars, pos, t, dt));
186 auto const porosity =
188 .template value<double>(vars, pos, t, dt);
192 .template value<double>(vars, pos, t, dt);
196 double const dSw_dpc =
198 .template dValue<double>(
205 .template value<double>(vars, pos, t, dt);
207 auto const drhow_dp =
210 .template dValue<double>(
216 .template value<double>(vars, pos, t, dt);
217 double const mass_mat_coeff = storage * Sw +
218 porosity * Sw / rho_w * drhow_dp -
221 local_M.noalias() += mass_mat_coeff *
_ip_data[ip].mass_operator;
226 relative_permeability)
227 .template value<double>(vars, pos, t, dt);
231 .template value<double>(vars, pos, t, dt);
232 local_K.noalias() +=
_ip_data[ip].dNdx.transpose() * permeability *
234 _ip_data[ip].integration_weight * (k_rel / mu);
239 assert(body_force.size() == GlobalDim);
241 _ip_data[ip].dNdx.transpose() * permeability * body_force *
243 local_b.noalias() += (k_rel / mu) * rho_w * gravity_operator;
248 for (
int idx_ml = 0; idx_ml < local_M.cols(); idx_ml++)
250 double const mass_lump_val = local_M.col(idx_ml).sum();
251 local_M.col(idx_ml).setZero();
252 local_M(idx_ml, idx_ml) = mass_lump_val;
258 const unsigned integration_point)
const override
260 auto const& N =
_ip_data[integration_point].N;
263 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
268 std::vector<GlobalVector*>
const& ,
269 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
270 std::vector<double>& )
const override
278 std::vector<GlobalVector*>
const& x,
279 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
280 std::vector<double>& cache)
const override
284 double const dt = std::numeric_limits<double>::quiet_NaN();
286 constexpr int process_id = 0;
289 assert(!indices.empty());
290 auto const local_x = x[process_id]->get(indices);
297 auto const& liquid_phase =
305 .template value<double>(vars, pos, t, dt);
307 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
308 &local_x[0], ShapeFunction::NPOINTS);
310 unsigned const n_integration_points =
315 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
316 cache, GlobalDim, n_integration_points);
318 for (
unsigned ip = 0; ip < n_integration_points; ++ip)
320 double p_int_pt = 0.0;
323 pos = {std::nullopt,
_element.getID(),
336 .template value<double>(vars, pos, t, dt);
339 auto const permeability =
342 .value(vars, pos, t, dt));
347 relative_permeability)
348 .template value<double>(vars, pos, t, dt);
351 .template value<double>(vars, pos, t, dt);
352 auto const K_mat_coeff = permeability * (k_rel / mu);
353 cache_vec.col(ip).noalias() =
354 -K_mat_coeff *
_ip_data[ip].dNdx * p_nodal_values;
360 .template value<double>(vars, pos, t, dt);
362 assert(body_force.size() == GlobalDim);
365 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
double gas_phase_pressure
double capillary_pressure
double liquid_phase_pressure
void setElementID(std::size_t element_id)
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< double > _saturation
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
std::vector< double > const & getIntPtSaturation(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
RichardsFlowProcessData const & _process_data
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)
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
MeshLib::Element const & _element
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim > LocalAssemblerTraits
typename ShapeMatricesType::NodalVectorType NodalVectorType
virtual std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual 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 =0
constexpr 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::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
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)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
const unsigned NUM_NODAL_DOF
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
GlobalDimNodalMatrixType const dNdx
double const integration_weight
NodalRowVectorType const N
IntegrationPointData(NodalRowVectorType const &N_, GlobalDimNodalMatrixType const &dNdx_, double const &integration_weight_, NodalMatrixType const mass_operator_)
NodalMatrixType const mass_operator