OGS
ThermoRichardsMechanicsFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
20#include "NumLib/DOF/LocalDOF.h"
28
29namespace ProcessLib
30{
31namespace ThermoRichardsMechanics
32{
33namespace MPL = MaterialPropertyLib;
34
35template <typename ShapeFunctionDisplacement, typename ShapeFunction,
36 int DisplacementDim, typename ConstitutiveTraits>
38 : public LocalAssemblerInterface<DisplacementDim, ConstitutiveTraits>
39{
40 static constexpr int temperature_index = 0;
41 static constexpr int temperature_size = ShapeFunction::NPOINTS;
42 static constexpr int pressure_index = temperature_size;
43 static constexpr int pressure_size = ShapeFunction::NPOINTS;
44 static constexpr int displacement_index = 2 * ShapeFunction::NPOINTS;
45 static constexpr int displacement_size =
46 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
47
48public:
51 // Note: temperature variable uses the same shape functions as that are used
52 // by pressure variable.
55
58
62
64 ShapeMatricesType, DisplacementDim,
65 ShapeFunctionDisplacement::NPOINTS>;
66
67 static int const KelvinVectorSize =
70
71 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
72
73 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
74 DisplacementDim,
76
81
83 MeshLib::Element const& e,
84 std::size_t const /*local_matrix_size*/,
85 NumLib::GenericIntegrationMethod const& integration_method,
86 bool const is_axially_symmetric,
88 process_data);
89
90 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
91 double const t,
92 int const process_id) override;
93
95 {
97
98 static auto constexpr local_matrix_dim =
100
101 template <Eigen::Index rows, Eigen::Index cols>
102 using Mat =
103 typename ShapeMatricesTypeDisplacement::template MatrixType<rows,
104 cols>;
105 using Vec = typename ShapeMatricesTypeDisplacement::template VectorType<
107
108 public:
109 void setZero()
110 {
111 M_TT = NodalMatrix::Zero(temperature_size, temperature_size);
112 M_Tp = NodalMatrix::Zero(temperature_size, pressure_size);
113 K_TT = NodalMatrix::Zero(temperature_size, temperature_size);
114 K_Tp = NodalMatrix::Zero(temperature_size, pressure_size);
115 dK_TT_dp = NodalMatrix::Zero(temperature_size, pressure_size);
116
119
120 M_pT = NodalMatrix::Zero(pressure_size, temperature_size);
121
122 K_pp = NodalMatrix::Zero(pressure_size, pressure_size);
123 K_pT = NodalMatrix::Zero(pressure_size, temperature_size);
124
125 storage_p_a_p = NodalMatrix::Zero(pressure_size, pressure_size);
126 storage_p_a_S_Jpp = NodalMatrix::Zero(pressure_size, pressure_size);
127 storage_p_a_S = NodalMatrix::Zero(pressure_size, pressure_size);
128
131 res = Vec::Zero(local_matrix_dim);
132 }
133
135 {
136 M_TT += other.M_TT;
137 M_Tp += other.M_Tp;
138 K_TT += other.K_TT;
139 K_Tp += other.K_Tp;
140 dK_TT_dp += other.dK_TT_dp;
141
142 M_pu += other.M_pu;
143
144 M_pT += other.M_pT;
145
146 K_pp += other.K_pp;
147 K_pT += other.K_pT;
148
152
153 Jac += other.Jac;
154 res += other.res;
155
156 return *this;
157 }
158
159 LocalMatrices& operator*=(double const a)
160 {
161 M_TT *= a;
162 M_Tp *= a;
163 K_TT *= a;
164 K_Tp *= a;
165 dK_TT_dp *= a;
166
167 M_pu *= a;
168
169 M_pT *= a;
170
171 K_pp *= a;
172 K_pT *= a;
173
174 storage_p_a_p *= a;
176 storage_p_a_S *= a;
177
178 Jac *= a;
179 res *= a;
180
181 return *this;
182 }
183
184 template <typename OStream>
185 friend OStream& operator<<(OStream& os, LocalMatrices const& loc_mat)
186 {
187 os << "- M_TT:\n"
188 << loc_mat.M_TT //
189 << "\n- M_Tp:\n"
190 << loc_mat.M_Tp //
191 << "\n- K_TT:\n"
192 << loc_mat.K_TT //
193 << "\n- K_Tp:\n"
194 << loc_mat.K_Tp //
195 << "\n- dK_TT_dp:\n"
196 << loc_mat.dK_TT_dp //
197 << "\n- M_pu:\n"
198 << loc_mat.M_pu //
199 << "\n- M_pT:\n"
200 << loc_mat.M_pT //
201 << "\n- K_pp:\n"
202 << loc_mat.K_pp //
203 << "\n- K_pT:\n"
204 << loc_mat.K_pT //
205 << "\n- storage_p_a_p:\n"
206 << loc_mat.storage_p_a_p //
207 << "\n- storage_p_a_S_Jpp:\n"
208 << loc_mat.storage_p_a_S_Jpp //
209 << "\n- storage_p_a_S:\n"
210 << loc_mat.storage_p_a_S //
211 << "\n- Jac:\n"
212 << loc_mat.Jac //
213 << "\n- res:\n"
214 << loc_mat.res;
215 return os;
216 }
217
223
225
227
230
234
238
242 };
243
244 void assembleWithJacobian(double const t, double const dt,
245 std::vector<double> const& local_x,
246 std::vector<double> const& local_x_prev,
247 std::vector<double>& local_rhs_data,
248 std::vector<double>& local_Jac_data) override;
249
250private:
252 double const t, double const dt,
253 ParameterLib::SpatialPosition const& x_position,
254 std::vector<double> const& local_x,
255 std::vector<double> const& local_x_prev, IpData const& ip_data,
256 typename ConstitutiveTraits::ConstitutiveSetting& CS,
258 typename ConstitutiveTraits::StatefulData& current_state,
259 typename ConstitutiveTraits::StatefulDataPrev const& prev_state,
261 typename ConstitutiveTraits::OutputData& output_data) const;
262
263 void addToLocalMatrixData(double const dt,
264 std::vector<double> const& local_x,
265 std::vector<double> const& local_x_prev,
266 LocalMatrices const& loc_mat,
267 std::vector<double>& local_rhs_data,
268 std::vector<double>& local_Jac_data) const;
269
270 void massLumping(LocalMatrices& loc_mat) const;
271
272public:
273 void initializeConcrete() override
274 {
275 unsigned const n_integration_points =
277 auto const time_independent = std::numeric_limits<double>::quiet_NaN();
278 auto const& medium =
279 *this->process_data_.media_map.getMedium(this->element_.getID());
280
281 for (unsigned ip = 0; ip < n_integration_points; ip++)
282 {
283 ParameterLib::SpatialPosition const x_position{
284 std::nullopt, this->element_.getID(), ip,
286 ShapeFunctionDisplacement,
288 this->element_, this->ip_data_[ip].N_u))};
289 auto& current_state = this->current_states_[ip];
290
291 // Set initial stress from parameter.
292 if (this->process_data_.initial_stress.value)
293 {
294 ConstitutiveTraits::ConstitutiveSetting::statefulStress(
295 current_state) =
297 DisplacementDim>(
298 (*this->process_data_.initial_stress.value)(
299 time_independent, x_position));
300 }
301
302 if (*this->process_data_.initialize_porosity_from_medium_property)
303 {
304 // Initial porosity. Could be read from integration point data
305 // or mesh.
306 std::get<PorosityData>(current_state).phi =
307 medium.property(MPL::porosity)
308 .template initialValue<double>(x_position,
309 time_independent);
310
311 if (medium.hasProperty(MPL::PropertyType::transport_porosity))
312 {
313 std::get<TransportPorosityData>(current_state).phi =
314 medium.property(MPL::transport_porosity)
315 .template initialValue<double>(x_position,
316 time_independent);
317 }
318 else
319 {
320 std::get<TransportPorosityData>(current_state).phi =
321 std::get<PorosityData>(current_state).phi;
322 }
323 }
324
325 double const t = 0; // TODO (naumov) pass t from top
326 this->solid_material_.initializeInternalStateVariables(
327 t, x_position,
328 *this->material_states_[ip].material_state_variables);
329 }
330
331 for (unsigned ip = 0; ip < n_integration_points; ip++)
332 {
333 this->material_states_[ip].pushBackState();
334 }
335
336 for (unsigned ip = 0; ip < n_integration_points; ip++)
337 {
338 this->prev_states_[ip] = this->current_states_[ip];
339 }
340 }
341
343 double const t, double const dt, Eigen::VectorXd const& local_x,
344 Eigen::VectorXd const& local_x_prev) override;
345
346 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
347 const unsigned integration_point) const override
348 {
349 auto const& N_u = ip_data_[integration_point].N_u;
350
351 // assumes N is stored contiguously in memory
352 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
353 }
354
355private:
356 std::vector<IpData> ip_data_;
357
358 static constexpr auto localDOF(std::vector<double> const& x)
359 {
360 return NumLib::localDOF<
361 ShapeFunction, ShapeFunction,
363 }
364
365 static auto block_uu(auto& mat)
366 {
367 return mat.template block<displacement_size, displacement_size>(
369 }
370 static auto block_up(auto& mat)
371 {
372 return mat.template block<displacement_size, pressure_size>(
374 }
375 static auto block_uT(auto& mat)
376 {
377 return mat.template block<displacement_size, temperature_size>(
379 }
380 static auto block_pu(auto& mat)
381 {
382 return mat.template block<pressure_size, displacement_size>(
384 }
385 static auto block_pp(auto& mat)
386 {
387 return mat.template block<pressure_size, pressure_size>(pressure_index,
389 }
390 static auto block_pT(auto& mat)
391 {
392 return mat.template block<pressure_size, temperature_size>(
394 }
395 static auto block_Tp(auto& mat)
396 {
397 return mat.template block<temperature_size, pressure_size>(
399 }
400 static auto block_TT(auto& mat)
401 {
402 return mat.template block<temperature_size, temperature_size>(
404 }
405
406 static auto block_u(auto& vec)
407 {
408 return vec.template segment<displacement_size>(displacement_index);
409 }
410 static auto block_p(auto& vec)
411 {
412 return vec.template segment<pressure_size>(pressure_index);
413 }
414 static auto block_T(auto& vec)
415 {
416 return vec.template segment<temperature_size>(temperature_index);
417 }
418
426 unsigned const ip, double const t,
427 ParameterLib::SpatialPosition const x_position,
428 MaterialPropertyLib::Medium const& medium,
429 MPL::VariableArray const& variables, double const p_at_ip);
430};
431
432} // namespace ThermoRichardsMechanics
433} // namespace ProcessLib
434
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
VectorType< _kelvin_vector_size > KelvinVectorType
typename ShapeMatricesTypeDisplacement::template MatrixType< rows, cols > Mat
typename ShapeMatricesTypeDisplacement::template VectorType< local_matrix_dim > Vec
void addToLocalMatrixData(double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, LocalMatrices const &loc_mat, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) const
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
void assembleWithJacobianSingleIP(double const t, double const dt, ParameterLib::SpatialPosition const &x_position, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, IpData const &ip_data, typename ConstitutiveTraits::ConstitutiveSetting &CS, MaterialPropertyLib::Medium &medium, LocalMatrices &out, typename ConstitutiveTraits::StatefulData &current_state, typename ConstitutiveTraits::StatefulDataPrev const &prev_state, MaterialStateData< DisplacementDim > &mat_state, typename ConstitutiveTraits::OutputData &output_data) const
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
ThermoRichardsMechanicsLocalAssembler(ThermoRichardsMechanicsLocalAssembler &&)=delete
void convertInitialStressType(unsigned const ip, double const t, ParameterLib::SpatialPosition const x_position, MaterialPropertyLib::Medium const &medium, MPL::VariableArray const &variables, double const p_at_ip)
ThermoRichardsMechanicsLocalAssembler(ThermoRichardsMechanicsLocalAssembler const &)=delete
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
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
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:64
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data_
std::vector< typename ConstitutiveTraits::StatefulDataPrev > prev_states_
std::vector< MaterialStateData< DisplacementDim > > material_states_
ConstitutiveTraits::SolidConstitutiveRelation const & solid_material_
std::vector< typename ConstitutiveTraits::StatefulData > current_states_