OGS
HydroMechanicsLocalAssemblerFracture-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
4#pragma once
5
6#include <variant>
7
15
16namespace ProcessLib
17{
18namespace LIE
19{
20namespace HydroMechanics
21{
22namespace MPL = MaterialPropertyLib;
23
24template <int DisplacementDim, typename RotationMatrix>
25Eigen::Matrix<double, DisplacementDim, DisplacementDim> createRotatedTensor(
26 RotationMatrix const& R, double const value)
27{
28 using M = Eigen::Matrix<double, DisplacementDim, DisplacementDim>;
29 M tensor = M::Zero();
30 tensor.diagonal().head(DisplacementDim - 1).setConstant(value);
31 return (R.transpose() * tensor * R).eval();
32}
33
34template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
35 int DisplacementDim>
36HydroMechanicsLocalAssemblerFracture<ShapeFunctionDisplacement,
37 ShapeFunctionPressure, DisplacementDim>::
38 HydroMechanicsLocalAssemblerFracture(
39 MeshLib::Element const& e,
40 std::size_t const /*local_matrix_size*/,
41 std::vector<unsigned> const& dofIndex_to_localIndex,
42 NumLib::GenericIntegrationMethod const& integration_method,
43 bool const is_axially_symmetric,
46 e, is_axially_symmetric, integration_method,
47 ShapeFunctionDisplacement::NPOINTS * DisplacementDim +
48 ShapeFunctionPressure::NPOINTS,
49 dofIndex_to_localIndex),
50 _process_data(process_data)
51{
52 assert(e.getDimension() == DisplacementDim - 1);
53
54 unsigned const n_integration_points =
55 integration_method.getNumberOfPoints();
56
57 _ip_data.reserve(n_integration_points);
58 _secondary_data.N.resize(n_integration_points);
59
60 auto const shape_matrices_u =
61 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
63 DisplacementDim>(e, is_axially_symmetric,
64 integration_method);
65
66 auto const shape_matrices_p =
67 NumLib::initShapeMatrices<ShapeFunctionPressure,
68 ShapeMatricesTypePressure, DisplacementDim>(
69 e, is_axially_symmetric, integration_method);
70
71 auto mat_id = (*_process_data.mesh_prop_materialIDs)[e.getID()];
72 auto frac_id = _process_data.map_materialID_to_fractureID[mat_id];
73 _fracture_property = &_process_data.fracture_properties[frac_id];
74
75 // Get element nodes for aperture0 interpolation from nodes to integration
76 // point. The aperture0 parameter is time-independent.
78 aperture0_node_values =
79 _fracture_property->aperture0.getNodalValuesOnElement(
80 e, /*time independent*/ 0);
81
82 for (unsigned ip = 0; ip < n_integration_points; ip++)
83 {
84 _ip_data.emplace_back(*_process_data.fracture_model);
85 auto const& sm_u = shape_matrices_u[ip];
86 auto const& sm_p = shape_matrices_p[ip];
88 std::nullopt, _element.getID(),
90 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
92 _element, sm_u.N))};
93
94 auto& ip_data = _ip_data[ip];
95 ip_data.integration_weight =
96 sm_u.detJ * sm_u.integralMeasure *
97 integration_method.getWeightedPoint(ip).getWeight();
98
99 ip_data.H_u.setZero(
100 DisplacementDim,
101 ShapeFunctionDisplacement::NPOINTS * DisplacementDim);
103 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
105 HMatrixType>(sm_u.N, ip_data.H_u);
106 ip_data.N_p = sm_p.N;
107 ip_data.dNdx_p = sm_p.dNdx;
108
109 _secondary_data.N[ip] = sm_u.N;
110
111 // Initialize current time step values
112 ip_data.w.setZero(DisplacementDim);
113 ip_data.sigma_eff.setZero(DisplacementDim);
114
115 // Previous time step values are not initialized and are set later.
116 ip_data.w_prev.resize(DisplacementDim);
117 ip_data.sigma_eff_prev.resize(DisplacementDim);
118
119 ip_data.C.resize(DisplacementDim, DisplacementDim);
120
121 ip_data.aperture0 = aperture0_node_values.dot(sm_u.N);
122 ip_data.aperture = ip_data.aperture0;
123
124 auto const initial_effective_stress =
125 _process_data.initial_fracture_effective_stress(0, x_position);
126 for (int i = 0; i < DisplacementDim; i++)
127 {
128 ip_data.sigma_eff[i] = initial_effective_stress[i];
129 ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
130 }
131 }
132}
133
134template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
135 int DisplacementDim>
137 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
138 assembleWithJacobianConcrete(double const t, double const dt,
139 Eigen::VectorXd const& local_x,
140 Eigen::VectorXd const& local_x_prev,
141 Eigen::VectorXd& local_b,
142 Eigen::MatrixXd& local_J)
143{
144 auto const p = local_x.segment(pressure_index, pressure_size);
145 auto const p_prev = local_x_prev.segment(pressure_index, pressure_size);
146 auto const g = local_x.segment(displacement_index, displacement_size);
147 auto const g_prev =
148 local_x_prev.segment(displacement_index, displacement_size);
149
150 auto rhs_p = local_b.segment(pressure_index, pressure_size);
151 auto rhs_g = local_b.segment(displacement_index, displacement_size);
152 auto J_pp = local_J.block(pressure_index, pressure_index, pressure_size,
154 auto J_pg = local_J.block(pressure_index, displacement_index, pressure_size,
156 auto J_gp = local_J.block(displacement_index, pressure_index,
158 auto J_gg = local_J.block(displacement_index, displacement_index,
160
161 assembleBlockMatricesWithJacobian(t, dt, p, p_prev, g, g_prev, rhs_p, rhs_g,
162 J_pp, J_pg, J_gg, J_gp);
163}
164
165template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
166 int DisplacementDim>
168 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
169 assembleBlockMatricesWithJacobian(
170 double const t, double const dt,
171 Eigen::Ref<const Eigen::VectorXd> const& p,
172 Eigen::Ref<const Eigen::VectorXd> const& p_prev,
173 Eigen::Ref<const Eigen::VectorXd> const& g,
174 Eigen::Ref<const Eigen::VectorXd> const& g_prev,
175 Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_g,
176 Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pg,
177 Eigen::Ref<Eigen::MatrixXd> J_gg, Eigen::Ref<Eigen::MatrixXd> J_gp)
178{
179 auto const& R = _fracture_property->R;
180
181 // the index of a normal (normal to a fracture plane) component
182 // in a displacement vector
183 auto constexpr index_normal = DisplacementDim - 1;
184
186 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
188
190 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
192
193 typename ShapeMatricesTypeDisplacement::template MatrixType<
195 Kgp = ShapeMatricesTypeDisplacement::template MatrixType<
198
199 // Projection of the body force vector at the element.
200 Eigen::MatrixXd const global2local_rotation =
201 R.template topLeftCorner<ShapeFunctionPressure::DIM, DisplacementDim>();
202
203 GlobalDimVectorType const gravity_vec = global2local_rotation.transpose() *
204 global2local_rotation *
205 _process_data.specific_body_force;
206
208 x_position.setElementID(_element.getID());
209
210 MPL::VariableArray variables;
211 auto const& medium = _process_data.media_map.getMedium(_element.getID());
212 auto const& liquid_phase = medium->phase("AqueousLiquid");
213 auto const T_ref =
215 .template value<double>(variables, x_position, t, dt);
216 variables.temperature = T_ref;
217
218 unsigned const n_integration_points = _ip_data.size();
219 for (unsigned ip = 0; ip < n_integration_points; ip++)
220 {
221 auto& ip_data = _ip_data[ip];
222 auto const& ip_w = ip_data.integration_weight;
223 auto const& N_p = ip_data.N_p;
224 auto const& dNdx_p = ip_data.dNdx_p;
225 auto const& H_g = ip_data.H_u;
226 auto const& identity2 =
228
229 x_position = {
230 std::nullopt, _element.getID(),
232 NumLib::interpolateCoordinates<ShapeFunctionPressure,
234 _element, N_p))};
235
236 auto& mat = ip_data.fracture_material;
237 auto& effective_stress = ip_data.sigma_eff;
238 auto const& effective_stress_prev = ip_data.sigma_eff_prev;
239 auto& w = ip_data.w;
240 auto const& w_prev = ip_data.w_prev;
241 auto& C = ip_data.C;
242 auto& state = *ip_data.material_state_variables;
243 auto& b_m = ip_data.aperture;
244
245 auto const rho_fr =
246 liquid_phase.property(MPL::PropertyType::density)
247 .template value<double>(variables, x_position, t, dt);
248 variables.density = rho_fr;
249
250 auto const alpha =
252 .template value<double>(variables, x_position, t, dt);
253
254 double const S =
255 medium->property(MPL::PropertyType::storage)
256 .template value<double>(variables, x_position, t, dt);
257
258 auto const mu =
259 liquid_phase.property(MPL::PropertyType::viscosity)
260 .template value<double>(variables, x_position, t, dt);
261
262 // displacement jumps in local coordinates
263 w.noalias() = R * H_g * g;
264
265 // aperture
266 b_m = ip_data.aperture0 + w[index_normal];
267 if (b_m < 0.0)
268 {
269 DBUG(
270 "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it must "
271 "be "
272 "non-negative. Setting it to zero.",
273 _element.getID(), ip, b_m);
274 b_m = 0;
275 }
276
277 auto const initial_effective_stress =
278 _process_data.initial_fracture_effective_stress(0, x_position);
279
280 Eigen::Map<typename HMatricesType::ForceVectorType const> const stress0(
281 initial_effective_stress.data(), initial_effective_stress.size());
282
283 // local C, local stress
284 mat.computeConstitutiveRelation(
285 t, x_position, ip_data.aperture0, stress0, w_prev, w,
286 effective_stress_prev, effective_stress, C, state);
287
288 //
289 // displacement equation, displacement jump part
290 //
291 rhs_g.noalias() -=
292 H_g.transpose() * R.transpose() * effective_stress * ip_w;
293 J_gg.noalias() += H_g.transpose() * R.transpose() * C * R * H_g * ip_w;
294
295 //
296 // displacement equation, pressure part
297 //
298 Kgp.noalias() +=
299 H_g.transpose() * R.transpose() * alpha * identity2 * N_p * ip_w;
300
301 //
302 // pressure equation, pressure part.
303 //
304
305 variables.fracture_aperture = b_m;
306 // Assume that the fracture permeability is isotropic
307 auto const permeability =
308 medium->property(MPL::PropertyType::permeability)
309 .value(variables, x_position, t, dt);
310
311 auto& k = ip_data.permeability;
312 k = std::get<double>(permeability);
313 double const k_over_mu = k / mu;
314 storage_p.noalias() += N_p.transpose() * b_m * S * N_p * ip_w;
315 laplace_p.noalias() +=
316 dNdx_p.transpose() * b_m * k_over_mu * dNdx_p * ip_w;
317 rhs_p.noalias() +=
318 dNdx_p.transpose() * b_m * k_over_mu * rho_fr * gravity_vec * ip_w;
319
320 //
321 // pressure equation, displacement jump part.
322 //
323 GlobalDimVectorType const grad_head = dNdx_p * p + rho_fr * gravity_vec;
324 Eigen::Matrix<double, 1, displacement_size> const mT_R_Hg =
325 identity2.transpose() * R * H_g;
326 // velocity in global coordinates
327 ip_data.darcy_velocity = -k_over_mu * grad_head;
328 J_pg.noalias() +=
329 N_p.transpose() * S * N_p * (p - p_prev) / dt * mT_R_Hg * ip_w;
330
331 // derivative of permeability with respect to aperture
332 double const dk_db_over_mu =
333 medium->property(MPL::PropertyType::permeability)
334 .template dValue<double>(variables,
336 x_position, t, dt) /
337 mu;
338 J_pg.noalias() +=
339 dNdx_p.transpose() * k_over_mu * grad_head * mT_R_Hg * ip_w;
340 J_pg.noalias() += dNdx_p.transpose() * b_m * dk_db_over_mu * grad_head *
341 mT_R_Hg * ip_w;
342 }
343
344 // displacement equation, pressure part
345 J_gp.noalias() -= Kgp;
346
347 // pressure equation, pressure part.
348 J_pp.noalias() += laplace_p + storage_p / dt;
349
350 // pressure equation, displacement jump part.
351 J_pg.noalias() += Kgp.transpose() / dt;
352
353 // pressure equation
354 rhs_p.noalias() -= laplace_p * p + storage_p * (p - p_prev) / dt +
355 Kgp.transpose() * (g - g_prev) / dt;
356
357 // displacement equation
358 rhs_g.noalias() -= -Kgp * p;
359}
360
361template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
362 int DisplacementDim>
364 ShapeFunctionDisplacement, ShapeFunctionPressure,
365 DisplacementDim>::postTimestepConcreteWithVector(const double t,
366 double const /*dt*/,
367 Eigen::VectorXd const&
368 local_x)
369{
370 auto const nodal_g = local_x.segment(displacement_index, displacement_size);
371
372 auto const& R = _fracture_property->R;
373 // the index of a normal (normal to a fracture plane) component
374 // in a displacement vector
375 auto constexpr index_normal = DisplacementDim - 1;
376
378 auto const e_id = _element.getID();
379 x_position.setElementID(e_id);
380
381 unsigned const n_integration_points = _ip_data.size();
382 for (unsigned ip = 0; ip < n_integration_points; ip++)
383 {
384 auto& ip_data = _ip_data[ip];
385 auto const& H_g = ip_data.H_u;
386 auto& mat = ip_data.fracture_material;
387 auto& effective_stress = ip_data.sigma_eff;
388 auto const& effective_stress_prev = ip_data.sigma_eff_prev;
389 auto& w = ip_data.w;
390 auto const& w_prev = ip_data.w_prev;
391 auto& C = ip_data.C;
392 auto& state = *ip_data.material_state_variables;
393 auto& b_m = ip_data.aperture;
394
395 x_position = {
396 std::nullopt, e_id,
398 NumLib::interpolateCoordinates<ShapeFunctionPressure,
400 _element, ip_data.N_p))};
401
402 // displacement jumps in local coordinates
403 w.noalias() = R * H_g * nodal_g;
404
405 // aperture
406 b_m = ip_data.aperture0 + w[index_normal];
407 if (b_m < 0.0)
408 {
409 DBUG(
410 "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it is "
411 "expected to be non-negative. Setting it to zero now.",
412 _element.getID(), ip, b_m);
413 b_m = 0;
414 }
415
416 auto const initial_effective_stress =
417 _process_data.initial_fracture_effective_stress(0, x_position);
418
419 Eigen::Map<typename HMatricesType::ForceVectorType const> const stress0(
420 initial_effective_stress.data(), initial_effective_stress.size());
421
422 // local C, local stress
423 mat.computeConstitutiveRelation(
424 t, x_position, ip_data.aperture0, stress0, w_prev, w,
425 effective_stress_prev, effective_stress, C, state);
426 }
427
428 double ele_b = 0;
429 double ele_k = 0;
430 typename HMatricesType::ForceVectorType ele_sigma_eff =
431 HMatricesType::ForceVectorType::Zero(DisplacementDim);
432 typename HMatricesType::ForceVectorType ele_w =
433 HMatricesType::ForceVectorType::Zero(DisplacementDim);
434 GlobalDimVectorType ele_velocity = GlobalDimVectorType::Zero();
435
436 double ele_Fs = -std::numeric_limits<double>::max();
437 for (auto const& ip : _ip_data)
438 {
439 ele_b += ip.aperture;
440 ele_k += ip.permeability;
441 ele_w += ip.w;
442 ele_sigma_eff += ip.sigma_eff;
443 ele_velocity += ip.darcy_velocity;
444 ele_Fs = std::max(
445 ele_Fs, ip.material_state_variables->getShearYieldFunctionValue());
446 }
447 ele_b /= static_cast<double>(n_integration_points);
448 ele_k /= static_cast<double>(n_integration_points);
449 ele_w /= static_cast<double>(n_integration_points);
450 ele_sigma_eff /= static_cast<double>(n_integration_points);
451 ele_velocity /= static_cast<double>(n_integration_points);
452 auto const element_id = _element.getID();
453 (*_process_data.mesh_prop_b)[element_id] = ele_b;
454 (*_process_data.mesh_prop_k_f)[element_id] = ele_k;
455
456 Eigen::Map<GlobalDimVectorType>(
457 &(*_process_data.element_fracture_stresses)[e_id * DisplacementDim]) =
458 ele_sigma_eff;
459
460 Eigen::Map<GlobalDimVectorType>(
461 &(*_process_data.element_fracture_velocities)[e_id * DisplacementDim]) =
462 ele_velocity;
463
464 Eigen::Map<GlobalDimVectorType>(
465 &(*_process_data.element_local_jumps)[e_id * DisplacementDim]) = ele_w;
466
467 (*_process_data.mesh_prop_fracture_shear_failure)[element_id] = ele_Fs;
468}
469
470template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
471 int DisplacementDim>
472std::vector<double> const& HydroMechanicsLocalAssemblerFracture<
473 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
474 getIntPtFractureVelocity(
475 const double /*t*/,
476 std::vector<GlobalVector*> const& /*x*/,
477 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
478 std::vector<double>& cache) const
479{
480 unsigned const n_integration_points = _ip_data.size();
481 cache.clear();
482 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
483 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
484 cache, DisplacementDim, n_integration_points);
485
486 for (unsigned ip = 0; ip < n_integration_points; ip++)
487 {
488 cache_matrix.col(ip).noalias() = _ip_data[ip].darcy_velocity;
489 }
490
491 return cache;
492}
493
494template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
495 int DisplacementDim>
496std::vector<double> const& HydroMechanicsLocalAssemblerFracture<
497 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
498 getIntPtFractureStress(
499 const double /*t*/,
500 std::vector<GlobalVector*> const& /*x*/,
501 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
502 std::vector<double>& cache) const
503{
504 unsigned const n_integration_points = _ip_data.size();
505 cache.clear();
506 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
507 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
508 cache, DisplacementDim, n_integration_points);
509
510 for (unsigned ip = 0; ip < n_integration_points; ip++)
511 {
512 cache_matrix.col(ip).noalias() = _ip_data[ip].sigma_eff;
513 }
514
515 return cache;
516}
517
518template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
519 int DisplacementDim>
520std::vector<double> const& HydroMechanicsLocalAssemblerFracture<
521 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
522 getIntPtFractureAperture(
523 const double /*t*/,
524 std::vector<GlobalVector*> const& /*x*/,
525 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
526 std::vector<double>& cache) const
527{
530}
531
532template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
533 int DisplacementDim>
534std::vector<double> const& HydroMechanicsLocalAssemblerFracture<
535 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
536 getIntPtFracturePermeability(
537 const double /*t*/,
538 std::vector<GlobalVector*> const& /*x*/,
539 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
540 std::vector<double>& cache) const
541{
544}
545
546} // namespace HydroMechanics
547} // namespace LIE
548} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
constexpr double getWeight() const
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
void postTimestepConcreteWithVector(const double t, double const dt, Eigen::VectorXd const &local_x) override
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assembleBlockMatricesWithJacobian(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_prev, Eigen::Ref< const Eigen::VectorXd > const &g, Eigen::Ref< const Eigen::VectorXd > const &g_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_g, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pg, Eigen::Ref< Eigen::MatrixXd > J_gg, Eigen::Ref< Eigen::MatrixXd > J_gp)
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
HydroMechanicsLocalAssemblerInterface(MeshLib::Element const &element, bool const is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
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)
Eigen::Matrix< double, DisplacementDim, DisplacementDim > createRotatedTensor(RotationMatrix const &R, double const value)
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
void computeHMatrix(N_Type const &N, HMatrixType &H)
Fills a H-matrix based on given shape function.
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
static MATERIALLIB_FRACTUREMODELS_EXPORT VectorType const value