OGS
HydroMechanicsFEM.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 <memory>
7#include <vector>
8
25
26namespace ProcessLib
27{
28namespace HydroMechanics
29{
30namespace MPL = MaterialPropertyLib;
31
32template <typename BMatricesType, typename ShapeMatricesTypeDisplacement,
33 typename ShapeMatricesTypePressure, int DisplacementDim, int NPoints>
35{
44
45 typename BMatricesType::KelvinVectorType sigma_eff, sigma_eff_prev;
46 typename BMatricesType::KelvinVectorType eps, eps_prev;
47
48 typename ShapeMatricesTypeDisplacement::NodalRowVectorType N_u;
49 typename ShapeMatricesTypeDisplacement::GlobalDimNodalMatrixType dNdx_u;
50
51 typename ShapeMatricesTypePressure::NodalRowVectorType N_p;
52 typename ShapeMatricesTypePressure::GlobalDimNodalMatrixType dNdx_p;
53
55 std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
56 DisplacementDim>::MaterialStateVariables>
59
60 // previous pressure rate for the fixed stress splitting
61 // approach in the staggered scheme.
62 // TODO: disable in monolithic scheme to save memory.
64
66 {
67 eps_prev = eps;
69 material_state_variables->pushBackState();
70 }
71
74 double const t,
75 ParameterLib::SpatialPosition const& x_position,
76 double const dt,
77 double const temperature)
78 {
79 namespace MPL = MaterialPropertyLib;
80
81 MPL::VariableArray variable_array;
82 MPL::VariableArray variable_array_prev;
83
84 auto const null_state = solid_material.createMaterialStateVariables();
85 solid_material.initializeInternalStateVariables(t, x_position,
86 *null_state);
87
89
90 variable_array.stress.emplace<KV>(KV::Zero());
91 variable_array.mechanical_strain.emplace<KV>(KV::Zero());
92 variable_array.temperature = temperature;
93
94 variable_array_prev.stress.emplace<KV>(KV::Zero());
95 variable_array_prev.mechanical_strain.emplace<KV>(KV::Zero());
96 variable_array_prev.temperature = temperature;
97
98 auto&& solution =
99 solid_material.integrateStress(variable_array_prev, variable_array,
100 t, x_position, dt, *null_state);
101
102 if (!solution)
103 {
104 OGS_FATAL("Computation of elastic tangent stiffness failed.");
105 }
106
108 std::move(std::get<2>(*solution));
109
110 return C;
111 }
112
113 template <typename DisplacementVectorType>
114 typename BMatricesType::KelvinMatrixType updateConstitutiveRelation(
115 MPL::VariableArray const& variable_array,
116 double const t,
117 ParameterLib::SpatialPosition const& x_position,
118 double const dt,
119 DisplacementVectorType const& /*u*/,
120 double const T)
121 {
122 MaterialPropertyLib::VariableArray variable_array_prev;
123 variable_array_prev.stress
126 variable_array_prev.mechanical_strain
128 eps_prev);
129 variable_array_prev.temperature = T;
130
131 auto&& solution = solid_material.integrateStress(
132 variable_array_prev, variable_array, t, x_position, dt,
134
135 if (!solution)
136 {
137 OGS_FATAL("Computation of local constitutive relation failed.");
138 }
139
141 std::tie(sigma_eff, material_state_variables, C) = std::move(*solution);
142
143 return C;
144 }
145
147};
148
151template <typename ShapeMatrixType>
153{
154 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
155};
156
157template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
158 int DisplacementDim>
160 : public LocalAssemblerInterface<DisplacementDim>
161{
162public:
165
166 // Types for pressure.
169
170 static int const KelvinVectorSize =
173
174 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
175
176 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
177 DisplacementDim,
179
182
184 MeshLib::Element const& e,
185 std::size_t const /*local_matrix_size*/,
186 NumLib::GenericIntegrationMethod const& integration_method,
187 bool const is_axially_symmetric,
189
191 std::size_t setIPDataInitialConditions(
192 std::string_view const name,
193 double const* values,
194 int const integration_order) override;
195
196 void assemble(double const /*t*/, double const /*dt*/,
197 std::vector<double> const& /*local_x*/,
198 std::vector<double> const& /*local_x_prev*/,
199 std::vector<double>& /*local_M_data*/,
200 std::vector<double>& /*local_K_data*/,
201 std::vector<double>& /*local_rhs_data*/) override
202 {
203 OGS_FATAL(
204 "HydroMechanicsLocalAssembler: assembly without jacobian is not "
205 "implemented.");
206 }
207
208 void assembleWithJacobian(double const t, double const dt,
209 std::vector<double> const& local_x,
210 std::vector<double> const& local_x_prev,
211 std::vector<double>& local_rhs_data,
212 std::vector<double>& local_Jac_data) override;
213
215 const double t, double const dt, Eigen::VectorXd const& local_x,
216 Eigen::VectorXd const& local_x_prev, int const process_id,
217 std::vector<double>& local_b_data,
218 std::vector<double>& local_Jac_data) override;
219
220 void initializeConcrete() override
221 {
222 unsigned const n_integration_points =
223 _integration_method.getNumberOfPoints();
224
225 for (unsigned ip = 0; ip < n_integration_points; ip++)
226 {
227 auto& ip_data = _ip_data[ip];
228
229 ParameterLib::SpatialPosition const x_position{
230 std::nullopt, _element.getID(),
232 NumLib::interpolateCoordinates<ShapeFunctionPressure,
234 _element, ip_data.N_p))};
235
237 if (_process_data.initial_stress.value)
238 {
239 ip_data.sigma_eff =
241 DisplacementDim>((*_process_data.initial_stress.value)(
242 std::numeric_limits<
243 double>::quiet_NaN() /* time independent */,
244 x_position));
245 }
246
247 double const t = 0; // TODO (naumov) pass t from top
248 ip_data.solid_material.initializeInternalStateVariables(
249 t, x_position, *ip_data.material_state_variables);
250
251 ip_data.pushBackState();
252 }
253 }
254
255 void postTimestepConcrete(Eigen::VectorXd const& local_x,
256 Eigen::VectorXd const& local_x_prev,
257 double const t, double const dt,
258 int const process_id) override;
259
261 double const t, double const dt, Eigen::VectorXd const& local_xs,
262 Eigen::VectorXd const& local_x_prev) override;
263
264 void postNonLinearSolverConcrete(Eigen::VectorXd const& local_x,
265 Eigen::VectorXd const& local_x_prev,
266 double const t, double const dt,
267 int const process_id) override;
268
269 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
270 double const t,
271 int const process_id) override;
272
273 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
274 const unsigned integration_point) const override
275 {
276 auto const& N_u = _secondary_data.N_u[integration_point];
277
278 // assumes N is stored contiguously in memory
279 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
280 }
281
282 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
283 // the ordering of the cache_mat.
284 // There should be only one.
285 std::vector<double> getSigma() const override;
286
287 std::vector<double> getEpsilon() const override;
288
289 std::vector<double> getStrainRateVariable() const override;
290
291 std::vector<double> const& getIntPtDarcyVelocity(
292 const double t,
293 std::vector<GlobalVector*> const& x,
294 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
295 std::vector<double>& cache) const override;
296
297 std::vector<double> const& getIntPtSigma(
298 const double /*t*/,
299 std::vector<GlobalVector*> const& /*x*/,
300 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
301 std::vector<double>& cache) const override
302 {
304 _ip_data, &IpData::sigma_eff, cache);
305 }
306
307 std::vector<double> const& getIntPtEpsilon(
308 const double /*t*/,
309 std::vector<GlobalVector*> const& /*x*/,
310 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
311 std::vector<double>& cache) const override
312 {
314 _ip_data, &IpData::eps, cache);
315 }
316
318 {
319 return displacement_size;
320 }
321
322private:
342 const double t, double const dt, Eigen::VectorXd const& local_x,
343 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
344
369 const double t, double const dt, Eigen::VectorXd const& local_x,
370 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
371 std::vector<double>& local_Jac_data);
372
373 unsigned getNumberOfIntegrationPoints() const override;
374
375 int getMaterialID() const override;
376
378 DisplacementDim>::MaterialStateVariables const&
379 getMaterialStateVariablesAt(unsigned integration_point) const override;
380
381private:
383
386 using IpData =
388 ShapeMatricesTypePressure, DisplacementDim,
389 ShapeFunctionDisplacement::NPOINTS>;
390 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
391
396 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
398
399 static const int pressure_index = 0;
400 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
401 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
402 static const int displacement_size =
403 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
404};
405
406} // namespace HydroMechanics
407} // namespace ProcessLib
408
#define OGS_FATAL(...)
Definition Error.h:19
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, 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.
NumLib::GenericIntegrationMethod const & _integration_method
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
HydroMechanicsLocalAssembler(HydroMechanicsLocalAssembler &&)=delete
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
void postNonLinearSolverConcrete(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const process_id) override
void postTimestepConcrete(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const process_id) override
void assembleWithJacobianForPressureEquations(const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
Returns number of read integration points.
IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS > IpData
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_xs, Eigen::VectorXd const &local_x_prev) override
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
void assembleWithJacobianForStaggeredScheme(const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
void assembleWithJacobianForDeformationEquations(const double t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
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
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
HydroMechanicsLocalAssembler(HydroMechanicsLocalAssembler const &)=delete
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
HydroMechanicsProcessData< DisplacementDim > & _process_data
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
MathLib::KelvinVector::Invariants< KelvinVectorSize > Invariants
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
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)
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
BMatricesType::KelvinMatrixType updateConstitutiveRelation(MPL::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x_position, double const dt, DisplacementVectorType const &, double const T)
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > computeElasticTangentStiffness(double const t, ParameterLib::SpatialPosition const &x_position, double const dt, double const temperature)
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u