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