OGS
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
class ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >

Definition at line 30 of file HydroMechanicsLocalAssemblerFracture.h.

#include <HydroMechanicsLocalAssemblerFracture.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]

Public Member Functions

 HydroMechanicsLocalAssemblerFracture (HydroMechanicsLocalAssemblerFracture const &)=delete
 
 HydroMechanicsLocalAssemblerFracture (HydroMechanicsLocalAssemblerFracture &&)=delete
 
 HydroMechanicsLocalAssemblerFracture (MeshLib::Element const &e, std::size_t const local_matrix_size, std::vector< unsigned > const &dofIndex_to_localIndex, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HydroMechanicsProcessData< DisplacementDim > &process_data)
 
void preTimestepConcrete (std::vector< double > const &, double const, double const) override
 
void postTimestepConcreteWithVector (const double t, double const dt, Eigen::VectorXd const &local_x) override
 
std::vector< double > const & getIntPtSigma (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilon (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtDarcyVelocity (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtFractureVelocity (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 & getIntPtFractureStress (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 & getIntPtFractureAperture (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 & getIntPtFracturePermeability (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
- Public Member Functions inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface
 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)
 
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
 
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &local_x_prev_, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
 
void postTimestepConcrete (Eigen::VectorXd const &local_x_, Eigen::VectorXd const &, const double t, double const dt, int const) override
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
 
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
 
virtual void preAssemble (double const, double const, std::vector< double > const &)
 
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
 
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
 
virtual std::optional< VectorSegmentgetVectorDeformationSegment () const
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Types

using ShapeMatricesTypeDisplacement
 
using HMatricesType
 
using HMatrixType = typename HMatricesType::HMatrixType
 
using ShapeMatricesTypePressure
 
using IntegrationPointDataType
 
using GlobalDimVectorType = Eigen::Matrix<double, DisplacementDim, 1>
 

Private Member Functions

void assembleWithJacobianConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
 
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)
 

Private Attributes

HydroMechanicsProcessData< DisplacementDim > & _process_data
 
FractureProperty const * _fracture_property = nullptr
 
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
 
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
 

Static Private Attributes

static const int pressure_index = 0
 
static const int pressure_size = ShapeFunctionPressure::NPOINTS
 
static const int displacement_index = ShapeFunctionPressure::NPOINTS
 
static const int displacement_size
 

Additional Inherited Members

- Protected Member Functions inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface
- Protected Attributes inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface
MeshLib::Element const & _element
 
bool const _is_axially_symmetric
 
NumLib::GenericIntegrationMethod const & _integration_method
 

Member Typedef Documentation

◆ GlobalDimVectorType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::GlobalDimVectorType = Eigen::Matrix<double, DisplacementDim, 1>
private

Definition at line 157 of file HydroMechanicsLocalAssemblerFracture.h.

◆ HMatricesType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HMatricesType
private
Initial value:
HMatrixPolicyType<ShapeFunctionDisplacement, DisplacementDim>

Definition at line 144 of file HydroMechanicsLocalAssemblerFracture.h.

◆ HMatrixType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HMatrixType = typename HMatricesType::HMatrixType
private

Definition at line 146 of file HydroMechanicsLocalAssemblerFracture.h.

◆ IntegrationPointDataType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::IntegrationPointDataType
private
Initial value:
IntegrationPointDataFracture<
DisplacementDim>
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
HMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > HMatricesType
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure

Definition at line 153 of file HydroMechanicsLocalAssemblerFracture.h.

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypeDisplacement
private

◆ ShapeMatricesTypePressure

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypePressure
private

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerFracture() [1/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HydroMechanicsLocalAssemblerFracture ( HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > const & )
delete

◆ HydroMechanicsLocalAssemblerFracture() [2/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HydroMechanicsLocalAssemblerFracture ( HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > && )
delete

◆ HydroMechanicsLocalAssemblerFracture() [3/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HydroMechanicsLocalAssemblerFracture ( MeshLib::Element const & e,
std::size_t const local_matrix_size,
std::vector< unsigned > const & dofIndex_to_localIndex,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
HydroMechanicsProcessData< DisplacementDim > & process_data )

Definition at line 44 of file HydroMechanicsLocalAssemblerFracture-impl.h.

53 e, is_axially_symmetric, integration_method,
54 ShapeFunctionDisplacement::NPOINTS * DisplacementDim +
55 ShapeFunctionPressure::NPOINTS,
56 dofIndex_to_localIndex),
57 _process_data(process_data)
58{
59 assert(e.getDimension() == DisplacementDim - 1);
60
61 unsigned const n_integration_points =
62 integration_method.getNumberOfPoints();
63
64 _ip_data.reserve(n_integration_points);
65 _secondary_data.N.resize(n_integration_points);
66
67 auto const shape_matrices_u =
68 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
70 DisplacementDim>(e, is_axially_symmetric,
71 integration_method);
72
73 auto const shape_matrices_p =
74 NumLib::initShapeMatrices<ShapeFunctionPressure,
75 ShapeMatricesTypePressure, DisplacementDim>(
76 e, is_axially_symmetric, integration_method);
77
78 auto mat_id = (*_process_data.mesh_prop_materialIDs)[e.getID()];
79 auto frac_id = _process_data.map_materialID_to_fractureID[mat_id];
80 _fracture_property = &_process_data.fracture_properties[frac_id];
81
82 // Get element nodes for aperture0 interpolation from nodes to integration
83 // point. The aperture0 parameter is time-independent.
85 aperture0_node_values =
87 e, /*time independent*/ 0);
88
89 for (unsigned ip = 0; ip < n_integration_points; ip++)
90 {
91 _ip_data.emplace_back(*_process_data.fracture_model);
92 auto const& sm_u = shape_matrices_u[ip];
93 auto const& sm_p = shape_matrices_p[ip];
95 std::nullopt, _element.getID(),
97 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
99 _element, sm_u.N))};
100
101 auto& ip_data = _ip_data[ip];
102 ip_data.integration_weight =
103 sm_u.detJ * sm_u.integralMeasure *
104 integration_method.getWeightedPoint(ip).getWeight();
105
106 ip_data.H_u.setZero(
107 DisplacementDim,
108 ShapeFunctionDisplacement::NPOINTS * DisplacementDim);
110 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
112 HMatrixType>(sm_u.N, ip_data.H_u);
113 ip_data.N_p = sm_p.N;
114 ip_data.dNdx_p = sm_p.dNdx;
115
116 _secondary_data.N[ip] = sm_u.N;
117
118 // Initialize current time step values
119 ip_data.w.setZero(DisplacementDim);
120 ip_data.sigma_eff.setZero(DisplacementDim);
121
122 // Previous time step values are not initialized and are set later.
123 ip_data.w_prev.resize(DisplacementDim);
124 ip_data.sigma_eff_prev.resize(DisplacementDim);
125
126 ip_data.C.resize(DisplacementDim, DisplacementDim);
127
128 ip_data.aperture0 = aperture0_node_values.dot(sm_u.N);
129 ip_data.aperture = ip_data.aperture0;
130
131 auto const initial_effective_stress =
132 _process_data.initial_fracture_effective_stress(0, x_position);
133 for (int i = 0; i < DisplacementDim; i++)
134 {
135 ip_data.sigma_eff[i] = initial_effective_stress[i];
136 ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
137 }
138 }
139}
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:91
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
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)
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)
void computeHMatrix(N_Type const &N, HMatrixType &H)
Fills a H-matrix based on given shape function.
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
virtual Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > getNodalValuesOnElement(MeshLib::Element const &element, double const t) const
Returns a matrix of values for all nodes of the given element.
Definition Parameter.h:164
ParameterLib::Parameter< double > const & aperture0
Initial aperture.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_fracture_property, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_process_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_secondary_data, ProcessLib::LIE::FractureProperty::aperture0, ProcessLib::computeHMatrix(), MeshLib::Element::getDimension(), MeshLib::Element::getID(), ParameterLib::Parameter< T >::getNodalValuesOnElement(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), NumLib::interpolateCoordinates(), and ProcessLib::LIE::HydroMechanics::SecondaryData< ShapeMatrixType >::N.

Member Function Documentation

◆ assembleBlockMatricesWithJacobian()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::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 )
private

Definition at line 175 of file HydroMechanicsLocalAssemblerFracture-impl.h.

185{
186 auto const& R = _fracture_property->R;
187
188 // the index of a normal (normal to a fracture plane) component
189 // in a displacement vector
190 auto constexpr index_normal = DisplacementDim - 1;
191
193 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
195
197 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
199
200 typename ShapeMatricesTypeDisplacement::template MatrixType<
202 Kgp = ShapeMatricesTypeDisplacement::template MatrixType<
205
206 // Projection of the body force vector at the element.
207 Eigen::MatrixXd const global2local_rotation =
208 R.template topLeftCorner<ShapeFunctionPressure::DIM, DisplacementDim>();
209
210 GlobalDimVectorType const gravity_vec = global2local_rotation.transpose() *
211 global2local_rotation *
212 _process_data.specific_body_force;
213
215 x_position.setElementID(_element.getID());
216
217 MPL::VariableArray variables;
218 auto const& medium = _process_data.media_map.getMedium(_element.getID());
219 auto const& liquid_phase = medium->phase("AqueousLiquid");
220 auto const T_ref =
221 medium->property(MPL::PropertyType::reference_temperature)
222 .template value<double>(variables, x_position, t, dt);
223 variables.temperature = T_ref;
224
225 unsigned const n_integration_points = _ip_data.size();
226 for (unsigned ip = 0; ip < n_integration_points; ip++)
227 {
228 auto& ip_data = _ip_data[ip];
229 auto const& ip_w = ip_data.integration_weight;
230 auto const& N_p = ip_data.N_p;
231 auto const& dNdx_p = ip_data.dNdx_p;
232 auto const& H_g = ip_data.H_u;
233 auto const& identity2 =
235
236 x_position = {
237 std::nullopt, _element.getID(),
239 NumLib::interpolateCoordinates<ShapeFunctionPressure,
241 _element, N_p))};
242
243 auto& mat = ip_data.fracture_material;
244 auto& effective_stress = ip_data.sigma_eff;
245 auto const& effective_stress_prev = ip_data.sigma_eff_prev;
246 auto& w = ip_data.w;
247 auto const& w_prev = ip_data.w_prev;
248 auto& C = ip_data.C;
249 auto& state = *ip_data.material_state_variables;
250 auto& b_m = ip_data.aperture;
251
252 auto const rho_fr =
253 liquid_phase.property(MPL::PropertyType::density)
254 .template value<double>(variables, x_position, t, dt);
255 variables.density = rho_fr;
256
257 auto const alpha =
258 medium->property(MPL::PropertyType::biot_coefficient)
259 .template value<double>(variables, x_position, t, dt);
260
261 double const S =
262 medium->property(MPL::PropertyType::storage)
263 .template value<double>(variables, x_position, t, dt);
264
265 auto const mu =
266 liquid_phase.property(MPL::PropertyType::viscosity)
267 .template value<double>(variables, x_position, t, dt);
268
269 // displacement jumps in local coordinates
270 w.noalias() = R * H_g * g;
271
272 // aperture
273 b_m = ip_data.aperture0 + w[index_normal];
274 if (b_m < 0.0)
275 {
276 DBUG(
277 "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it must "
278 "be "
279 "non-negative. Setting it to zero.",
280 _element.getID(), ip, b_m);
281 b_m = 0;
282 }
283
284 auto const initial_effective_stress =
285 _process_data.initial_fracture_effective_stress(0, x_position);
286
287 Eigen::Map<typename HMatricesType::ForceVectorType const> const stress0(
288 initial_effective_stress.data(), initial_effective_stress.size());
289
290 // local C, local stress
291 mat.computeConstitutiveRelation(
292 t, x_position, ip_data.aperture0, stress0, w_prev, w,
293 effective_stress_prev, effective_stress, C, state);
294
295 //
296 // displacement equation, displacement jump part
297 //
298 rhs_g.noalias() -=
299 H_g.transpose() * R.transpose() * effective_stress * ip_w;
300 J_gg.noalias() += H_g.transpose() * R.transpose() * C * R * H_g * ip_w;
301
302 //
303 // displacement equation, pressure part
304 //
305 Kgp.noalias() +=
306 H_g.transpose() * R.transpose() * alpha * identity2 * N_p * ip_w;
307
308 //
309 // pressure equation, pressure part.
310 //
311
312 variables.fracture_aperture = b_m;
313 // Assume that the fracture permeability is isotropic
314 auto const permeability =
315 medium->property(MPL::PropertyType::permeability)
316 .value(variables, x_position, t, dt);
317
318 auto& k = ip_data.permeability;
319 k = std::get<double>(permeability);
320 double const k_over_mu = k / mu;
321 storage_p.noalias() += N_p.transpose() * b_m * S * N_p * ip_w;
322 laplace_p.noalias() +=
323 dNdx_p.transpose() * b_m * k_over_mu * dNdx_p * ip_w;
324 rhs_p.noalias() +=
325 dNdx_p.transpose() * b_m * k_over_mu * rho_fr * gravity_vec * ip_w;
326
327 //
328 // pressure equation, displacement jump part.
329 //
330 GlobalDimVectorType const grad_head = dNdx_p * p + rho_fr * gravity_vec;
331 Eigen::Matrix<double, 1, displacement_size> const mT_R_Hg =
332 identity2.transpose() * R * H_g;
333 // velocity in global coordinates
334 ip_data.darcy_velocity = -k_over_mu * grad_head;
335 J_pg.noalias() +=
336 N_p.transpose() * S * N_p * (p - p_prev) / dt * mT_R_Hg * ip_w;
337
338 // derivative of permeability with respect to aperture
339 double const dk_db_over_mu =
340 medium->property(MPL::PropertyType::permeability)
341 .template dValue<double>(variables,
342 MPL::Variable::fracture_aperture,
343 x_position, t, dt) /
344 mu;
345 J_pg.noalias() +=
346 dNdx_p.transpose() * k_over_mu * grad_head * mT_R_Hg * ip_w;
347 J_pg.noalias() += dNdx_p.transpose() * b_m * dk_db_over_mu * grad_head *
348 mT_R_Hg * ip_w;
349 }
350
351 // displacement equation, pressure part
352 J_gp.noalias() -= Kgp;
353
354 // pressure equation, pressure part.
355 J_pp.noalias() += laplace_p + storage_p / dt;
356
357 // pressure equation, displacement jump part.
358 J_pg.noalias() += Kgp.transpose() / dt;
359
360 // pressure equation
361 rhs_p.noalias() -= laplace_p * p + storage_p * (p - p_prev) / dt +
362 Kgp.transpose() * (g - g_prev) / dt;
363
364 // displacement equation
365 rhs_g.noalias() -= -Kgp * p;
366}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void setElementID(std::size_t element_id)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
Eigen::MatrixXd R
Rotation matrix from global to local coordinates.

References DBUG(), MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::VariableArray::fracture_aperture, NumLib::interpolateCoordinates(), ParameterLib::SpatialPosition::setElementID(), and MaterialPropertyLib::VariableArray::temperature.

◆ assembleWithJacobianConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianConcrete ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
Eigen::VectorXd & local_b,
Eigen::MatrixXd & local_J )
overrideprivatevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 144 of file HydroMechanicsLocalAssemblerFracture-impl.h.

150{
151 auto const p = local_x.segment(pressure_index, pressure_size);
152 auto const p_prev = local_x_prev.segment(pressure_index, pressure_size);
153 auto const g = local_x.segment(displacement_index, displacement_size);
154 auto const g_prev =
155 local_x_prev.segment(displacement_index, displacement_size);
156
157 auto rhs_p = local_b.segment(pressure_index, pressure_size);
158 auto rhs_g = local_b.segment(displacement_index, displacement_size);
159 auto J_pp = local_J.block(pressure_index, pressure_index, pressure_size,
161 auto J_pg = local_J.block(pressure_index, displacement_index, pressure_size,
163 auto J_gp = local_J.block(displacement_index, pressure_index,
165 auto J_gg = local_J.block(displacement_index, displacement_index,
167
168 assembleBlockMatricesWithJacobian(t, dt, p, p_prev, g, g_prev, rhs_p, rhs_g,
169 J_pp, J_pg, J_gg, J_gp);
170}
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)

◆ getIntPtDarcyVelocity()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtDarcyVelocity ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 81 of file HydroMechanicsLocalAssemblerFracture.h.

86 {
87 cache.resize(0);
88 return cache;
89 }

◆ getIntPtEpsilon()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilon ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 71 of file HydroMechanicsLocalAssemblerFracture.h.

76 {
77 cache.resize(0);
78 return cache;
79 }

◆ getIntPtFractureAperture()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFractureAperture ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 528 of file HydroMechanicsLocalAssemblerFracture-impl.h.

534{
537}
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)

References ProcessLib::getIntegrationPointScalarData().

◆ getIntPtFracturePermeability()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFracturePermeability ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

◆ getIntPtFractureStress()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFractureStress ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 504 of file HydroMechanicsLocalAssemblerFracture-impl.h.

510{
511 unsigned const n_integration_points = _ip_data.size();
512 cache.clear();
513 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
514 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
515 cache, DisplacementDim, n_integration_points);
516
517 for (unsigned ip = 0; ip < n_integration_points; ip++)
518 {
519 cache_matrix.col(ip).noalias() = _ip_data[ip].sigma_eff;
520 }
521
522 return cache;
523}
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)

References MathLib::createZeroedMatrix().

◆ getIntPtFractureVelocity()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFractureVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 480 of file HydroMechanicsLocalAssemblerFracture-impl.h.

486{
487 unsigned const n_integration_points = _ip_data.size();
488 cache.clear();
489 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
490 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
491 cache, DisplacementDim, n_integration_points);
492
493 for (unsigned ip = 0; ip < n_integration_points; ip++)
494 {
495 cache_matrix.col(ip).noalias() = _ip_data[ip].darcy_velocity;
496 }
497
498 return cache;
499}

References MathLib::createZeroedMatrix().

◆ getIntPtSigma()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtSigma ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 61 of file HydroMechanicsLocalAssemblerFracture.h.

66 {
67 cache.resize(0);
68 return cache;
69 }

◆ getShapeMatrix()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 115 of file HydroMechanicsLocalAssemblerFracture.h.

117 {
118 auto const& N = _secondary_data.N[integration_point];
119
120 // assumes N is stored contiguously in memory
121 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
122 }

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_secondary_data, and ProcessLib::LIE::HydroMechanics::SecondaryData< ShapeMatrixType >::N.

◆ postTimestepConcreteWithVector()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcreteWithVector ( const double t,
double const dt,
Eigen::VectorXd const & local_x )
overridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 372 of file HydroMechanicsLocalAssemblerFracture-impl.h.

376{
377 auto const nodal_g = local_x.segment(displacement_index, displacement_size);
378
379 auto const& R = _fracture_property->R;
380 // the index of a normal (normal to a fracture plane) component
381 // in a displacement vector
382 auto constexpr index_normal = DisplacementDim - 1;
383
385 auto const e_id = _element.getID();
386 x_position.setElementID(e_id);
387
388 unsigned const n_integration_points = _ip_data.size();
389 for (unsigned ip = 0; ip < n_integration_points; ip++)
390 {
391 auto& ip_data = _ip_data[ip];
392 auto const& H_g = ip_data.H_u;
393 auto& mat = ip_data.fracture_material;
394 auto& effective_stress = ip_data.sigma_eff;
395 auto const& effective_stress_prev = ip_data.sigma_eff_prev;
396 auto& w = ip_data.w;
397 auto const& w_prev = ip_data.w_prev;
398 auto& C = ip_data.C;
399 auto& state = *ip_data.material_state_variables;
400 auto& b_m = ip_data.aperture;
401
402 x_position = {
403 std::nullopt, e_id,
405 NumLib::interpolateCoordinates<ShapeFunctionPressure,
407 _element, ip_data.N_p))};
408
409 // displacement jumps in local coordinates
410 w.noalias() = R * H_g * nodal_g;
411
412 // aperture
413 b_m = ip_data.aperture0 + w[index_normal];
414 if (b_m < 0.0)
415 {
416 DBUG(
417 "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it is "
418 "expected to be non-negative. Setting it to zero now.",
419 _element.getID(), ip, b_m);
420 b_m = 0;
421 }
422
423 auto const initial_effective_stress =
424 _process_data.initial_fracture_effective_stress(0, x_position);
425
426 Eigen::Map<typename HMatricesType::ForceVectorType const> const stress0(
427 initial_effective_stress.data(), initial_effective_stress.size());
428
429 // local C, local stress
430 mat.computeConstitutiveRelation(
431 t, x_position, ip_data.aperture0, stress0, w_prev, w,
432 effective_stress_prev, effective_stress, C, state);
433 }
434
435 double ele_b = 0;
436 double ele_k = 0;
437 typename HMatricesType::ForceVectorType ele_sigma_eff =
438 HMatricesType::ForceVectorType::Zero(DisplacementDim);
439 typename HMatricesType::ForceVectorType ele_w =
440 HMatricesType::ForceVectorType::Zero(DisplacementDim);
441 GlobalDimVectorType ele_velocity = GlobalDimVectorType::Zero();
442
443 double ele_Fs = -std::numeric_limits<double>::max();
444 for (auto const& ip : _ip_data)
445 {
446 ele_b += ip.aperture;
447 ele_k += ip.permeability;
448 ele_w += ip.w;
449 ele_sigma_eff += ip.sigma_eff;
450 ele_velocity += ip.darcy_velocity;
451 ele_Fs = std::max(
452 ele_Fs, ip.material_state_variables->getShearYieldFunctionValue());
453 }
454 ele_b /= static_cast<double>(n_integration_points);
455 ele_k /= static_cast<double>(n_integration_points);
456 ele_w /= static_cast<double>(n_integration_points);
457 ele_sigma_eff /= static_cast<double>(n_integration_points);
458 ele_velocity /= static_cast<double>(n_integration_points);
459 auto const element_id = _element.getID();
460 (*_process_data.mesh_prop_b)[element_id] = ele_b;
461 (*_process_data.mesh_prop_k_f)[element_id] = ele_k;
462
463 Eigen::Map<GlobalDimVectorType>(
464 &(*_process_data.element_fracture_stresses)[e_id * DisplacementDim]) =
465 ele_sigma_eff;
466
467 Eigen::Map<GlobalDimVectorType>(
468 &(*_process_data.element_fracture_velocities)[e_id * DisplacementDim]) =
469 ele_velocity;
470
471 Eigen::Map<GlobalDimVectorType>(
472 &(*_process_data.element_local_jumps)[e_id * DisplacementDim]) = ele_w;
473
474 (*_process_data.mesh_prop_fracture_shear_failure)[element_id] = ele_Fs;
475}
VectorType< DisplacementDim > ForceVectorType

References DBUG(), NumLib::interpolateCoordinates(), and ParameterLib::SpatialPosition::setElementID().

◆ preTimestepConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::preTimestepConcrete ( std::vector< double > const & ,
double const ,
double const  )
inlineoverridevirtual

Member Data Documentation

◆ _fracture_property

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
FractureProperty const* ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_fracture_property = nullptr
private

◆ _ip_data

◆ _process_data

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
HydroMechanicsProcessData<DisplacementDim>& ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_process_data
private

◆ _secondary_data

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType> ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_secondary_data
private

◆ displacement_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_index = ShapeFunctionPressure::NPOINTS
staticprivate

Definition at line 169 of file HydroMechanicsLocalAssemblerFracture.h.

◆ displacement_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_size
staticprivate
Initial value:
=
ShapeFunctionDisplacement::NPOINTS * DisplacementDim

Definition at line 170 of file HydroMechanicsLocalAssemblerFracture.h.

◆ pressure_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_index = 0
staticprivate

Definition at line 167 of file HydroMechanicsLocalAssemblerFracture.h.

◆ pressure_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_size = ShapeFunctionPressure::NPOINTS
staticprivate

Definition at line 168 of file HydroMechanicsLocalAssemblerFracture.h.


The documentation for this class was generated from the following files: