OGS
HydroMechanicsLocalAssemblerMatrix.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <optional>
7#include <vector>
8
16
17namespace ProcessLib
18{
19namespace LIE
20{
21namespace HydroMechanics
22{
23namespace MPL = MaterialPropertyLib;
24
25template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
26 int DisplacementDim>
29{
30public:
34 delete;
35
37 MeshLib::Element const& e,
38 std::size_t const n_variables,
39 std::size_t const local_matrix_size,
40 std::vector<unsigned> const& dofIndex_to_localIndex,
41 NumLib::GenericIntegrationMethod const& integration_method,
42 bool const is_axially_symmetric,
44
45 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
46 double const /*t*/,
47 double const /*delta_t*/) override
48 {
49 for (auto& data : _ip_data)
50 {
51 data.pushBackState();
52 }
53 }
54
55 std::vector<double> const& getIntPtSigma(
56 const double t,
57 std::vector<GlobalVector*> const& x,
58 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
59 std::vector<double>& cache) const override;
60
61 std::vector<double> const& getIntPtEpsilon(
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 std::vector<double> const& getIntPtDarcyVelocity(
68 const double t,
69 std::vector<GlobalVector*> const& x,
70 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
71 std::vector<double>& cache) const override;
72
73 std::vector<double> const& getIntPtFractureVelocity(
74 const double /*t*/,
75 std::vector<GlobalVector*> const& /*x*/,
76 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
77 std::vector<double>& cache) const override
78 {
79 cache.resize(0);
80 return cache;
81 }
82
83 std::vector<double> const& getIntPtFractureStress(
84 const double /*t*/,
85 std::vector<GlobalVector*> const& /*x*/,
86 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
87 std::vector<double>& cache) const override
88 {
89 cache.resize(0);
90 return cache;
91 }
92
93 std::vector<double> const& getIntPtFractureAperture(
94 const double /*t*/,
95 std::vector<GlobalVector*> const& /*x*/,
96 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
97 std::vector<double>& cache) const override
98 {
99 cache.resize(0);
100 return cache;
101 }
102
103 std::vector<double> const& getIntPtFracturePermeability(
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 cache.resize(0);
110 return cache;
111 }
112
113 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
114 const unsigned integration_point) const override
115 {
116 auto const& N = _secondary_data.N[integration_point];
117
118 // assumes N is stored contiguously in memory
119 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
120 }
121
122protected:
123 void assembleWithJacobianConcrete(double const t, double const dt,
124 Eigen::VectorXd const& local_x,
125 Eigen::VectorXd const& local_x_prev,
126 Eigen::VectorXd& local_rhs,
127 Eigen::MatrixXd& local_Jac) override;
128
130 double const t, double const dt,
131 Eigen::Ref<const Eigen::VectorXd> const& p,
132 Eigen::Ref<const Eigen::VectorXd> const& p_prev,
133 Eigen::Ref<const Eigen::VectorXd> const& u,
134 Eigen::Ref<const Eigen::VectorXd> const& u_prev,
135 Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_u,
136 Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pu,
137 Eigen::Ref<Eigen::MatrixXd> J_uu, Eigen::Ref<Eigen::MatrixXd> J_up);
138
140 double const t, double const dt,
141 Eigen::VectorXd const& local_x) override;
142
144 double const t, double const dt,
145 Eigen::Ref<const Eigen::VectorXd> const& p,
146 Eigen::Ref<const Eigen::VectorXd> const& u);
147
149 double const t, Eigen::Ref<Eigen::VectorXd> p);
150
151 // Types for displacement.
156
158
159 // Types for pressure.
162
165 ShapeMatricesTypePressure, DisplacementDim,
166 ShapeFunctionDisplacement::NPOINTS>;
169
170 using DisplacementDimVector = Eigen::Matrix<double, DisplacementDim, 1>;
171
172 std::optional<BBarMatrixType> getDilatationalBBarMatrix() const
173 {
174 if (!(_process_data.use_b_bar))
175 {
176 return std::nullopt;
177 }
178
180 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
181 ShapeFunctionDisplacement, BBarMatrixType,
185 }
186
188
189 std::vector<IntegrationPointDataType,
190 Eigen::aligned_allocator<IntegrationPointDataType>>
192
193 static const int pressure_index = 0;
194 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
195 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
196 static const int displacement_size =
197 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
198 static const int kelvin_vector_size =
200
202 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
204};
205
206} // namespace HydroMechanics
207} // namespace LIE
208} // namespace ProcessLib
209
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
MatrixType< 3, ShapeFunctionDisplacement::NPOINTS > BBarMatrixType
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< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
void postTimestepConcreteWithVector(double const t, double const dt, Eigen::VectorXd const &local_x) override
std::vector< double > const & getIntPtFracturePermeability(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
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
void postTimestepConcreteWithBlockVectors(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)
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, DisplacementDim > ShapeMatricesTypePressure
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
std::vector< double > const & getIntPtFractureAperture(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureStress(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 & getIntPtFractureVelocity(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, 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.
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
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)
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix &&)=delete
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix const &)=delete
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
IntegrationPointDataMatrix< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS > IntegrationPointDataType
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)
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
Used for the extrapolation of the integration point values.