OGS
ThermoMechanicsFEM.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
22
23namespace ProcessLib
24{
25namespace ThermoMechanics
26{
27namespace MPL = MaterialPropertyLib;
28
29template <typename BMatricesType, typename ShapeMatricesType,
30 int DisplacementDim>
32{
41
42 // total stress
43 typename BMatricesType::KelvinVectorType sigma, sigma_prev;
44
45 // total strain
46 typename BMatricesType::KelvinVectorType eps, eps_prev;
47
48 // mechanical strain
49 typename BMatricesType::KelvinVectorType eps_m, eps_m_prev;
50
52 std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
53 DisplacementDim>::MaterialStateVariables>
55
57 typename ShapeMatricesType::NodalRowVectorType N;
58 typename ShapeMatricesType::GlobalDimNodalMatrixType dNdx;
59
61 {
62 eps_prev = eps;
65 material_state_variables->pushBackState();
66 }
67
69};
70
73template <typename ShapeMatrixType>
75{
76 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
77};
78
118
119template <typename ShapeFunction, int DisplacementDim>
121 : public ThermoMechanicsLocalAssemblerInterface<DisplacementDim>
122{
123public:
126
127 // Types for displacement.
128 // (Higher order elements = ShapeFunction).
132
133 static int const KelvinVectorSize =
136
138 using RhsVector = typename ShapeMatricesType::template VectorType<
139 ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim>;
140 using JacobianMatrix = typename ShapeMatricesType::template MatrixType<
141 ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim,
142 ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim>;
143 using IpData =
145
147 delete;
149
151 MeshLib::Element const& e,
152 std::size_t const /*local_matrix_size*/,
153 NumLib::GenericIntegrationMethod const& integration_method,
154 bool const is_axially_symmetric,
156
158 std::size_t setIPDataInitialConditions(
159 std::string_view const name,
160 double const* values,
161 int const integration_order) override;
162
163 void assemble(double const /*t*/, double const /*dt*/,
164 std::vector<double> const& /*local_x*/,
165 std::vector<double> const& /*local_x_prev*/,
166 std::vector<double>& /*local_M_data*/,
167 std::vector<double>& /*local_K_data*/,
168 std::vector<double>& /*local_rhs_data*/) override
169 {
170 OGS_FATAL(
171 "ThermoMechanicsLocalAssembler: assembly without jacobian is not "
172 "implemented.");
173 }
174
176 double const t, double const dt, Eigen::VectorXd const& local_x,
177 Eigen::VectorXd const& local_x_prev, int const process_id,
178 std::vector<double>& local_b_data,
179 std::vector<double>& local_Jac_data) override;
180
181 void assembleWithJacobian(double const t, double const dt,
182 std::vector<double> const& local_x,
183 std::vector<double> const& local_x_prev,
184 std::vector<double>& local_rhs_data,
185 std::vector<double>& local_Jac_data) override;
186
187 void initializeConcrete() override
188 {
189 unsigned const n_integration_points =
190 _integration_method.getNumberOfPoints();
191 for (unsigned ip = 0; ip < n_integration_points; ip++)
192 {
193 auto& ip_data = _ip_data[ip];
194
195 ParameterLib::SpatialPosition const x_position{
196 std::nullopt, _element.getID(),
200 _element, ip_data.N))};
201
203 if (_process_data.initial_stress != nullptr)
204 {
205 ip_data.sigma =
207 DisplacementDim>((*_process_data.initial_stress)(
208 std::numeric_limits<
209 double>::quiet_NaN() /* time independent */,
210 x_position));
211 }
212
213 double const t = 0; // TODO (naumov) pass t from top
214 ip_data.solid_material.initializeInternalStateVariables(
215 t, x_position, *ip_data.material_state_variables);
216
217 ip_data.pushBackState();
218 }
219 }
220
221 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
222 double const t,
223 int const process_id) override;
224
225 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
226 Eigen::VectorXd const& /*local_x_prev*/,
227 double const /*t*/, double const /*dt*/,
228 int const /*process_id*/) override
229 {
230 unsigned const n_integration_points =
231 _integration_method.getNumberOfPoints();
232
233 for (unsigned ip = 0; ip < n_integration_points; ip++)
234 {
235 _ip_data[ip].pushBackState();
236 }
237 }
238
239 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
240 const unsigned integration_point) const override
241 {
242 auto const& N = _secondary_data.N[integration_point];
243
244 // assumes N is stored contiguously in memory
245 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
246 }
247
249 {
250 return displacement_size;
251 }
252
253private:
275
277 const double t, double const dt, Eigen::VectorXd const& local_x,
278 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
279 std::vector<double>& local_Jac_data);
280
302 const double t, double const dt, Eigen::VectorXd const& local_x,
303 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_b_data,
304 std::vector<double>& local_Jac_data);
305
306 std::size_t setSigma(double const* values);
307
308 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
309 // the ordering of the cache_mat.
310 // There should be only one.
311 std::vector<double> getSigma() const override;
312
313 std::vector<double> const& getIntPtSigma(
314 const double t,
315 std::vector<GlobalVector*> const& x,
316 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
317 std::vector<double>& cache) const override;
318
319 std::size_t setEpsilon(double const* values);
320
321 std::vector<double> getEpsilon() const override;
322
323 std::vector<double> const& getIntPtEpsilon(
324 const double t,
325 std::vector<GlobalVector*> const& x,
326 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
327 std::vector<double>& cache) const override;
328
329 std::vector<double> const& getIntPtEpsilonMechanical(
330 const double t,
331 std::vector<GlobalVector*> const& x,
332 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
333 std::vector<double>& cache) const override;
334
335 std::size_t setEpsilonMechanical(double const* values);
336
337 std::vector<double> getEpsilonMechanical() const override;
338
339 unsigned getNumberOfIntegrationPoints() const override;
340
341 int getMaterialID() const override;
342
344 DisplacementDim>::MaterialStateVariables const&
345 getMaterialStateVariablesAt(unsigned integration_point) const override;
346
347private:
349
350 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
351
356
357 static const int temperature_index = 0;
358 static const int temperature_size = ShapeFunction::NPOINTS;
359 static const int displacement_index = ShapeFunction::NPOINTS;
360 static const int displacement_size =
361 ShapeFunction::NPOINTS * DisplacementDim;
362};
363
364} // namespace ThermoMechanics
365} // namespace ProcessLib
366
#define OGS_FATAL(...)
Definition Error.h:19
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
std::vector< double > const & getIntPtEpsilonMechanical(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
typename BMatricesType::NodalForceVectorType NodalForceVectorType
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
Returns number of read integration points.
MathLib::KelvinVector::Invariants< KelvinVectorSize > Invariants
typename ShapeMatricesType::template MatrixType< ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim, ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim > JacobianMatrix
BMatrixPolicyType< ShapeFunction, DisplacementDim > BMatricesType
std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
ThermoMechanicsLocalAssembler(ThermoMechanicsLocalAssembler const &)=delete
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
typename ShapeMatricesType::template VectorType< ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim > RhsVector
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim > IpData
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
ThermoMechanicsProcessData< DisplacementDim > & _process_data
NumLib::GenericIntegrationMethod const & _integration_method
void assembleWithJacobianForHeatConductionEquations(const double 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 assembleWithJacobianForDeformationEquations(const double 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 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
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
ThermoMechanicsLocalAssembler(ThermoMechanicsLocalAssembler &&)=delete
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
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)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material
ShapeMatricesType::NodalRowVectorType N
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
ShapeMatricesType::GlobalDimNodalMatrixType dNdx
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N