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