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 23 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 int getNumberOfVectorElementsForDeformation () 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 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 150 of file HydroMechanicsLocalAssemblerFracture.h.

◆ HMatricesType

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

◆ HMatrixType

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

Definition at line 139 of file HydroMechanicsLocalAssemblerFracture.h.

◆ IntegrationPointDataType

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

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypeDisplacement
private
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 135 of file HydroMechanicsLocalAssemblerFracture.h.

◆ 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 37 of file HydroMechanicsLocalAssemblerFracture-impl.h.

51{
52 assert(e.getDimension() == DisplacementDim - 1);
53
54 unsigned const n_integration_points =
55 integration_method.getNumberOfPoints();
56
59
60 auto const shape_matrices_u =
65
66 auto const shape_matrices_p =
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.
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(),
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(
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
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}
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)
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::HydroMechanicsLocalAssemblerInterface(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element, _fracture_property, _ip_data, _process_data, _secondary_data, ProcessLib::computeHMatrix(), MeshLib::Element::getDimension(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), and NumLib::interpolateCoordinates().

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 168 of file HydroMechanicsLocalAssemblerFracture-impl.h.

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
188
192
198
199 // Projection of the body force vector at the element.
202
205 _process_data.specific_body_force;
206
208 x_position.setElementID(_element.getID());
209
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(),
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 =
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 =
256 .template value<double>(variables, x_position, t, dt);
257
258 auto const mu =
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
282
283 // local C, local stress
284 mat.computeConstitutiveRelation(
285 t, x_position, ip_data.aperture0, stress0, w_prev, w,
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 =
309 .value(variables, x_position, t, dt);
310
311 auto& k = ip_data.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 //
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 =
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}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element, _fracture_property, _ip_data, _process_data, MaterialPropertyLib::biot_coefficient, DBUG(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, displacement_size, MaterialPropertyLib::fracture_aperture, MaterialPropertyLib::VariableArray::fracture_aperture, NumLib::interpolateCoordinates(), MaterialPropertyLib::permeability, pressure_size, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::storage, MaterialPropertyLib::VariableArray::temperature, MaterialLib::Fracture::FractureIdentity2< DisplacementDim >::value, and MaterialPropertyLib::viscosity.

Referenced by assembleWithJacobianConcrete().

◆ 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 137 of file HydroMechanicsLocalAssemblerFracture-impl.h.

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 =
149
150 auto rhs_p = local_b.segment(pressure_index, pressure_size);
160
162 J_pp, J_pg, J_gg, J_gp);
163}
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)

References assembleBlockMatricesWithJacobian(), displacement_index, displacement_size, pressure_index, and pressure_size.

◆ 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 74 of file HydroMechanicsLocalAssemblerFracture.h.

79 {
80 cache.resize(0);
81 return cache;
82 }

◆ 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 64 of file HydroMechanicsLocalAssemblerFracture.h.

69 {
70 cache.resize(0);
71 return cache;
72 }

◆ 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

◆ getIntPtFracturePermeability()

◆ 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 497 of file HydroMechanicsLocalAssemblerFracture-impl.h.

503{
504 unsigned const n_integration_points = _ip_data.size();
505 cache.clear();
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}

References _ip_data, and 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 473 of file HydroMechanicsLocalAssemblerFracture-impl.h.

479{
480 unsigned const n_integration_points = _ip_data.size();
481 cache.clear();
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}

References _ip_data, and 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 54 of file HydroMechanicsLocalAssemblerFracture.h.

59 {
60 cache.resize(0);
61 return cache;
62 }

◆ 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 108 of file HydroMechanicsLocalAssemblerFracture.h.

110 {
111 auto const& N = _secondary_data.N[integration_point];
112
113 // assumes N is stored contiguously in memory
114 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
115 }

References _secondary_data.

◆ 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 365 of file HydroMechanicsLocalAssemblerFracture-impl.h.

369{
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 = {
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
421
422 // local C, local stress
423 mat.computeConstitutiveRelation(
424 t, x_position, ip_data.aperture0, stress0, w_prev, w,
426 }
427
428 double ele_b = 0;
429 double ele_k = 0;
435
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;
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
457 &(*_process_data.element_fracture_stresses)[e_id * DisplacementDim]) =
459
461 &(*_process_data.element_fracture_velocities)[e_id * DisplacementDim]) =
463
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}

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element, _fracture_property, _ip_data, _process_data, DBUG(), displacement_index, displacement_size, NumLib::interpolateCoordinates(), postTimestepConcreteWithVector(), and ParameterLib::SpatialPosition::setElementID().

Referenced by postTimestepConcreteWithVector().

◆ 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

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 40 of file HydroMechanicsLocalAssemblerFracture.h.

43 {
44 for (auto& data : _ip_data)
45 {
46 data.pushBackState();
47 }
48 }

References _ip_data.

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

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector<IntegrationPointDataType, Eigen::aligned_allocator<IntegrationPointDataType> > ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data
private

◆ _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

◆ 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 163 of file HydroMechanicsLocalAssemblerFracture.h.

Referenced by assembleBlockMatricesWithJacobian(), assembleWithJacobianConcrete(), and postTimestepConcreteWithVector().

◆ pressure_index

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

◆ pressure_size

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

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