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

Detailed Description

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

Definition at line 30 of file HydroMechanicsLocalAssemblerFracture.h.

#include <HydroMechanicsLocalAssemblerFracture.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >:
[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< GlobalDim > &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.
 
- 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, GlobalDim, 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< GlobalDim > & _process_data
 
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 GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::GlobalDimVectorType = Eigen::Matrix<double, GlobalDim, 1>
private

Definition at line 158 of file HydroMechanicsLocalAssemblerFracture.h.

◆ HMatricesType

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

Definition at line 144 of file HydroMechanicsLocalAssemblerFracture.h.

◆ HMatrixType

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

Definition at line 146 of file HydroMechanicsLocalAssemblerFracture.h.

◆ IntegrationPointDataType

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

Definition at line 153 of file HydroMechanicsLocalAssemblerFracture.h.

◆ ShapeMatricesTypeDisplacement

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

◆ ShapeMatricesTypePressure

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

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerFracture() [1/3]

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

◆ HydroMechanicsLocalAssemblerFracture() [2/3]

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

◆ HydroMechanicsLocalAssemblerFracture() [3/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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< GlobalDim > & process_data )

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

53 e, is_axially_symmetric, integration_method,
54 ShapeFunctionDisplacement::NPOINTS * GlobalDim +
55 ShapeFunctionPressure::NPOINTS,
56 dofIndex_to_localIndex),
57 _process_data(process_data)
58{
59 assert(e.getDimension() == GlobalDim - 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 e, is_axially_symmetric, integration_method);
71
72 auto const shape_matrices_p =
73 NumLib::initShapeMatrices<ShapeFunctionPressure,
74 ShapeMatricesTypePressure, GlobalDim>(
75 e, is_axially_symmetric, integration_method);
76
77 auto const& frac_prop = *_process_data.fracture_property;
78
79 // Get element nodes for aperture0 interpolation from nodes to integration
80 // point. The aperture0 parameter is time-independent.
82 aperture0_node_values = frac_prop.aperture0.getNodalValuesOnElement(
83 e, /*time independent*/ 0);
84
86 x_position.setElementID(e.getID());
87 for (unsigned ip = 0; ip < n_integration_points; ip++)
88 {
89 x_position.setIntegrationPoint(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];
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(GlobalDim,
100 ShapeFunctionDisplacement::NPOINTS * GlobalDim);
102 GlobalDim, ShapeFunctionDisplacement::NPOINTS,
104 HMatrixType>(sm_u.N, ip_data.H_u);
105 ip_data.N_p = sm_p.N;
106 ip_data.dNdx_p = sm_p.dNdx;
107
108 _secondary_data.N[ip] = sm_u.N;
109
110 // Initialize current time step values
111 ip_data.w.setZero(GlobalDim);
112 ip_data.sigma_eff.setZero(GlobalDim);
113
114 // Previous time step values are not initialized and are set later.
115 ip_data.w_prev.resize(GlobalDim);
116 ip_data.sigma_eff_prev.resize(GlobalDim);
117
118 ip_data.C.resize(GlobalDim, GlobalDim);
119
120 ip_data.aperture0 = aperture0_node_values.dot(sm_u.N);
121 ip_data.aperture = ip_data.aperture0;
122
123 auto const initial_effective_stress =
124 _process_data.initial_fracture_effective_stress(0, x_position);
125 for (int i = 0; i < GlobalDim; i++)
126 {
127 ip_data.sigma_eff[i] = initial_effective_stress[i];
128 ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
129 }
130 }
131}
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
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)
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
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_ip_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_process_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_secondary_data, ProcessLib::computeHMatrix(), MeshLib::Element::getDimension(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::LIE::HydroMechanics::SecondaryData< ShapeMatrixType >::N, ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

Member Function Documentation

◆ assembleBlockMatricesWithJacobian()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 167 of file HydroMechanicsLocalAssemblerFracture-impl.h.

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

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

◆ assembleWithJacobianConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 136 of file HydroMechanicsLocalAssemblerFracture-impl.h.

142{
143 auto const p = local_x.segment(pressure_index, pressure_size);
144 auto const p_prev = local_x_prev.segment(pressure_index, pressure_size);
145 auto const g = local_x.segment(displacement_index, displacement_size);
146 auto const g_prev =
147 local_x_prev.segment(displacement_index, displacement_size);
148
149 auto rhs_p = local_b.segment(pressure_index, pressure_size);
150 auto rhs_g = local_b.segment(displacement_index, displacement_size);
151 auto J_pp = local_J.block(pressure_index, pressure_index, pressure_size,
153 auto J_pg = local_J.block(pressure_index, displacement_index, pressure_size,
155 auto J_gp = local_J.block(displacement_index, pressure_index,
157 auto J_gg = local_J.block(displacement_index, displacement_index,
159
160 assembleBlockMatricesWithJacobian(t, dt, p, p_prev, g, g_prev, rhs_p, rhs_g,
161 J_pp, J_pg, J_gg, J_gp);
162}
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 GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 511 of file HydroMechanicsLocalAssemblerFracture-impl.h.

517{
520}
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 GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 487 of file HydroMechanicsLocalAssemblerFracture-impl.h.

493{
494 unsigned const n_integration_points = _ip_data.size();
495 cache.clear();
496 auto cache_matrix = MathLib::createZeroedMatrix<
497 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
498 cache, GlobalDim, n_integration_points);
499
500 for (unsigned ip = 0; ip < n_integration_points; ip++)
501 {
502 cache_matrix.col(ip).noalias() = _ip_data[ip].sigma_eff;
503 }
504
505 return cache;
506}
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 GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 463 of file HydroMechanicsLocalAssemblerFracture-impl.h.

469{
470 unsigned const n_integration_points = _ip_data.size();
471 cache.clear();
472 auto cache_matrix = MathLib::createZeroedMatrix<
473 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
474 cache, GlobalDim, n_integration_points);
475
476 for (unsigned ip = 0; ip < n_integration_points; ip++)
477 {
478 cache_matrix.col(ip).noalias() = _ip_data[ip].darcy_velocity;
479 }
480
481 return cache;
482}

References MathLib::createZeroedMatrix().

◆ getIntPtSigma()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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, GlobalDim >::_secondary_data, and ProcessLib::LIE::HydroMechanics::SecondaryData< ShapeMatrixType >::N.

◆ postTimestepConcreteWithVector()

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

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

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

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

References DBUG(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ preTimestepConcrete()

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

Member Data Documentation

◆ _ip_data

◆ _process_data

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

◆ _secondary_data

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

◆ displacement_index

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

Definition at line 168 of file HydroMechanicsLocalAssemblerFracture.h.

◆ displacement_size

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

Definition at line 169 of file HydroMechanicsLocalAssemblerFracture.h.

◆ pressure_index

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

Definition at line 166 of file HydroMechanicsLocalAssemblerFracture.h.

◆ pressure_size

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

Definition at line 167 of file HydroMechanicsLocalAssemblerFracture.h.


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