Loading [MathJax]/extensions/tex2jax.js
OGS
TH2MFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
35
36namespace ProcessLib
37{
38namespace TH2M
39{
42template <typename ShapeMatrixType>
44{
45 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
46};
47
48template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
49 int DisplacementDim>
50class TH2MLocalAssembler : public LocalAssemblerInterface<DisplacementDim>
51{
52public:
55
58
59 template <int N>
60 using VectorType =
61 typename ShapeMatricesTypePressure::template VectorType<N>;
62
63 template <int M, int N>
64 using MatrixType =
65 typename ShapeMatricesTypePressure::template MatrixType<M, N>;
66
69
72
73 static int const KelvinVectorSize =
75 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
76
77 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
78 DisplacementDim,
80
82
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
93private:
96 std::string_view const name,
97 double const* values,
98 int const integration_order) override;
99
100 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
101 double const t,
102 int const process_id) override;
103
104 void assemble(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_M_data*/,
108 std::vector<double>& /*local_K_data*/,
109 std::vector<double>& /*local_rhs_data*/) override;
110
111 void assembleWithJacobian(double const t, double const dt,
112 std::vector<double> const& local_x,
113 std::vector<double> const& local_x_prev,
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 auto const time_independent = std::numeric_limits<double>::quiet_NaN();
122 auto const& medium =
123 *this->process_data_.media_map.getMedium(this->element_.getID());
124
125 for (unsigned ip = 0; ip < n_integration_points; ip++)
126 {
127 auto& ip_data = _ip_data[ip];
128
129 ParameterLib::SpatialPosition const x_position{
130 std::nullopt, this->element_.getID(),
132 ShapeFunctionDisplacement,
134 ip_data.N_u))};
135 auto& current_state = this->current_states_[ip];
136
138 if (this->process_data_.initial_stress.value)
139 {
140 current_state.eff_stress_data.sigma_eff.noalias() =
142 DisplacementDim>(
143 (*this->process_data_.initial_stress.value)(
144 std::numeric_limits<
145 double>::quiet_NaN() /* time independent */,
146 x_position));
147 }
148
149 if (*this->process_data_.initialize_porosity_from_medium_property)
150 {
151 // Initial porosity. Could be read from integration point data
152 // or mesh.
153 current_state.porosity_data.phi =
154 medium.property(MaterialPropertyLib::porosity)
155 .template initialValue<double>(x_position,
156 time_independent);
157
158 if (medium.hasProperty(
160 {
161 current_state.transport_porosity_data.phi =
163 .template initialValue<double>(x_position,
164 time_independent);
165 }
166 else
167 {
168 current_state.transport_porosity_data.phi =
169 current_state.porosity_data.phi;
170 }
171 }
172
173 double const t = 0; // TODO (naumov) pass t from top
174 auto& material_state = this->material_states_[ip];
175 this->solid_material_.initializeInternalStateVariables(
176 t, x_position, *material_state.material_state_variables);
177
178 material_state.pushBackState();
179 }
180
181 for (unsigned ip = 0; ip < n_integration_points; ip++)
182 {
183 this->prev_states_[ip] = this->current_states_[ip];
184 }
185 }
186
187 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
188 Eigen::VectorXd const& /*local_x_prev*/,
189 double const /*t*/, double const /*dt*/,
190 int const /*process_id*/) override
191 {
192 unsigned const n_integration_points =
194
195 for (unsigned ip = 0; ip < n_integration_points; ip++)
196 {
197 this->material_states_[ip].pushBackState();
198 }
199
200 for (unsigned ip = 0; ip < n_integration_points; ip++)
201 {
202 this->prev_states_[ip] = this->current_states_[ip];
203 }
204 }
205
207 double const t, double const dt, Eigen::VectorXd const& local_x,
208 Eigen::VectorXd const& local_x_prev) override;
209
210 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
211 const unsigned integration_point) const override
212 {
213 auto const& N_u = _secondary_data.N_u[integration_point];
214
215 // assumes N is stored contiguously in memory
216 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
217 }
218
219 std::tuple<
220 std::vector<ConstitutiveRelations::ConstitutiveData<DisplacementDim>>,
221 std::vector<
224 Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_prev,
225 double const t, double const dt,
227 models);
228
229 std::vector<ConstitutiveRelations::DerivativesData<DisplacementDim>>
231 Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_prev,
232 double const t, double const dt,
233 std::vector<
235 ip_constitutive_data,
236 std::vector<
238 ip_constitutive_variables,
240 models);
241
242 virtual std::optional<VectorSegment> getVectorDeformationSegment()
243 const override
244 {
245 return std::optional<VectorSegment>{
247 }
248
249private:
254 std::vector<IpData> _ip_data;
255
257 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
259
260 // The shape function of pressure has the same form with the shape function
261 // of temperature
262 static const int gas_pressure_index = 0;
263 static const int gas_pressure_size = ShapeFunctionPressure::NPOINTS;
264 static const int capillary_pressure_index = ShapeFunctionPressure::NPOINTS;
265 static const int capillary_pressure_size = ShapeFunctionPressure::NPOINTS;
266 static const int temperature_index = 2 * ShapeFunctionPressure::NPOINTS;
267 static const int temperature_size = ShapeFunctionPressure::NPOINTS;
268 static const int displacement_index = ShapeFunctionPressure::NPOINTS * 3;
269 static const int displacement_size =
270 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
271
272 static const int C_index = 0;
273 static const int C_size = ShapeFunctionPressure::NPOINTS;
274 static const int W_index = ShapeFunctionPressure::NPOINTS;
275 static const int W_size = ShapeFunctionPressure::NPOINTS;
276};
277
278} // namespace TH2M
279} // namespace ProcessLib
280
281#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:77
static const int capillary_pressure_index
Definition TH2MFEM.h:264
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
Definition TH2MFEM.h:53
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
Definition TH2MFEM.h:75
typename ShapeMatricesTypePressure::template VectorType< N > VectorType
Definition TH2MFEM.h:60
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)
virtual std::optional< VectorSegment > getVectorDeformationSegment() const override
Definition TH2MFEM.h:242
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Definition TH2MFEM.h:258
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
Definition TH2MFEM.h:187
static const int capillary_pressure_size
Definition TH2MFEM.h:265
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:64
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
Definition TH2MFEM.h:210
static const int gas_pressure_index
Definition TH2MFEM.h:262
static int const KelvinVectorSize
Definition TH2MFEM.h:73
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:70
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
Definition TH2MFEM.h:67
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
static const int displacement_index
Definition TH2MFEM.h:268
std::vector< IpData > _ip_data
Definition TH2MFEM.h:254
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
Data that is needed for the equation system assembly.
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:45