OGS
RichardsFlowFEM.h
Go to the documentation of this file.
1
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 =
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,
130 integration_factor,
131 sm.N.transpose() * sm.N * integration_factor *
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
148 auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
149 local_M_data, local_matrix_size, local_matrix_size);
150 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
151 local_K_data, local_matrix_size, local_matrix_size);
152 auto local_b = MathLib::createZeroedVector<NodalVectorType>(
153 local_b_data, local_matrix_size);
154
155 unsigned const n_integration_points =
159
160 auto const& medium =
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 pos.setIntegrationPoint(ip);
173 double p_int_pt = 0.0;
174 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, p_int_pt);
175
176 vars.liquid_phase_pressure = p_int_pt;
177 vars.capillary_pressure = -p_int_pt;
178 // setting pG to 1 atm
179 // TODO : rewrite equations s.t. p_L = pG-p_cap
180 vars.gas_phase_pressure = 1.0e5;
181
182 auto const permeability =
183 MaterialPropertyLib::formEigenTensor<GlobalDim>(
184 medium.property(MaterialPropertyLib::permeability)
185 .value(vars, pos, t, dt));
186
187 auto const porosity =
189 .template value<double>(vars, pos, t, dt);
190
191 double const Sw =
193 .template value<double>(vars, pos, t, dt);
194 _saturation[ip] = Sw;
195 vars.liquid_saturation = Sw;
196
197 double const dSw_dpc =
199 .template dValue<double>(
201 pos, t, dt);
202
203 auto const drhow_dp =
204 liquid_phase
206 .template dValue<double>(
207 vars,
209 pos, t, dt);
210 auto const storage =
212 .template value<double>(vars, pos, t, dt);
213 double const mass_mat_coeff =
214 storage * Sw + porosity * Sw * drhow_dp - porosity * dSw_dpc;
215
216 local_M.noalias() += mass_mat_coeff * _ip_data[ip].mass_operator;
217
218 double const k_rel =
219 medium
221 relative_permeability)
222 .template value<double>(vars, pos, t, dt);
223
224 auto const mu =
225 liquid_phase.property(MaterialPropertyLib::viscosity)
226 .template value<double>(vars, pos, t, dt);
227 local_K.noalias() += _ip_data[ip].dNdx.transpose() * permeability *
228 _ip_data[ip].dNdx *
229 _ip_data[ip].integration_weight * (k_rel / mu);
230
232 {
233 auto const rho_w =
234 liquid_phase
236 .template value<double>(vars, pos, t, dt);
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 }
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
293
294 auto const& medium =
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 =
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 vars.liquid_phase_pressure = p_int_pt;
321 vars.capillary_pressure = -p_int_pt;
322 // setting pG to 1 atm
323 // TODO : rewrite equations s.t. p_L = pG-p_cap
324 vars.gas_phase_pressure = 1.0e5;
325
326 double const Sw =
328 .template value<double>(vars, pos, t, dt);
329 vars.liquid_saturation = Sw;
330
331 auto const permeability =
332 MaterialPropertyLib::formEigenTensor<GlobalDim>(
333 medium.property(MaterialPropertyLib::permeability)
334 .value(vars, pos, t, dt));
335
336 double const k_rel =
337 medium
339 relative_permeability)
340 .template value<double>(vars, pos, t, dt);
341 auto const mu =
342 liquid_phase.property(MaterialPropertyLib::viscosity)
343 .template value<double>(vars, pos, t, dt);
344 auto const K_mat_coeff = permeability * (k_rel / mu);
345 cache_vec.col(ip).noalias() =
346 -K_mat_coeff * _ip_data[ip].dNdx * p_nodal_values;
348 {
349 auto const rho_w =
350 liquid_phase
352 .template value<double>(vars, pos, t, dt);
353 auto const& body_force = _process_data.specific_body_force;
354 assert(body_force.size() == GlobalDim);
355 // here it is assumed that the vector body_force is directed
356 // 'downwards'
357 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
358 }
359 }
360
361 return cache;
362 }
363
364private:
367
369 std::vector<
372 Eigen::aligned_allocator<IntegrationPointData<
375 std::vector<double> _saturation;
376};
377
378} // namespace RichardsFlow
379} // namespace ProcessLib
Definition of the PiecewiseLinearInterpolation class.
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
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
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::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)
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_)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map