OGS
HydroMechanicsLocalAssemblerMatrix.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <optional>
14#include <vector>
15
23
24namespace ProcessLib
25{
26namespace LIE
27{
28namespace HydroMechanics
29{
30namespace MPL = MaterialPropertyLib;
31
32template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
33 int GlobalDim>
36{
37public:
41 delete;
42
44 MeshLib::Element const& e,
45 std::size_t const n_variables,
46 std::size_t const local_matrix_size,
47 std::vector<unsigned> const& dofIndex_to_localIndex,
48 NumLib::GenericIntegrationMethod const& integration_method,
49 bool const is_axially_symmetric,
51
52 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
53 double const /*t*/,
54 double const /*delta_t*/) override
55 {
56 for (auto& data : _ip_data)
57 {
58 data.pushBackState();
59 }
60 }
61
62 std::vector<double> const& getIntPtSigma(
63 const double t,
64 std::vector<GlobalVector*> const& x,
65 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
66 std::vector<double>& cache) const override;
67
68 std::vector<double> const& getIntPtEpsilon(
69 const double t,
70 std::vector<GlobalVector*> const& x,
71 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
72 std::vector<double>& cache) const override;
73
74 std::vector<double> const& getIntPtDarcyVelocity(
75 const double t,
76 std::vector<GlobalVector*> const& x,
77 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
78 std::vector<double>& cache) const override;
79
80 std::vector<double> const& getIntPtFractureVelocity(
81 const double /*t*/,
82 std::vector<GlobalVector*> const& /*x*/,
83 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
84 std::vector<double>& cache) const override
85 {
86 cache.resize(0);
87 return cache;
88 }
89
90 std::vector<double> const& getIntPtFractureStress(
91 const double /*t*/,
92 std::vector<GlobalVector*> const& /*x*/,
93 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
94 std::vector<double>& cache) const override
95 {
96 cache.resize(0);
97 return cache;
98 }
99
100 std::vector<double> const& getIntPtFractureAperture(
101 const double /*t*/,
102 std::vector<GlobalVector*> const& /*x*/,
103 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
104 std::vector<double>& cache) const override
105 {
106 cache.resize(0);
107 return cache;
108 }
109
110 std::vector<double> const& getIntPtFracturePermeability(
111 const double /*t*/,
112 std::vector<GlobalVector*> const& /*x*/,
113 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
114 std::vector<double>& cache) const override
115 {
116 cache.resize(0);
117 return cache;
118 }
119
120 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
121 const unsigned integration_point) const override
122 {
123 auto const& N = _secondary_data.N[integration_point];
124
125 // assumes N is stored contiguously in memory
126 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
127 }
128
129protected:
130 void assembleWithJacobianConcrete(double const t, double const dt,
131 Eigen::VectorXd const& local_x,
132 Eigen::VectorXd const& local_x_prev,
133 Eigen::VectorXd& local_rhs,
134 Eigen::MatrixXd& local_Jac) override;
135
137 double const t, double const dt,
138 Eigen::Ref<const Eigen::VectorXd> const& p,
139 Eigen::Ref<const Eigen::VectorXd> const& p_prev,
140 Eigen::Ref<const Eigen::VectorXd> const& u,
141 Eigen::Ref<const Eigen::VectorXd> const& u_prev,
142 Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_u,
143 Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pu,
144 Eigen::Ref<Eigen::MatrixXd> J_uu, Eigen::Ref<Eigen::MatrixXd> J_up);
145
147 double const t, double const dt,
148 Eigen::VectorXd const& local_x) override;
149
151 double const t, double const dt,
152 Eigen::Ref<const Eigen::VectorXd> const& p,
153 Eigen::Ref<const Eigen::VectorXd> const& u);
154
156 double const t, Eigen::Ref<Eigen::VectorXd> p);
157
158 // Types for displacement.
163
165
166 // Types for pressure.
169
172 ShapeMatricesTypePressure, GlobalDim,
173 ShapeFunctionDisplacement::NPOINTS>;
174
175 using GlobalDimVector = Eigen::Matrix<double, GlobalDim, 1>;
176
177 std::optional<BBarMatrixType> getDilatationalBBarMatrix() const
178 {
179 if (!(_process_data.use_b_bar))
180 {
181 return std::nullopt;
182 }
183
185 GlobalDim, ShapeFunctionDisplacement::NPOINTS,
186 ShapeFunctionDisplacement, BBarMatrixType,
190 }
191
193
194 std::vector<IntegrationPointDataType,
195 Eigen::aligned_allocator<IntegrationPointDataType>>
197
198 static const int pressure_index = 0;
199 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
200 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
201 static const int displacement_size =
202 ShapeFunctionDisplacement::NPOINTS * GlobalDim;
203 static const int kelvin_vector_size =
205
207 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
209};
210
211} // namespace HydroMechanics
212} // namespace LIE
213} // namespace ProcessLib
214
MatrixType< 3, ShapeFunction::NPOINTS > BBarMatrixType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > ShapeMatricesTypeDisplacement
BMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > BMatricesType
std::vector< double > const & getIntPtFractureVelocity(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 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, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
IntegrationPointDataMatrix< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, GlobalDim, ShapeFunctionDisplacement::NPOINTS > IntegrationPointDataType
std::vector< double > const & getIntPtFractureStress(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void postTimestepConcreteWithBlockVectors(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assembleWithJacobianConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, Eigen::VectorXd &local_rhs, Eigen::MatrixXd &local_Jac) override
std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim > ShapeMatricesTypePressure
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 &u, Eigen::Ref< const Eigen::VectorXd > const &u_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_u, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pu, Eigen::Ref< Eigen::MatrixXd > J_uu, Eigen::Ref< Eigen::MatrixXd > J_up)
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix const &)=delete
void postTimestepConcreteWithVector(double const t, double const dt, Eigen::VectorXd const &local_x) override
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix &&)=delete
std::vector< double > const & getIntPtFractureAperture(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
BBarMatrixType computeDilatationalBbar(std::vector< IpData, Eigen::aligned_allocator< IpData > > const &ip_data, MeshLib::Element const &element, NumLib::GenericIntegrationMethod const &integration_method, const bool is_axially_symmetric)
Used for the extrapolation of the integration point values.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N