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