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 =
158 vars.temperature =
159 medium
160 .property(
162 .template value<double>(vars, pos, t, dt);
163
164 for (unsigned ip = 0; ip < n_integration_points; ip++)
165 {
166 double p_int_pt = 0.0;
167 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
168
169 pos = {std::nullopt, _element.getID(),
173 _element, _ip_data[ip].N))};
174
175 vars.liquid_phase_pressure = p_int_pt;
176 vars.capillary_pressure = -p_int_pt;
177 // setting pG to 1 atm
178 // TODO : rewrite equations s.t. p_L = pG-p_cap
179 vars.gas_phase_pressure = 1.0e5;
180
181 auto const permeability =
183 medium.property(MaterialPropertyLib::permeability)
184 .value(vars, pos, t, dt));
185
186 auto const porosity =
188 .template value<double>(vars, pos, t, dt);
189
190 double const Sw =
192 .template value<double>(vars, pos, t, dt);
193 _saturation[ip] = Sw;
194 vars.liquid_saturation = Sw;
195
196 double const dSw_dpc =
198 .template dValue<double>(
200 pos, t, dt);
201
202 auto const rho_w =
203 liquid_phase
205 .template value<double>(vars, pos, t, dt);
206
207 auto const drhow_dp =
208 liquid_phase
210 .template dValue<double>(
211 vars,
213 pos, t, dt);
214 auto const storage =
216 .template value<double>(vars, pos, t, dt);
217 double const mass_mat_coeff = storage * Sw +
218 porosity * Sw / rho_w * drhow_dp -
219 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
236 if (_process_data.has_gravity)
237 {
238 auto const& body_force = _process_data.specific_body_force;
239 assert(body_force.size() == GlobalDim);
240 NodalVectorType gravity_operator =
241 _ip_data[ip].dNdx.transpose() * permeability * body_force *
242 _ip_data[ip].integration_weight;
243 local_b.noalias() += (k_rel / mu) * rho_w * gravity_operator;
244 }
245 }
246 if (_process_data.has_mass_lumping)
247 {
248 for (int idx_ml = 0; idx_ml < local_M.cols(); idx_ml++)
249 {
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;
253 }
254 } // end of mass lumping
255 }
256
257 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
258 const unsigned integration_point) const override
259 {
260 auto const& N = _ip_data[integration_point].N;
261
262 // assumes N is stored contiguously in memory
263 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
264 }
265
266 std::vector<double> const& getIntPtSaturation(
267 const double /*t*/,
268 std::vector<GlobalVector*> const& /*x*/,
269 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
270 std::vector<double>& /*cache*/) const override
271 {
272 assert(!_saturation.empty());
273 return _saturation;
274 }
275
276 std::vector<double> const& getIntPtDarcyVelocity(
277 const double t,
278 std::vector<GlobalVector*> const& x,
279 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
280 std::vector<double>& cache) const override
281 {
282 // TODO (tf) Temporary value not used by current material models. Need
283 // extension of secondary variable interface
284 double const dt = std::numeric_limits<double>::quiet_NaN();
285
286 constexpr int process_id = 0; // monolithic scheme.
287 auto const indices =
288 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
289 assert(!indices.empty());
290 auto const local_x = x[process_id]->get(indices);
291
293 pos.setElementID(_element.getID());
294
295 auto const& medium =
296 *_process_data.media_map.getMedium(_element.getID());
297 auto const& liquid_phase =
299
301 vars.temperature =
302 medium
303 .property(
305 .template value<double>(vars, pos, t, dt);
306
307 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
308 &local_x[0], ShapeFunction::NPOINTS);
309
310 unsigned const n_integration_points =
311 _integration_method.getNumberOfPoints();
312
313 cache.clear();
314 auto cache_vec = MathLib::createZeroedMatrix<
315 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
316 cache, GlobalDim, n_integration_points);
317
318 for (unsigned ip = 0; ip < n_integration_points; ++ip)
319 {
320 double p_int_pt = 0.0;
321 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
322
323 pos = {std::nullopt, _element.getID(),
327 _element, _ip_data[ip].N))};
328 vars.liquid_phase_pressure = p_int_pt;
329 vars.capillary_pressure = -p_int_pt;
330 // setting pG to 1 atm
331 // TODO : rewrite equations s.t. p_L = pG-p_cap
332 vars.gas_phase_pressure = 1.0e5;
333
334 double const Sw =
336 .template value<double>(vars, pos, t, dt);
337 vars.liquid_saturation = Sw;
338
339 auto const permeability =
341 medium.property(MaterialPropertyLib::permeability)
342 .value(vars, pos, t, dt));
343
344 double const k_rel =
345 medium
347 relative_permeability)
348 .template value<double>(vars, pos, t, dt);
349 auto const mu =
350 liquid_phase.property(MaterialPropertyLib::viscosity)
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;
355 if (_process_data.has_gravity)
356 {
357 auto const rho_w =
358 liquid_phase
360 .template value<double>(vars, pos, t, dt);
361 auto const& body_force = _process_data.specific_body_force;
362 assert(body_force.size() == GlobalDim);
363 // here it is assumed that the vector body_force is directed
364 // 'downwards'
365 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
366 }
367 }
368
369 return cache;
370 }
371
372private:
375
377 std::vector<
380 Eigen::aligned_allocator<IntegrationPointData<
383 std::vector<double> _saturation;
384};
385
386} // namespace RichardsFlow
387} // 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_)