OGS
TH2MFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
33
34namespace ProcessLib
35{
36namespace TH2M
37{
40template <typename ShapeMatrixType>
42{
43 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
44};
45
46template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
47 int DisplacementDim>
48class TH2MLocalAssembler : public LocalAssemblerInterface<DisplacementDim>
49{
50public:
53
56
57 template <int N>
58 using VectorType =
59 typename ShapeMatricesTypePressure::template VectorType<N>;
60
61 template <int M, int N>
62 using MatrixType =
63 typename ShapeMatricesTypePressure::template MatrixType<M, N>;
64
67
70
71 static int const KelvinVectorSize =
73 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
74
75 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
76 DisplacementDim,
78
80
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
91private:
94 std::string_view const name,
95 double const* values,
96 int const integration_order) override;
97
98 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
99 double const t,
100 int const process_id) override;
101
102 void assemble(double const /*t*/, double const /*dt*/,
103 std::vector<double> const& /*local_x*/,
104 std::vector<double> const& /*local_x_prev*/,
105 std::vector<double>& /*local_M_data*/,
106 std::vector<double>& /*local_K_data*/,
107 std::vector<double>& /*local_rhs_data*/) override;
108
109 void assembleWithJacobian(double const t, double const dt,
110 std::vector<double> const& local_x,
111 std::vector<double> const& local_x_prev,
112 std::vector<double>& /*local_M_data*/,
113 std::vector<double>& /*local_K_data*/,
114 std::vector<double>& local_rhs_data,
115 std::vector<double>& local_Jac_data) override;
116
117 void initializeConcrete() override
118 {
119 unsigned const n_integration_points =
121
122 for (unsigned ip = 0; ip < n_integration_points; ip++)
123 {
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 this->current_states_[ip].eff_stress_data.sigma.noalias() =
138 DisplacementDim>(
139 (*this->process_data_.initial_stress.value)(
140 std::numeric_limits<
141 double>::quiet_NaN() /* time independent */,
142 x_position));
143 }
144
145 double const t = 0; // TODO (naumov) pass t from top
146 auto& material_state = this->material_states_[ip];
147 this->solid_material_.initializeInternalStateVariables(
148 t, x_position, *material_state.material_state_variables);
149
150 material_state.pushBackState();
151 }
152
153 for (unsigned ip = 0; ip < n_integration_points; ip++)
154 {
155 this->prev_states_[ip] = this->current_states_[ip];
156 }
157 }
158
159 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
160 Eigen::VectorXd const& /*local_x_prev*/,
161 double const /*t*/, double const /*dt*/,
162 int const /*process_id*/) override
163 {
164 unsigned const n_integration_points =
166
167 for (unsigned ip = 0; ip < n_integration_points; ip++)
168 {
169 this->material_states_[ip].pushBackState();
170 }
171
172 for (unsigned ip = 0; ip < n_integration_points; ip++)
173 {
174 this->prev_states_[ip] = this->current_states_[ip];
175 }
176 }
177
179 double const t, double const dt, Eigen::VectorXd const& local_x,
180 Eigen::VectorXd const& local_x_prev) override;
181
182 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
183 const unsigned integration_point) const override
184 {
185 auto const& N_u = _secondary_data.N_u[integration_point];
186
187 // assumes N is stored contiguously in memory
188 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
189 }
190
191 std::tuple<
192 std::vector<ConstitutiveRelations::ConstitutiveData<DisplacementDim>>,
193 std::vector<
196 Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_prev,
197 double const t, double const dt,
199 models);
200
201private:
206 std::vector<IpData> _ip_data;
207
209 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
211
212 // The shape function of pressure has the same form with the shape function
213 // of temperature
214 static const int gas_pressure_index = 0;
215 static const int gas_pressure_size = ShapeFunctionPressure::NPOINTS;
216 static const int capillary_pressure_index = ShapeFunctionPressure::NPOINTS;
217 static const int capillary_pressure_size = ShapeFunctionPressure::NPOINTS;
218 static const int temperature_index = 2 * ShapeFunctionPressure::NPOINTS;
219 static const int temperature_size = ShapeFunctionPressure::NPOINTS;
220 static const int displacement_index = ShapeFunctionPressure::NPOINTS * 3;
221 static const int displacement_size =
222 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
223
224 static const int C_index = 0;
225 static const int C_size = ShapeFunctionPressure::NPOINTS;
226 static const int W_index = ShapeFunctionPressure::NPOINTS;
227 static const int W_size = ShapeFunctionPressure::NPOINTS;
228};
229
230} // namespace TH2M
231} // namespace ProcessLib
232
233#include "TH2MFEM-impl.h"
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
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:75
static const int capillary_pressure_index
Definition TH2MFEM.h:216
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
Definition TH2MFEM.h:51
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
Definition TH2MFEM.h:73
typename ShapeMatricesTypePressure::template VectorType< N > VectorType
Definition TH2MFEM.h:58
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Definition TH2MFEM.h:210
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
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
Definition TH2MFEM.h:159
static const int capillary_pressure_size
Definition TH2MFEM.h:217
TH2MLocalAssembler(TH2MLocalAssembler &&)=delete
typename ShapeMatricesTypePressure::template MatrixType< M, N > MatrixType
Definition TH2MFEM.h:62
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
Definition TH2MFEM.h:182
static const int gas_pressure_index
Definition TH2MFEM.h:214
static int const KelvinVectorSize
Definition TH2MFEM.h:71
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:68
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
Definition TH2MFEM.h:65
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
static const int displacement_index
Definition TH2MFEM.h:220
std::vector< IpData > _ip_data
Definition TH2MFEM.h:206
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
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
NumLib::GenericIntegrationMethod const & integration_method_
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:43