OGS
RichardsMechanicsFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
32
33namespace ProcessLib
34{
35namespace RichardsMechanics
36{
37namespace MPL = MaterialPropertyLib;
38
41template <typename ShapeMatrixType>
43{
44 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
45};
46
47template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
48 int DisplacementDim>
50 : public LocalAssemblerInterface<DisplacementDim>
51{
52public:
57
60
64
65 using IpData =
67 ShapeMatricesTypePressure, DisplacementDim,
68 ShapeFunctionDisplacement::NPOINTS>;
69
70 static int const KelvinVectorSize =
73
74 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
75
76 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
77 DisplacementDim,
79
81 delete;
83
85 MeshLib::Element const& e,
86 std::size_t const /*local_matrix_size*/,
87 NumLib::GenericIntegrationMethod const& integration_method,
88 bool const is_axially_symmetric,
90
93 std::string_view const name,
94 double const* values,
95 int const integration_order) override;
96
97 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
98 double const t,
99 int const process_id) override;
100
101 void assemble(double const t, double const dt,
102 std::vector<double> const& local_x,
103 std::vector<double> const& local_x_prev,
104 std::vector<double>& local_M_data,
105 std::vector<double>& local_K_data,
106 std::vector<double>& local_rhs_data) override;
107
108 void assembleWithJacobian(double const t, double const dt,
109 std::vector<double> const& local_x,
110 std::vector<double> const& local_x_prev,
111 std::vector<double>& /*local_M_data*/,
112 std::vector<double>& /*local_K_data*/,
113 std::vector<double>& local_rhs_data,
114 std::vector<double>& local_Jac_data) override;
115
117 double const t, double const dt, Eigen::VectorXd const& local_x,
118 Eigen::VectorXd const& local_x_prev, int const process_id,
119 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
120 std::vector<double>& local_b_data,
121 std::vector<double>& local_Jac_data) override;
122
123 void initializeConcrete() override
124 {
125 unsigned const n_integration_points =
127
128 for (unsigned ip = 0; ip < n_integration_points; ip++)
129 {
130 auto& ip_data = _ip_data[ip];
131
132 ParameterLib::SpatialPosition const x_position{
133 std::nullopt, _element.getID(), ip,
136 ShapeFunctionDisplacement,
138
140 if (_process_data.initial_stress != nullptr)
141 {
142 ip_data.sigma_eff =
144 DisplacementDim>((*_process_data.initial_stress)(
145 std::numeric_limits<
146 double>::quiet_NaN() /* time independent */,
147 x_position));
148 }
149
150 double const t = 0; // TODO (naumov) pass t from top
151 ip_data.solid_material.initializeInternalStateVariables(
152 t, x_position, *ip_data.material_state_variables);
153
154 ip_data.pushBackState();
155 }
156 }
157
158 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
159 Eigen::VectorXd const& /*local_x_prev*/,
160 double const /*t*/, double const /*dt*/,
161 int const /*process_id*/) override
162 {
163 unsigned const n_integration_points =
165
166 for (unsigned ip = 0; ip < n_integration_points; ip++)
167 {
168 _ip_data[ip].pushBackState();
169 }
170 }
171
173 double const t, double const dt, Eigen::VectorXd const& local_x,
174 Eigen::VectorXd const& local_x_prev) override;
175
176 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
177 const unsigned integration_point) const override
178 {
179 auto const& N_u = _secondary_data.N_u[integration_point];
180
181 // assumes N is stored contiguously in memory
182 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
183 }
184
185 std::vector<double> getSigma() const override;
186
187 std::vector<double> const& getIntPtDarcyVelocity(
188 const double t,
189 std::vector<GlobalVector*> const& x,
190 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
191 std::vector<double>& cache) const override;
192
193 std::vector<double> getSaturation() const override;
194 std::vector<double> const& getIntPtSaturation(
195 const double t,
196 std::vector<GlobalVector*> const& x,
197 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
198 std::vector<double>& cache) const override;
199
200 std::vector<double> getMicroSaturation() const override;
201 std::vector<double> const& getIntPtMicroSaturation(
202 const double t,
203 std::vector<GlobalVector*> const& x,
204 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
205 std::vector<double>& cache) const override;
206
207 std::vector<double> getMicroPressure() const override;
208 std::vector<double> const& getIntPtMicroPressure(
209 const double t,
210 std::vector<GlobalVector*> const& x,
211 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
212 std::vector<double>& cache) const override;
213
214 std::vector<double> getPorosity() const override;
215 std::vector<double> const& getIntPtPorosity(
216 const double t,
217 std::vector<GlobalVector*> const& x,
218 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
219 std::vector<double>& cache) const override;
220
221 std::vector<double> getTransportPorosity() const override;
222 std::vector<double> const& getIntPtTransportPorosity(
223 const double t,
224 std::vector<GlobalVector*> const& x,
225 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
226 std::vector<double>& cache) const override;
227
228 std::vector<double> const& getIntPtSigma(
229 const double t,
230 std::vector<GlobalVector*> const& x,
231 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
232 std::vector<double>& cache) const override;
233
234 std::vector<double> getSwellingStress() const override;
235 std::vector<double> const& getIntPtSwellingStress(
236 const double t,
237 std::vector<GlobalVector*> const& x,
238 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
239 std::vector<double>& cache) const override;
240
241 std::vector<double> getEpsilon() const override;
242 std::vector<double> const& getIntPtEpsilon(
243 const double t,
244 std::vector<GlobalVector*> const& x,
245 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
246 std::vector<double>& cache) const override;
247
248 int getMaterialID() const override;
249
250 std::vector<double> getMaterialStateVariableInternalState(
251 std::function<std::span<double>(
253 MaterialStateVariables&)> const& get_values_span,
254 int const& n_components) const override;
255
256 std::vector<double> const& getIntPtDryDensitySolid(
257 const double t,
258 std::vector<GlobalVector*> const& x,
259 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
260 std::vector<double>& cache) const override;
261
262private:
287 double const t, double const dt, Eigen::VectorXd const& local_x,
288 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_M_data,
289 std::vector<double>& local_K_data, std::vector<double>& local_b_data,
290 std::vector<double>& local_Jac_data);
291
319 double const t, double const dt, Eigen::VectorXd const& local_x,
320 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_M_data,
321 std::vector<double>& local_K_data, std::vector<double>& local_b_data,
322 std::vector<double>& local_Jac_data);
323
324 unsigned getNumberOfIntegrationPoints() const override;
325
327 DisplacementDim>::MaterialStateVariables const&
328 getMaterialStateVariablesAt(unsigned integration_point) const override;
329
330private:
332
333 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
334
339 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
341
342 static const int pressure_index = 0;
343 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
344 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
345 static const int displacement_size =
346 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
347};
348
349} // namespace RichardsMechanics
350} // namespace ProcessLib
351
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
VectorType< _kelvin_vector_size > KelvinVectorType
std::vector< double > const & getIntPtPorosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
void assembleWithJacobianForDeformationEquations(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
std::vector< double > const & getIntPtMicroPressure(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
void assembleWithJacobianForPressureEquations(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
RichardsMechanicsLocalAssembler(RichardsMechanicsLocalAssembler &&)=delete
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &, std::vector< double > &, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) 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
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
RichardsMechanicsLocalAssembler(RichardsMechanicsLocalAssembler const &)=delete
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
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
std::vector< double > const & getIntPtDryDensitySolid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data) override
std::vector< double > const & getIntPtTransportPorosity(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 & getIntPtMicroSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
RichardsMechanicsProcessData< DisplacementDim > & _process_data
std::vector< double > getMaterialStateVariableInternalState(std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const override
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
std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
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_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
std::vector< double > const & getIntPtSwellingStress(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u