OGS
RichardsMechanicsFEM.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
19#include "NumLib/DOF/LocalDOF.h"
28
29namespace ProcessLib
30{
31namespace RichardsMechanics
32{
33namespace MPL = MaterialPropertyLib;
34
37template <typename ShapeMatrixType>
39{
40 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
41};
42
43template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
44 int DisplacementDim>
46 : public LocalAssemblerInterface<DisplacementDim>
47{
48public:
53
56
60
61 using IpData =
63 ShapeMatricesTypePressure, DisplacementDim,
64 ShapeFunctionDisplacement::NPOINTS>;
65
66 static int const KelvinVectorSize =
69
70 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
71
72 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
73 DisplacementDim,
75
77 delete;
79
81 MeshLib::Element const& e,
82 std::size_t const /*local_matrix_size*/,
83 NumLib::GenericIntegrationMethod const& integration_method,
84 bool const is_axially_symmetric,
86
87 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
88 double const t,
89 int const process_id) override;
90
91 void assemble(double const t, double const dt,
92 std::vector<double> const& local_x,
93 std::vector<double> const& local_x_prev,
94 std::vector<double>& local_M_data,
95 std::vector<double>& local_K_data,
96 std::vector<double>& local_rhs_data) override;
97
98 void assembleWithJacobian(double const t, double const dt,
99 std::vector<double> const& local_x,
100 std::vector<double> const& local_x_prev,
101 std::vector<double>& local_rhs_data,
102 std::vector<double>& local_Jac_data) override;
103
105 double const t, double const dt, Eigen::VectorXd const& local_x,
106 Eigen::VectorXd const& local_x_prev, int const process_id,
107 std::vector<double>& local_b_data,
108 std::vector<double>& local_Jac_data) override;
109
110 void initializeConcrete() override
111 {
112 unsigned const n_integration_points =
113 this->integration_method_.getNumberOfPoints();
114
115 for (unsigned ip = 0; ip < n_integration_points; ip++)
116 {
117 auto& SD = this->current_states_[ip];
118 auto& ip_data = ip_data_[ip];
119
120 ParameterLib::SpatialPosition const x_position{
121 std::nullopt, this->element_.getID(),
123 ShapeFunctionDisplacement,
125 ip_data.N_u))};
126
128 if (this->process_data_.initial_stress.value)
129 {
131 DisplacementDim>>(SD)
132 .sigma_eff =
134 DisplacementDim>(
135 // The data in process_data_.initial_stress.value can
136 // be total stress or effective stress.
137 (*this->process_data_.initial_stress.value)(
138 std::numeric_limits<
139 double>::quiet_NaN() /* time independent */,
140 x_position));
141 }
142
143 double const t = 0; // TODO (naumov) pass t from top
144 this->solid_material_.initializeInternalStateVariables(
145 t, x_position,
146 *this->material_states_[ip].material_state_variables);
147
148 this->material_states_[ip].pushBackState();
149
150 this->prev_states_[ip] = SD;
151 }
152 }
153
155 double const t, double const dt, Eigen::VectorXd const& local_x,
156 Eigen::VectorXd const& local_x_prev) override;
157
158 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
159 const unsigned integration_point) const override
160 {
161 auto const& N_u = secondary_data_.N_u[integration_point];
162
163 // assumes N is stored contiguously in memory
164 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
165 }
166
168 {
169 return displacement_size;
170 }
171
172private:
192 double const t, double const dt, Eigen::VectorXd const& local_x,
193 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
194 std::vector<double>& local_Jac_data);
195
218 double const t, double const dt, Eigen::VectorXd const& local_x,
219 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
220 std::vector<double>& local_Jac_data);
221
222private:
224 double const t, double const dt,
225 ParameterLib::SpatialPosition const& x_position, IpData& ip_data,
226 MPL::VariableArray& variables, MPL::VariableArray& variables_prev,
227 MPL::Medium const* const medium, TemperatureData const T_data,
232 std::optional<MicroPorosityParameters> const& micro_porosity_parameters,
234 solid_material,
236 material_state_data);
237
238 static constexpr auto localDOF(auto const& x)
239 {
240 return NumLib::localDOF<
241 ShapeFunctionPressure,
243 }
244
245 std::vector<IpData, Eigen::aligned_allocator<IpData>> ip_data_;
246
248 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
250
251 static const int pressure_index = 0;
252 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
253 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
254 static const int displacement_size =
255 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
256};
257
258} // namespace RichardsMechanics
259} // namespace ProcessLib
260
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS > IpData
void assembleWithJacobianForPressureEquations(double const 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)
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
MathLib::KelvinVector::Invariants< KelvinVectorSize > Invariants
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > secondary_data_
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
RichardsMechanicsLocalAssembler(RichardsMechanicsLocalAssembler &&)=delete
static void assembleWithJacobianEvalConstitutiveSetting(double const t, double const dt, ParameterLib::SpatialPosition const &x_position, IpData &ip_data, MPL::VariableArray &variables, MPL::VariableArray &variables_prev, MPL::Medium const *const medium, TemperatureData const T_data, CapillaryPressureData< DisplacementDim > const &p_cap_data, ConstitutiveData< DisplacementDim > &CD, StatefulData< DisplacementDim > &SD, StatefulDataPrev< DisplacementDim > const &SD_prev, std::optional< MicroPorosityParameters > const &micro_porosity_parameters, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material, ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > &material_state_data)
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.
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
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_b_data, std::vector< double > &local_Jac_data) override
RichardsMechanicsLocalAssembler(RichardsMechanicsLocalAssembler const &)=delete
void assembleWithJacobianForDeformationEquations(double const 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::vector< IpData, Eigen::aligned_allocator< IpData > > ip_data_
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
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
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)
auto localDOF(ElementDOFVector const &x)
Definition LocalDOF.h:56
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BaseLib::StrongType< double, struct TemperatureDataTag > TemperatureData
ProcessLib::ConstitutiveRelations::PrevStateOf< StatefulData< DisplacementDim > > StatefulDataPrev
std::tuple< StrainData< DisplacementDim >, ProcessLib::ConstitutiveRelations::EffectiveStressData< DisplacementDim >, ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature:: SwellingDataStateful< DisplacementDim >, ProcessLib::ConstitutiveRelations::MechanicalStrainData< DisplacementDim >, ProcessLib::ThermoRichardsMechanics::SaturationData, ProcessLib::ThermoRichardsMechanics::PorosityData, ProcessLib::ThermoRichardsMechanics::TransportPorosityData, MicroPressure, MicroSaturation > StatefulData
Data whose state must be tracked by the process.
std::tuple< StiffnessTensor< DisplacementDim >, ProcessLib::ThermoRichardsMechanics::PorosityData, Density, LiquidDensity, ProcessLib::ThermoRichardsMechanics::BiotData, ProcessLib::ThermoRichardsMechanics::SaturationDataDeriv, ProcessLib::ThermoRichardsMechanics::LiquidViscosityData, ProcessLib::ThermoRichardsMechanics::SolidCompressibilityData, ProcessLib::ThermoRichardsMechanics::BishopsData, PrevState< ProcessLib::ThermoRichardsMechanics::BishopsData >, ProcessLib::ThermoRichardsMechanics::PermeabilityData< DisplacementDim >, SaturationSecantDerivative > ConstitutiveData
Data that is needed for the equation system assembly.
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, RichardsMechanicsProcessData< DisplacementDim > &process_data)
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_
std::vector< ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > > material_states_
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u