OGS
HydroMechanicsLocalAssemblerFracture.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
20
21namespace ProcessLib
22{
23namespace LIE
24{
25namespace HydroMechanics
26{
27template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
28
29 int GlobalDim>
32{
33public:
38
40 MeshLib::Element const& e,
41 std::size_t const local_matrix_size,
42 std::vector<unsigned> const& dofIndex_to_localIndex,
43 NumLib::GenericIntegrationMethod const& integration_method,
44 bool const is_axially_symmetric,
46
47 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
48 double const /*t*/,
49 double const /*delta_t*/) override
50 {
51 for (auto& data : _ip_data)
52 {
53 data.pushBackState();
54 }
55 }
56
58 const double t, double const dt,
59 Eigen::VectorXd const& local_x) override;
60
61 std::vector<double> const& getIntPtSigma(
62 const double /*t*/,
63 std::vector<GlobalVector*> const& /*x*/,
64 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
65 std::vector<double>& cache) const override
66 {
67 cache.resize(0);
68 return cache;
69 }
70
71 std::vector<double> const& getIntPtEpsilon(
72 const double /*t*/,
73 std::vector<GlobalVector*> const& /*x*/,
74 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
75 std::vector<double>& cache) const override
76 {
77 cache.resize(0);
78 return cache;
79 }
80
81 std::vector<double> const& getIntPtDarcyVelocity(
82 const double /*t*/,
83 std::vector<GlobalVector*> const& /*x*/,
84 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
85 std::vector<double>& cache) const override
86 {
87 cache.resize(0);
88 return cache;
89 }
90
91 std::vector<double> const& getIntPtFractureVelocity(
92 const double t,
93 std::vector<GlobalVector*> const& x,
94 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
95 std::vector<double>& cache) const override;
96
97 std::vector<double> const& getIntPtFractureStress(
98 const double t,
99 std::vector<GlobalVector*> const& x,
100 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
101 std::vector<double>& cache) const override;
102
103 std::vector<double> const& getIntPtFractureAperture(
104 const double t,
105 std::vector<GlobalVector*> const& x,
106 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
107 std::vector<double>& cache) const override;
108
109 std::vector<double> const& getIntPtFracturePermeability(
110 const double t,
111 std::vector<GlobalVector*> const& x,
112 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
113 std::vector<double>& cache) const override;
114
115 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
116 const unsigned integration_point) const override
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 }
123
124private:
125 void assembleWithJacobianConcrete(double const t, double const dt,
126 Eigen::VectorXd const& local_x,
127 Eigen::VectorXd const& local_x_prev,
128 Eigen::VectorXd& local_b,
129 Eigen::MatrixXd& local_J) override;
130
132 double const t, double const dt,
133 Eigen::Ref<const Eigen::VectorXd> const& p,
134 Eigen::Ref<const Eigen::VectorXd> const& p_prev,
135 Eigen::Ref<const Eigen::VectorXd> const& g,
136 Eigen::Ref<const Eigen::VectorXd> const& g_prev,
137 Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_g,
138 Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pg,
139 Eigen::Ref<Eigen::MatrixXd> J_gg, Eigen::Ref<Eigen::MatrixXd> J_gp);
140
141 // Types for displacement.
147
148 // Types for pressure.
151
152 // Types for the integration point data
156 ShapeMatricesTypePressure, GlobalDim>;
157
158 using GlobalDimVectorType = Eigen::Matrix<double, GlobalDim, 1>;
159
161
162 std::vector<IntegrationPointDataType,
163 Eigen::aligned_allocator<IntegrationPointDataType>>
165
166 static const int pressure_index = 0;
167 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
168 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
169 static const int displacement_size =
170 ShapeFunctionDisplacement::NPOINTS * GlobalDim;
171
173 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
175};
176
177} // namespace HydroMechanics
178} // namespace LIE
179} // namespace ProcessLib
180
MatrixType< DisplacementDim, _number_of_dof > HMatrixType
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 & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
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
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)
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
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
std::vector< double > const & getIntPtDarcyVelocity(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
HydroMechanicsLocalAssemblerFracture(HydroMechanicsLocalAssemblerFracture const &)=delete
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
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > ShapeMatricesTypeDisplacement
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 postTimestepConcreteWithVector(const double t, double const dt, Eigen::VectorXd const &local_x) override
ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim > ShapeMatricesTypePressure
IntegrationPointDataFracture< HMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, GlobalDim > IntegrationPointDataType
HydroMechanicsLocalAssemblerFracture(HydroMechanicsLocalAssemblerFracture &&)=delete
Used for the extrapolation of the integration point values.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N