OGS
RichardsComponentTransportFEM-impl.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
10
11namespace ProcessLib
12{
14{
15template <typename ShapeFunction, int GlobalDim>
17 MeshLib::Element const& element,
18 std::size_t const local_matrix_size,
19 NumLib::GenericIntegrationMethod const& integration_method,
20 bool is_axially_symmetric,
21 RichardsComponentTransportProcessData const& process_data,
22 ProcessVariable const& transport_process_variable)
23 : _element(element),
24 _process_data(process_data),
25 _integration_method(integration_method),
26 _transport_process_variable(transport_process_variable)
27{
28 // This assertion is valid only if all nodal d.o.f. use the same shape
29 // matrices.
30 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
31 (void)local_matrix_size;
32
33 unsigned const n_integration_points =
34 _integration_method.getNumberOfPoints();
35 _ip_data.reserve(n_integration_points);
36
37 auto const shape_matrices =
39 element, is_axially_symmetric, _integration_method);
40
41 for (unsigned ip = 0; ip < n_integration_points; ip++)
42 {
43 auto const& sm = shape_matrices[ip];
44 double const integration_factor = sm.integralMeasure * sm.detJ;
45 _ip_data.emplace_back(
46 sm.N, sm.dNdx,
47 _integration_method.getWeightedPoint(ip).getWeight() *
48 integration_factor,
49 sm.N.transpose() * sm.N * integration_factor *
50 _integration_method.getWeightedPoint(ip).getWeight());
51 }
52}
53
54template <typename ShapeFunction, int GlobalDim>
56 double const t, double const dt, std::vector<double> const& local_x,
57 std::vector<double> const& /*local_x_prev*/,
58 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
59 std::vector<double>& local_b_data)
60{
61 auto const local_matrix_size = local_x.size();
62 // This assertion is valid only if all nodal d.o.f. use the same shape
63 // matrices.
64 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
65
67 local_M_data, local_matrix_size, local_matrix_size);
69 local_K_data, local_matrix_size, local_matrix_size);
71 local_b_data, local_matrix_size);
72
73 unsigned const n_integration_points =
74 _integration_method.getNumberOfPoints();
75
76 auto p_nodal_values = Eigen::Map<const NodalVectorType>(
77 &local_x[pressure_index], pressure_size);
78
79 auto const& b = _process_data.specific_body_force;
80
82
83 GlobalDimMatrixType const& I(
84 GlobalDimMatrixType::Identity(GlobalDim, GlobalDim));
85
86 // Get material properties
87 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
88 auto const& phase =
90 auto const& component =
91 phase.component(_transport_process_variable.getName());
92
93 auto KCC = local_K.template block<concentration_size, concentration_size>(
95 auto MCC = local_M.template block<concentration_size, concentration_size>(
97 auto Kpp = local_K.template block<pressure_size, pressure_size>(
99 auto Mpp = local_M.template block<pressure_size, pressure_size>(
101 auto Bp = local_b.template block<pressure_size, 1>(pressure_index, 0);
102
103 for (std::size_t ip(0); ip < n_integration_points; ++ip)
104 {
105 auto const& ip_data = _ip_data[ip];
106 auto const& N = ip_data.N;
107 auto const& dNdx = ip_data.dNdx;
108 auto const& w = ip_data.integration_weight;
109
111 std::nullopt, _element.getID(),
114 _element, N))};
115
116 double C_int_pt = 0.0;
117 double p_int_pt = 0.0;
118 // Order matters: First C, then p!
119 NumLib::shapeFunctionInterpolate(local_x, N, C_int_pt, p_int_pt);
120
121 vars.capillary_pressure = -p_int_pt;
123 .template value<double>(vars, pos, t, dt);
124
125 double const dSw_dpc =
127 .template dValue<double>(
129 pos, t, dt);
130
131 vars.concentration = C_int_pt;
132 vars.liquid_phase_pressure = p_int_pt;
133 // setting pG to 1 atm
134 // TODO : rewrite equations s.t. p_L = pG-p_cap
135 vars.gas_phase_pressure = 1.0e5;
136
137 // \todo the argument to getValue() has to be changed for non
138 // constant storage model
139 auto const specific_storage =
141 .template value<double>(vars, pos, t, dt);
142 // \todo the first argument has to be changed for non constant
143 // porosity model
144 auto const porosity =
146 .template value<double>(vars, pos, t, dt);
147
148 auto const retardation_factor =
150 .template value<double>(vars, pos, t, dt);
151
152 auto const solute_dispersivity_transverse =
154 .template value<double>(vars, pos, t, dt);
155 auto const solute_dispersivity_longitudinal =
157 .template value<double>(vars, pos, t, dt);
158
159 // Use the fluid density model to compute the density
160 auto const density = phase[MaterialPropertyLib::PropertyType::density]
161 .template value<double>(vars, pos, t, dt);
162 vars.density = density;
163
164 auto const decay_rate =
166 .template value<double>(vars, pos, t, dt);
167 auto const pore_diffusion_coefficient =
170 .value(vars, pos, t, dt));
171
174 vars, pos, t, dt));
175 vars.liquid_saturation = Sw;
176 auto const k_rel =
178 .template value<double>(vars, pos, t, dt);
179 // Use the viscosity model to compute the viscosity
181 .template value<double>(vars, pos, t, dt);
182 auto const K_times_k_rel_over_mu = K * (k_rel / mu);
183
184 GlobalDimVectorType const velocity =
185 _process_data.has_gravity
186 ? GlobalDimVectorType(-K_times_k_rel_over_mu *
187 (dNdx * p_nodal_values - density * b))
188 : GlobalDimVectorType(-K_times_k_rel_over_mu * dNdx *
189 p_nodal_values);
190
191 double const velocity_magnitude = velocity.norm();
192 GlobalDimMatrixType const hydrodynamic_dispersion =
193 velocity_magnitude != 0.0
195 porosity * pore_diffusion_coefficient +
196 solute_dispersivity_transverse * velocity_magnitude * I +
197 (solute_dispersivity_longitudinal -
198 solute_dispersivity_transverse) /
199 velocity_magnitude * velocity * velocity.transpose())
200 : GlobalDimMatrixType(porosity * pore_diffusion_coefficient +
201 solute_dispersivity_transverse *
202 velocity_magnitude * I);
203
204 // matrix assembly
205 KCC.noalias() +=
206 (dNdx.transpose() * hydrodynamic_dispersion * dNdx +
207 N.transpose() * velocity.transpose() * dNdx +
208 N.transpose() * decay_rate * porosity * retardation_factor * N) *
209 w;
210 MCC.noalias() += w * N.transpose() * porosity * retardation_factor * N;
211 Kpp.noalias() += w * dNdx.transpose() * K_times_k_rel_over_mu * dNdx;
212 // \TODO Extend to pressure dependent density.
213 double const drhow_dp(0.0);
214 Mpp.noalias() += (specific_storage * Sw + porosity * Sw * drhow_dp -
215 porosity * dSw_dpc) *
216 ip_data.mass_operator;
217
218 if (_process_data.has_gravity)
219 {
220 Bp += w * density * dNdx.transpose() * K_times_k_rel_over_mu * b;
221 }
222 /* with Oberbeck-Boussing assumption density difference only exists
223 * in buoyancy effects */
224 }
225}
226
227template <typename ShapeFunction, int GlobalDim>
228std::vector<double> const&
230 const double t,
231 std::vector<GlobalVector*> const& x,
232 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
233 std::vector<double>& cache) const
234{
235 unsigned const n_integration_points =
236 _integration_method.getNumberOfPoints();
237
238 constexpr int process_id = 0; // monolithic scheme
239 auto const indices =
240 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
241 assert(!indices.empty());
242 auto const local_x = x[process_id]->get(indices);
243
244 cache.clear();
245 auto cache_mat = MathLib::createZeroedMatrix<
246 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
247 cache, GlobalDim, n_integration_points);
248
250
251 // Get material properties
252 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
253 auto const& phase =
255
256 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
257 &local_x[ShapeFunction::NPOINTS], ShapeFunction::NPOINTS);
258
259 for (unsigned ip = 0; ip < n_integration_points; ++ip)
260 {
261 auto const& ip_data = _ip_data[ip];
262 auto const& N = ip_data.N;
263 auto const& dNdx = ip_data.dNdx;
264
266 std::nullopt, _element.getID(),
269 _element, N))};
270
271 auto const dt = std::numeric_limits<double>::quiet_NaN();
274 vars, pos, t, dt));
276 .template value<double>(vars, pos, t, dt);
277
278 double C_int_pt = 0.0;
279 double p_int_pt = 0.0;
280 NumLib::shapeFunctionInterpolate(local_x, N, C_int_pt, p_int_pt);
281
282 // saturation
283 vars.capillary_pressure = -p_int_pt;
285 .template value<double>(vars, pos, t, dt);
286
287 vars.liquid_saturation = Sw;
288 auto const k_rel =
290 .template value<double>(vars, pos, t, dt);
291
292 cache_mat.col(ip).noalias() = -dNdx * p_nodal_values;
293 if (_process_data.has_gravity)
294 {
295 vars.concentration = C_int_pt;
296 vars.liquid_phase_pressure = p_int_pt;
297 auto const rho_w = phase[MaterialPropertyLib::PropertyType::density]
298 .template value<double>(vars, pos, t, dt);
299 auto const b = _process_data.specific_body_force;
300 // here it is assumed that the vector b is directed 'downwards'
301 cache_mat.col(ip).noalias() += rho_w * b;
302 }
303 cache_mat.col(ip).noalias() = k_rel / mu * (K * cache_mat.col(ip));
304 }
305 return cache;
306}
307
308template <typename ShapeFunction, int GlobalDim>
309Eigen::Map<const Eigen::RowVectorXd>
311 const unsigned integration_point) const
312{
313 auto const& N = _ip_data[integration_point].N;
314
315 // assumes N is stored contiguously in memory
316 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
317}
318
319template <typename ShapeFunction, int GlobalDim>
320std::vector<double> const&
322 const double t,
323 std::vector<GlobalVector*> const& x,
324 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
325 std::vector<double>& cache) const
326{
328
329 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
330
331 unsigned const n_integration_points =
332 _integration_method.getNumberOfPoints();
333
334 constexpr int process_id = 0; // monolithic scheme
335 auto const indices =
336 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
337 assert(!indices.empty());
338 auto const local_x = x[process_id]->get(indices);
339
340 cache.clear();
341 auto cache_vec = MathLib::createZeroedVector<
342 Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor>>(
343 cache, n_integration_points);
344
345 for (unsigned ip = 0; ip < n_integration_points; ++ip)
346 {
347 auto const& ip_data = _ip_data[ip];
348 auto const& N = ip_data.N;
349
351 std::nullopt, _element.getID(),
354 _element, N))};
355
356 double C_int_pt = 0.0;
357 double p_int_pt = 0.0;
358 NumLib::shapeFunctionInterpolate(local_x, N, C_int_pt, p_int_pt);
359
360 // saturation
361 vars.capillary_pressure = -p_int_pt;
362 auto const dt = std::numeric_limits<double>::quiet_NaN();
364 .template value<double>(vars, pos, t, dt);
365 cache_vec[ip] = Sw;
366 }
367
368 return cache;
369}
370
371} // namespace RichardsComponentTransport
372} // namespace ProcessLib
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
LocalAssemblerData(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, RichardsComponentTransportProcessData const &process_data, ProcessVariable const &transport_process_variable)
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 override
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
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
@ longitudinal_dispersivity
used to compute the hydrodynamic dispersion tensor.
@ transversal_dispersivity
used to compute the hydrodynamic dispersion tensor.
@ retardation_factor
specify retardation factor used in component transport process.
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)