OGS
RichardsMechanicsFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
34
35namespace ProcessLib
36{
37namespace RichardsMechanics
38{
39namespace MPL = MaterialPropertyLib;
40
43template <typename ShapeMatrixType>
45{
46 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
47};
48
49template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
50 int DisplacementDim>
52 : public LocalAssemblerInterface<DisplacementDim>
53{
54public:
59
62
66
67 using IpData =
69 ShapeMatricesTypePressure, DisplacementDim,
70 ShapeFunctionDisplacement::NPOINTS>;
71
72 static int const KelvinVectorSize =
75
76 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
77
78 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
79 DisplacementDim,
81
83 delete;
85
87 MeshLib::Element const& e,
88 std::size_t const /*local_matrix_size*/,
89 NumLib::GenericIntegrationMethod const& integration_method,
90 bool const is_axially_symmetric,
92
93 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
94 double const t,
95 int const process_id) override;
96
97 void assemble(double const t, double const dt,
98 std::vector<double> const& local_x,
99 std::vector<double> const& local_x_prev,
100 std::vector<double>& local_M_data,
101 std::vector<double>& local_K_data,
102 std::vector<double>& local_rhs_data) override;
103
104 void assembleWithJacobian(double const t, double const dt,
105 std::vector<double> const& local_x,
106 std::vector<double> const& local_x_prev,
107 std::vector<double>& local_rhs_data,
108 std::vector<double>& local_Jac_data) override;
109
111 double const t, double const dt, Eigen::VectorXd const& local_x,
112 Eigen::VectorXd const& local_x_prev, int const process_id,
113 std::vector<double>& local_b_data,
114 std::vector<double>& local_Jac_data) override;
115
116 void initializeConcrete() override
117 {
118 unsigned const n_integration_points =
120
121 for (unsigned ip = 0; ip < n_integration_points; ip++)
122 {
123 auto& SD = this->current_states_[ip];
124 auto& ip_data = _ip_data[ip];
125
126 ParameterLib::SpatialPosition const x_position{
127 std::nullopt, this->element_.getID(), ip,
129 ShapeFunctionDisplacement,
131 ip_data.N_u))};
132
134 if (this->process_data_.initial_stress.value)
135 {
136 std::get<ProcessLib::ThermoRichardsMechanics::
137 ConstitutiveStress_StrainTemperature::
138 EffectiveStressData<DisplacementDim>>(SD)
139 .sigma_eff =
141 DisplacementDim>(
142 // The data in process_data_.initial_stress.value can
143 // be total stress or effective stress.
144 (*this->process_data_.initial_stress.value)(
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 this->solid_material_.initializeInternalStateVariables(
152 t, x_position,
153 *this->material_states_[ip].material_state_variables);
154
155 this->material_states_[ip].pushBackState();
156
157 this->prev_states_[ip] = SD;
158 }
159 }
160
162 double const t, double const dt, Eigen::VectorXd const& local_x,
163 Eigen::VectorXd const& local_x_prev) override;
164
165 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
166 const unsigned integration_point) const override
167 {
168 auto const& N_u = _secondary_data.N_u[integration_point];
169
170 // assumes N is stored contiguously in memory
171 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
172 }
173
174private:
194 double const t, double const dt, Eigen::VectorXd const& local_x,
195 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
196 std::vector<double>& local_Jac_data);
197
220 double const t, double const dt, Eigen::VectorXd const& local_x,
221 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
222 std::vector<double>& local_Jac_data);
223
224private:
226 double const t, double const dt,
227 ParameterLib::SpatialPosition const& x_position, IpData& ip_data,
228 MPL::VariableArray& variables, MPL::VariableArray& variables_prev,
229 MPL::Medium const* const medium, TemperatureData const T_data,
234 std::optional<MicroPorosityParameters> const& micro_porosity_parameters,
236 solid_material,
238 material_state_data);
239
240 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
241
243 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
245
246 static const int pressure_index = 0;
247 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
248 static const int displacement_index = ShapeFunctionPressure::NPOINTS;
249 static const int displacement_size =
250 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
251};
252
253} // namespace RichardsMechanics
254} // namespace ProcessLib
255
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
VectorType< _kelvin_vector_size > KelvinVectorType
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
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_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
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)
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
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)
std::tuple< StrainData< DisplacementDim >, ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature:: EffectiveStressData< DisplacementDim >, ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature:: SwellingDataStateful< DisplacementDim >, ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature:: MechanicalStrainData< DisplacementDim >, ProcessLib::ThermoRichardsMechanics::SaturationData, ProcessLib::ThermoRichardsMechanics::PorosityData, ProcessLib::ThermoRichardsMechanics::TransportPorosityData, MicroPressure, MicroSaturation > StatefulData
Data whose state must be tracked by the process.
ProcessLib::ConstitutiveRelations::PrevStateOf< StatefulData< DisplacementDim > > StatefulDataPrev
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
std::vector< StatefulData< DisplacementDim > > current_states_
std::vector< StatefulDataPrev< DisplacementDim > > prev_states_
RichardsMechanicsProcessData< DisplacementDim > & process_data_
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_
std::vector< ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > > material_states_
NumLib::GenericIntegrationMethod const & integration_method_
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u