OGS
HydroMechanicsLocalAssemblerMatrix.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
21
22namespace ProcessLib
23{
24namespace LIE
25{
26namespace HydroMechanics
27{
28namespace MPL = MaterialPropertyLib;
29
30template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
31 int GlobalDim>
34{
35public:
39 delete;
40
42 MeshLib::Element const& e,
43 std::size_t const n_variables,
44 std::size_t const local_matrix_size,
45 std::vector<unsigned> const& dofIndex_to_localIndex,
46 NumLib::GenericIntegrationMethod const& integration_method,
47 bool const is_axially_symmetric,
49
50 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
51 double const /*t*/,
52 double const /*delta_t*/) override
53 {
54 for (auto& data : _ip_data)
55 {
56 data.pushBackState();
57 }
58 }
59
60 std::vector<double> const& getIntPtSigma(
61 const double t,
62 std::vector<GlobalVector*> const& x,
63 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
64 std::vector<double>& cache) const override;
65
66 std::vector<double> const& getIntPtEpsilon(
67 const double t,
68 std::vector<GlobalVector*> const& x,
69 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
70 std::vector<double>& cache) const override;
71
72 std::vector<double> const& getIntPtDarcyVelocity(
73 const double t,
74 std::vector<GlobalVector*> const& x,
75 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
76 std::vector<double>& cache) const override;
77
78 std::vector<double> const& getIntPtFractureVelocity(
79 const double /*t*/,
80 std::vector<GlobalVector*> const& /*x*/,
81 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
82 std::vector<double>& cache) const override
83 {
84 cache.resize(0);
85 return cache;
86 }
87
88 std::vector<double> const& getIntPtFractureStress(
89 const double /*t*/,
90 std::vector<GlobalVector*> const& /*x*/,
91 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
92 std::vector<double>& cache) const override
93 {
94 cache.resize(0);
95 return cache;
96 }
97
98 std::vector<double> const& getIntPtFractureAperture(
99 const double /*t*/,
100 std::vector<GlobalVector*> const& /*x*/,
101 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
102 std::vector<double>& cache) const override
103 {
104 cache.resize(0);
105 return cache;
106 }
107
108 std::vector<double> const& getIntPtFracturePermeability(
109 const double /*t*/,
110 std::vector<GlobalVector*> const& /*x*/,
111 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
112 std::vector<double>& cache) const override
113 {
114 cache.resize(0);
115 return cache;
116 }
117
118 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
119 const unsigned integration_point) const override
120 {
121 auto const& N = _secondary_data.N[integration_point];
122
123 // assumes N is stored contiguously in memory
124 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
125 }
126
127protected:
128 void assembleWithJacobianConcrete(double const t, double const dt,
129 Eigen::VectorXd const& local_x,
130 Eigen::VectorXd const& local_x_prev,
131 Eigen::VectorXd& local_rhs,
132 Eigen::MatrixXd& local_Jac) override;
133
135 double const t, double const dt,
136 Eigen::Ref<const Eigen::VectorXd> const& p,
137 Eigen::Ref<const Eigen::VectorXd> const& p_prev,
138 Eigen::Ref<const Eigen::VectorXd> const& u,
139 Eigen::Ref<const Eigen::VectorXd> const& u_prev,
140 Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_u,
141 Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pu,
142 Eigen::Ref<Eigen::MatrixXd> J_uu, Eigen::Ref<Eigen::MatrixXd> J_up);
143
145 double const t, double const dt,
146 Eigen::VectorXd const& local_x) override;
147
149 double const t, double const dt,
150 Eigen::Ref<const Eigen::VectorXd> const& p,
151 Eigen::Ref<const Eigen::VectorXd> const& u);
152
154 double const t, Eigen::Ref<Eigen::VectorXd> p);
155
156 // Types for displacement.
161
162 // Types for pressure.
165
168 ShapeMatricesTypePressure, GlobalDim,
169 ShapeFunctionDisplacement::NPOINTS>;
170
171 using GlobalDimVector = Eigen::Matrix<double, GlobalDim, 1>;
172
174
175 std::vector<IntegrationPointDataType,
176 Eigen::aligned_allocator<IntegrationPointDataType>>
178
179 static const int pressure_index = 0;
180 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
181 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
182 static const int displacement_size =
183 ShapeFunctionDisplacement::NPOINTS * GlobalDim;
184 static const int kelvin_vector_size =
186
188 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
190};
191
192} // namespace HydroMechanics
193} // namespace LIE
194} // namespace ProcessLib
195
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.
Used for the extrapolation of the integration point values.