OGS
RichardsFlowFEM.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <vector>
7
25
26namespace ProcessLib
27{
28namespace RichardsFlow
29{
30template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType,
31 typename NodalMatrixType>
33{
34 IntegrationPointData(NodalRowVectorType const& N_,
35 GlobalDimNodalMatrixType const& dNdx_,
36 double const& integration_weight_,
37 NodalMatrixType const mass_operator_)
38 : N(N_),
39 dNdx(dNdx_),
40 integration_weight(integration_weight_),
41 mass_operator(mass_operator_)
42 {
43 }
44
45 NodalRowVectorType const N;
46 GlobalDimNodalMatrixType const dNdx;
47 double const integration_weight;
48 NodalMatrixType const mass_operator;
49
51};
52const unsigned NUM_NODAL_DOF = 1;
53
57{
58public:
59 virtual std::vector<double> const& getIntPtSaturation(
60 const double t,
61 std::vector<GlobalVector*> const& x,
62 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
63 std::vector<double>& cache) const = 0;
64
65 virtual std::vector<double> const& getIntPtDarcyVelocity(
66 const double t,
67 std::vector<GlobalVector*> const& x,
68 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
69 std::vector<double>& cache) const = 0;
70};
71
72template <typename ShapeFunction, int GlobalDim>
74{
77
79 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
80
88
89public:
91 MeshLib::Element const& element,
92 std::size_t const local_matrix_size,
93 NumLib::GenericIntegrationMethod const& integration_method,
94 bool is_axially_symmetric,
95 RichardsFlowProcessData const& process_data)
96 : _element(element),
97 _process_data(process_data),
98 _integration_method(integration_method),
100 std::vector<double>(_integration_method.getNumberOfPoints()))
101 {
102 // This assertion is valid only if all nodal d.o.f. use the same shape
103 // matrices.
104 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
105 (void)local_matrix_size;
106
107 unsigned const n_integration_points =
108 _integration_method.getNumberOfPoints();
109 _ip_data.reserve(n_integration_points);
110
111 auto const shape_matrices =
113 GlobalDim>(element, is_axially_symmetric,
115
116 for (unsigned ip = 0; ip < n_integration_points; ip++)
117 {
118 auto const& sm = shape_matrices[ip];
119 const double integration_factor = sm.integralMeasure * sm.detJ;
120 _ip_data.emplace_back(
121 sm.N, sm.dNdx,
122 _integration_method.getWeightedPoint(ip).getWeight() *
123 integration_factor,
124 sm.N.transpose() * sm.N * integration_factor *
125 _integration_method.getWeightedPoint(ip).getWeight());
126 }
127 }
128
129 void assemble(double const t, double const dt,
130 std::vector<double> const& local_x,
131 std::vector<double> const& /*local_x_prev*/,
132 std::vector<double>& local_M_data,
133 std::vector<double>& local_K_data,
134 std::vector<double>& local_b_data) override
135 {
136 auto const local_matrix_size = local_x.size();
137 // This assertion is valid only if all nodal d.o.f. use the same shape
138 // matrices.
139 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
140
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);
147
148 unsigned const n_integration_points =
149 _integration_method.getNumberOfPoints();
151 pos.setElementID(_element.getID());
152
153 auto const& medium =
154 *_process_data.media_map.getMedium(_element.getID());
155 auto const& liquid_phase = medium.phase("AqueousLiquid");
157 vars.temperature =
158 medium
159 .property(
161 .template value<double>(vars, pos, t, dt);
162
163 for (unsigned ip = 0; ip < n_integration_points; ip++)
164 {
165 double p_int_pt = 0.0;
166 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
167
168 pos = {std::nullopt, _element.getID(),
172 _element, _ip_data[ip].N))};
173
174 vars.liquid_phase_pressure = p_int_pt;
175 vars.capillary_pressure = -p_int_pt;
176 // setting pG to 1 atm
177 // TODO : rewrite equations s.t. p_L = pG-p_cap
178 vars.gas_phase_pressure = 1.0e5;
179
180 auto const permeability =
182 medium.property(MaterialPropertyLib::permeability)
183 .value(vars, pos, t, dt));
184
185 auto const porosity =
187 .template value<double>(vars, pos, t, dt);
188
189 double const Sw =
191 .template value<double>(vars, pos, t, dt);
192 _saturation[ip] = Sw;
193 vars.liquid_saturation = Sw;
194
195 double const dSw_dpc =
197 .template dValue<double>(
199 pos, t, dt);
200
201 auto const rho_w =
202 liquid_phase
204 .template value<double>(vars, pos, t, dt);
205
206 auto const drhow_dp =
207 liquid_phase
209 .template dValue<double>(
210 vars,
212 pos, t, dt);
213 auto const storage =
215 .template value<double>(vars, pos, t, dt);
216 double const mass_mat_coeff = storage * Sw +
217 porosity * Sw / rho_w * drhow_dp -
218 porosity * dSw_dpc;
219
220 local_M.noalias() += mass_mat_coeff * _ip_data[ip].mass_operator;
221
222 double const k_rel =
223 medium
225 relative_permeability)
226 .template value<double>(vars, pos, t, dt);
227
228 auto const mu =
229 liquid_phase.property(MaterialPropertyLib::viscosity)
230 .template value<double>(vars, pos, t, dt);
231 local_K.noalias() += _ip_data[ip].dNdx.transpose() * permeability *
232 _ip_data[ip].dNdx *
233 _ip_data[ip].integration_weight * (k_rel / mu);
234
235 if (_process_data.has_gravity)
236 {
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 }
245 if (_process_data.has_mass_lumping)
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 }
255
256 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
257 const unsigned integration_point) const override
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 }
264
265 std::vector<double> const& getIntPtSaturation(
266 const double /*t*/,
267 std::vector<GlobalVector*> const& /*x*/,
268 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
269 std::vector<double>& /*cache*/) const override
270 {
271 assert(!_saturation.empty());
272 return _saturation;
273 }
274
275 std::vector<double> const& getIntPtDarcyVelocity(
276 const double t,
277 std::vector<GlobalVector*> const& x,
278 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
279 std::vector<double>& cache) const override
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
292 pos.setElementID(_element.getID());
293
294 auto const& medium =
295 *_process_data.media_map.getMedium(_element.getID());
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 =
309 _integration_method.getNumberOfPoints();
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
321 pos = {std::nullopt, _element.getID(),
325 _element, _ip_data[ip].N))};
326 vars.liquid_phase_pressure = p_int_pt;
327 vars.capillary_pressure = -p_int_pt;
328 // setting pG to 1 atm
329 // TODO : rewrite equations s.t. p_L = pG-p_cap
330 vars.gas_phase_pressure = 1.0e5;
331
332 double const Sw =
334 .template value<double>(vars, pos, t, dt);
335 vars.liquid_saturation = Sw;
336
337 auto const permeability =
339 medium.property(MaterialPropertyLib::permeability)
340 .value(vars, pos, t, dt));
341
342 double const k_rel =
343 medium
345 relative_permeability)
346 .template value<double>(vars, pos, t, dt);
347 auto const mu =
348 liquid_phase.property(MaterialPropertyLib::viscosity)
349 .template value<double>(vars, pos, t, dt);
350 auto const K_mat_coeff = permeability * (k_rel / mu);
351 cache_vec.col(ip).noalias() =
352 -K_mat_coeff * _ip_data[ip].dNdx * p_nodal_values;
353 if (_process_data.has_gravity)
354 {
355 auto const rho_w =
356 liquid_phase
358 .template value<double>(vars, pos, t, dt);
359 auto const& body_force = _process_data.specific_body_force;
360 assert(body_force.size() == GlobalDim);
361 // here it is assumed that the vector body_force is directed
362 // 'downwards'
363 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
364 }
365 }
366
367 return cache;
368 }
369
370private:
373
375 std::vector<
378 Eigen::aligned_allocator<IntegrationPointData<
381 std::vector<double> _saturation;
382};
383
384} // namespace RichardsFlow
385} // namespace ProcessLib
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
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 > 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
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)
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
IntegrationPointData(NodalRowVectorType const &N_, GlobalDimNodalMatrixType const &dNdx_, double const &integration_weight_, NodalMatrixType const mass_operator_)