OGS
TH2MFEM.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
28
29namespace ProcessLib
30{
31namespace TH2M
32{
35template <typename ShapeMatrixType>
37{
38 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
39};
40
41template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
42 int DisplacementDim>
43class TH2MLocalAssembler : public LocalAssemblerInterface<DisplacementDim>
44{
45public:
48
51
52 template <int N>
53 using VectorType =
54 typename ShapeMatricesTypePressure::template VectorType<N>;
55
56 template <int M, int N>
57 using MatrixType =
58 typename ShapeMatricesTypePressure::template MatrixType<M, N>;
59
62
65
66 static int const KelvinVectorSize =
68 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
69
70 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
71 DisplacementDim,
73
75
78
80 MeshLib::Element const& e,
81 std::size_t const /*local_matrix_size*/,
82 NumLib::GenericIntegrationMethod const& integration_method,
83 bool const is_axially_symmetric,
85
86private:
89 std::string_view const name,
90 double const* values,
91 int const integration_order) override;
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
110 void initializeConcrete() override
111 {
112 unsigned const n_integration_points =
113 this->integration_method_.getNumberOfPoints();
114 auto const time_independent = std::numeric_limits<double>::quiet_NaN();
115 auto const& medium =
116 *this->process_data_.media_map.getMedium(this->element_.getID());
117
118 for (unsigned ip = 0; ip < n_integration_points; ip++)
119 {
120 auto& ip_data = _ip_data[ip];
121
122 ParameterLib::SpatialPosition const x_position{
123 std::nullopt, this->element_.getID(),
125 ShapeFunctionDisplacement,
127 ip_data.N_u))};
128 auto& current_state = this->current_states_[ip];
129
131 if (this->process_data_.initial_stress.value)
132 {
133 current_state.eff_stress_data.sigma_eff.noalias() =
135 DisplacementDim>(
136 (*this->process_data_.initial_stress.value)(
137 std::numeric_limits<
138 double>::quiet_NaN() /* time independent */,
139 x_position));
140 }
141
142 if (*this->process_data_.initialize_porosity_from_medium_property)
143 {
144 // Initial porosity. Could be read from integration point data
145 // or mesh.
146 current_state.porosity_data.phi =
147 medium.property(MaterialPropertyLib::porosity)
148 .template initialValue<double>(x_position,
149 time_independent);
150
151 if (medium.hasProperty(
153 {
154 current_state.transport_porosity_data.phi =
156 .template initialValue<double>(x_position,
157 time_independent);
158 }
159 else
160 {
161 current_state.transport_porosity_data.phi =
162 current_state.porosity_data.phi;
163 }
164 }
165
166 double const t = 0; // TODO (naumov) pass t from top
167 auto& material_state = this->material_states_[ip];
168 this->solid_material_.initializeInternalStateVariables(
169 t, x_position, *material_state.material_state_variables);
170
171 material_state.pushBackState();
172 }
173
174 for (unsigned ip = 0; ip < n_integration_points; ip++)
175 {
176 this->prev_states_[ip] = this->current_states_[ip];
177 }
178 }
179
180 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
181 Eigen::VectorXd const& /*local_x_prev*/,
182 double const /*t*/, double const /*dt*/,
183 int const /*process_id*/) override
184 {
185 unsigned const n_integration_points =
186 this->integration_method_.getNumberOfPoints();
187
188 for (unsigned ip = 0; ip < n_integration_points; ip++)
189 {
190 this->material_states_[ip].pushBackState();
191 }
192
193 for (unsigned ip = 0; ip < n_integration_points; ip++)
194 {
195 this->prev_states_[ip] = this->current_states_[ip];
196 }
197 }
198
200 double const t, double const dt, Eigen::VectorXd const& local_x,
201 Eigen::VectorXd const& local_x_prev) override;
202
203 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
204 const unsigned integration_point) const override
205 {
206 auto const& N_u = _secondary_data.N_u[integration_point];
207
208 // assumes N is stored contiguously in memory
209 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
210 }
211
212 std::tuple<
213 std::vector<ConstitutiveRelations::ConstitutiveData<DisplacementDim>>,
214 std::vector<
217 Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_prev,
218 double const t, double const dt,
220 models);
221
222 std::vector<ConstitutiveRelations::DerivativesData<DisplacementDim>>
224 Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_prev,
225 double const t, double const dt,
226 std::vector<
228 ip_constitutive_data,
229 std::vector<
231 ip_constitutive_variables,
233 models);
234
236 {
237 return displacement_size;
238 }
239
240private:
245 std::vector<IpData> _ip_data;
246
248 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
250
251 // The shape function of pressure has the same form with the shape function
252 // of temperature
253 static const int gas_pressure_index = 0;
254 static const int gas_pressure_size = ShapeFunctionPressure::NPOINTS;
255 static const int capillary_pressure_index = ShapeFunctionPressure::NPOINTS;
256 static const int capillary_pressure_size = ShapeFunctionPressure::NPOINTS;
257 static const int temperature_index = 2 * ShapeFunctionPressure::NPOINTS;
258 static const int temperature_size = ShapeFunctionPressure::NPOINTS;
259 static const int displacement_index = ShapeFunctionPressure::NPOINTS * 3;
260 static const int displacement_size =
261 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
262
263 static const int C_index = 0;
264 static const int C_size = ShapeFunctionPressure::NPOINTS;
265 static const int W_index = ShapeFunctionPressure::NPOINTS;
266 static const int W_size = ShapeFunctionPressure::NPOINTS;
267};
268
269} // namespace TH2M
270} // namespace ProcessLib
271
272#include "TH2MFEM-impl.h"
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
static constexpr auto & N_u_op
Definition TH2MFEM.h:70
static const int capillary_pressure_index
Definition TH2MFEM.h:255
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
Definition TH2MFEM.h:46
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
Definition TH2MFEM.h:68
typename ShapeMatricesTypePressure::template VectorType< N > VectorType
Definition TH2MFEM.h:53
std::vector< ConstitutiveRelations::DerivativesData< DisplacementDim > > updateConstitutiveVariablesDerivatives(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, std::vector< ConstitutiveRelations::ConstitutiveData< DisplacementDim > > const &ip_constitutive_data, std::vector< ConstitutiveRelations::ConstitutiveTempData< DisplacementDim > > const &ip_constitutive_variables, ConstitutiveRelations::ConstitutiveModels< DisplacementDim > const &models)
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Definition TH2MFEM.h:249
MathLib::KelvinVector::Invariants< KelvinVectorSize > Invariants
Definition TH2MFEM.h:74
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
Definition TH2MFEM.h:180
static const int capillary_pressure_size
Definition TH2MFEM.h:256
TH2MLocalAssembler(TH2MLocalAssembler &&)=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 > &local_rhs_data, std::vector< double > &local_Jac_data) override
typename ShapeMatricesTypePressure::template MatrixType< M, N > MatrixType
Definition TH2MFEM.h:57
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
Definition TH2MFEM.h:203
static const int gas_pressure_index
Definition TH2MFEM.h:253
static int const KelvinVectorSize
Definition TH2MFEM.h:66
int getNumberOfVectorElementsForDeformation() const override
Definition TH2MFEM.h:235
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
Definition TH2MFEM.h:241
IntegrationPointData< ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure > IpData
Definition TH2MFEM.h:243
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
Definition TH2MFEM.h:49
TH2MLocalAssembler(TH2MLocalAssembler const &)=delete
std::tuple< std::vector< ConstitutiveRelations::ConstitutiveData< DisplacementDim > >, std::vector< ConstitutiveRelations::ConstitutiveTempData< DisplacementDim > > > updateConstitutiveVariables(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, ConstitutiveRelations::ConstitutiveModels< DisplacementDim > const &models)
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
Definition TH2MFEM.h:63
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
Definition TH2MFEM.h:60
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
static const int displacement_index
Definition TH2MFEM.h:259
std::vector< IpData > _ip_data
Definition TH2MFEM.h:245
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
VectorType< GlobalDim > GlobalDimVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
NumLib::GenericIntegrationMethod const & integration_method_
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TH2MProcessData< DisplacementDim > &process_data)
TH2MProcessData< DisplacementDim > & process_data_
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< ConstitutiveRelations::MaterialStateData< DisplacementDim > > material_states_
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
Definition TH2MFEM.h:38