OGS
RichardsFlowFEM.h
Go to the documentation of this file.
1
10
11#pragma once
12
13#include <vector>
14
32
33namespace ProcessLib
34{
35namespace RichardsFlow
36{
37template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType,
38 typename NodalMatrixType>
40{
41 IntegrationPointData(NodalRowVectorType const& N_,
42 GlobalDimNodalMatrixType const& dNdx_,
43 double const& integration_weight_,
44 NodalMatrixType const mass_operator_)
45 : N(N_),
46 dNdx(dNdx_),
47 integration_weight(integration_weight_),
48 mass_operator(mass_operator_)
49 {
50 }
51
52 NodalRowVectorType const N;
53 GlobalDimNodalMatrixType const dNdx;
54 double const integration_weight;
55 NodalMatrixType const mass_operator;
56
58};
59const unsigned NUM_NODAL_DOF = 1;
60
64{
65public:
66 virtual std::vector<double> const& getIntPtSaturation(
67 const double t,
68 std::vector<GlobalVector*> const& x,
69 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
70 std::vector<double>& cache) const = 0;
71
72 virtual std::vector<double> const& getIntPtDarcyVelocity(
73 const double t,
74 std::vector<GlobalVector*> const& x,
75 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
76 std::vector<double>& cache) const = 0;
77};
78
79template <typename ShapeFunction, int GlobalDim>
81{
84
86 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
87
95
96public:
98 MeshLib::Element const& element,
99 std::size_t const local_matrix_size,
100 NumLib::GenericIntegrationMethod const& integration_method,
101 bool is_axially_symmetric,
102 RichardsFlowProcessData const& process_data)
103 : _element(element),
104 _process_data(process_data),
105 _integration_method(integration_method),
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.
111 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
112 (void)local_matrix_size;
113
114 unsigned const n_integration_points =
115 _integration_method.getNumberOfPoints();
116 _ip_data.reserve(n_integration_points);
117
118 auto const shape_matrices =
120 GlobalDim>(element, is_axially_symmetric,
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() *
130 integration_factor,
131 sm.N.transpose() * sm.N * integration_factor *
132 _integration_method.getWeightedPoint(ip).getWeight());
133 }
134 }
135
136 void assemble(double const t, double const dt,
137 std::vector<double> const& local_x,
138 std::vector<double> const& /*local_x_prev*/,
139 std::vector<double>& local_M_data,
140 std::vector<double>& local_K_data,
141 std::vector<double>& local_b_data) override
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.
146 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
147
149 local_M_data, local_matrix_size, local_matrix_size);
151 local_K_data, local_matrix_size, local_matrix_size);
153 local_b_data, local_matrix_size);
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;
173 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
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 =
189 medium.property(MaterialPropertyLib::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 =
209 liquid_phase
211 .template value<double>(vars, pos, t, dt);
212
213 auto const drhow_dp =
214 liquid_phase
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 -
225 porosity * dSw_dpc;
226
227 local_M.noalias() += mass_mat_coeff * _ip_data[ip].mass_operator;
228
229 double const k_rel =
230 medium
232 relative_permeability)
233 .template value<double>(vars, pos, t, dt);
234
235 auto const mu =
236 liquid_phase.property(MaterialPropertyLib::viscosity)
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);
246 NodalVectorType gravity_operator =
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();
258 local_M(idx_ml, idx_ml) = mass_lump_val;
259 }
260 } // end of mass lumping
261 }
262
263 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
264 const unsigned integration_point) const override
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 }
271
272 std::vector<double> const& getIntPtSaturation(
273 const double /*t*/,
274 std::vector<GlobalVector*> const& /*x*/,
275 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
276 std::vector<double>& /*cache*/) const override
277 {
278 assert(!_saturation.empty());
279 return _saturation;
280 }
281
282 std::vector<double> const& getIntPtDarcyVelocity(
283 const double t,
284 std::vector<GlobalVector*> const& x,
285 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
286 std::vector<double>& cache) const override
287 {
288 // TODO (tf) Temporary value not used by current material models. Need
289 // extension of secondary variable interface
290 double const dt = std::numeric_limits<double>::quiet_NaN();
291
292 constexpr int process_id = 0; // monolithic scheme.
293 auto const indices =
294 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
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
312 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
313 &local_x[0], ShapeFunction::NPOINTS);
314
315 unsigned const n_integration_points =
316 _integration_method.getNumberOfPoints();
317
318 cache.clear();
319 auto cache_vec = MathLib::createZeroedMatrix<
320 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
321 cache, GlobalDim, n_integration_points);
322
323 for (unsigned ip = 0; ip < n_integration_points; ++ip)
324 {
325 double p_int_pt = 0.0;
326 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
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 =
346 medium.property(MaterialPropertyLib::permeability)
347 .value(vars, pos, t, dt));
348
349 double const k_rel =
350 medium
352 relative_permeability)
353 .template value<double>(vars, pos, t, dt);
354 auto const mu =
355 liquid_phase.property(MaterialPropertyLib::viscosity)
356 .template value<double>(vars, pos, t, dt);
357 auto const K_mat_coeff = permeability * (k_rel / mu);
358 cache_vec.col(ip).noalias() =
359 -K_mat_coeff * _ip_data[ip].dNdx * p_nodal_values;
360 if (_process_data.has_gravity)
361 {
362 auto const rho_w =
363 liquid_phase
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 }
376
377private:
380
382 std::vector<
385 Eigen::aligned_allocator<IntegrationPointData<
388 std::vector<double> _saturation;
389};
390
391} // namespace RichardsFlow
392} // namespace ProcessLib
Definition of the PiecewiseLinearInterpolation class.
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
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_)