OGS
ThermoMechanicsFEM.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <memory>
14 #include <vector>
15 
18 #include "MathLib/KelvinVector.h"
24 #include "ParameterLib/Parameter.h"
28 
29 namespace ProcessLib
30 {
31 namespace ThermoMechanics
32 {
33 namespace MPL = MaterialPropertyLib;
34 
35 template <typename BMatricesType, typename ShapeMatricesType,
36  int DisplacementDim>
38 {
44  solid_material.createMaterialStateVariables())
45  {
46  }
47 
48  // total stress
50 
51  // total strain
53 
54  // mechanical strain
56 
58  std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
59  DisplacementDim>::MaterialStateVariables>
61 
63  typename ShapeMatricesType::NodalRowVectorType N;
64  typename ShapeMatricesType::GlobalDimNodalMatrixType dNdx;
65 
67  {
68  eps_prev = eps;
69  eps_m_prev = eps_m;
70  sigma_prev = sigma;
71  material_state_variables->pushBackState();
72  }
73 
75 };
76 
79 template <typename ShapeMatrixType>
81 {
82  std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
83 };
84 
125 template <typename ShapeFunction, typename IntegrationMethod,
126  int DisplacementDim>
128  : public ThermoMechanicsLocalAssemblerInterface<DisplacementDim>
129 {
130 public:
133 
134  // Types for displacement.
135  // (Higher order elements = ShapeFunction).
139 
140  static int const KelvinVectorSize =
143 
145  using RhsVector = typename ShapeMatricesType::template VectorType<
146  ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim>;
147  using JacobianMatrix = typename ShapeMatricesType::template MatrixType<
148  ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim,
149  ShapeFunction::NPOINTS + ShapeFunction::NPOINTS * DisplacementDim>;
150  using IpData =
152 
154  delete;
156 
158  MeshLib::Element const& e,
159  std::size_t const /*local_matrix_size*/,
160  bool const is_axially_symmetric,
161  unsigned const integration_order,
163 
165  std::size_t setIPDataInitialConditions(
166  std::string const& name,
167  double const* values,
168  int const integration_order) override;
169 
170  void assemble(double const /*t*/, double const /*dt*/,
171  std::vector<double> const& /*local_x*/,
172  std::vector<double> const& /*local_xdot*/,
173  std::vector<double>& /*local_M_data*/,
174  std::vector<double>& /*local_K_data*/,
175  std::vector<double>& /*local_rhs_data*/) override
176  {
177  OGS_FATAL(
178  "ThermoMechanicsLocalAssembler: assembly without jacobian is not "
179  "implemented.");
180  }
181 
183  double const t, double const dt, Eigen::VectorXd const& local_x,
184  Eigen::VectorXd const& local_xdot, const double dxdot_dx,
185  const double dx_dx, int const process_id,
186  std::vector<double>& local_M_data, std::vector<double>& local_K_data,
187  std::vector<double>& local_b_data,
188  std::vector<double>& local_Jac_data) override;
189 
190  void assembleWithJacobian(double const t, double const dt,
191  std::vector<double> const& local_x,
192  std::vector<double> const& local_xdot,
193  const double dxdot_dx, const double dx_dx,
194  std::vector<double>& local_M_data,
195  std::vector<double>& local_K_data,
196  std::vector<double>& local_rhs_data,
197  std::vector<double>& local_Jac_data) override;
198 
199  void initializeConcrete() override
200  {
201  unsigned const n_integration_points =
202  _integration_method.getNumberOfPoints();
203  for (unsigned ip = 0; ip < n_integration_points; ip++)
204  {
205  auto& ip_data = _ip_data[ip];
206 
208  if (_process_data.initial_stress != nullptr)
209  {
210  ParameterLib::SpatialPosition const x_position{
211  std::nullopt, _element.getID(), ip,
213  NumLib::interpolateCoordinates<ShapeFunction,
215  _element, ip_data.N))};
216 
217  ip_data.sigma =
219  DisplacementDim>((*_process_data.initial_stress)(
220  std::numeric_limits<
221  double>::quiet_NaN() /* time independent */,
222  x_position));
223  }
224 
225  ip_data.pushBackState();
226  }
227  }
228 
229  void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
230  double const /*t*/, double const /*dt*/) override
231  {
232  unsigned const n_integration_points =
233  _integration_method.getNumberOfPoints();
234 
235  for (unsigned ip = 0; ip < n_integration_points; ip++)
236  {
237  _ip_data[ip].pushBackState();
238  }
239  }
240 
241  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
242  const unsigned integration_point) const override
243  {
244  auto const& N = _secondary_data.N[integration_point];
245 
246  // assumes N is stored contiguously in memory
247  return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
248  }
249 
250 private:
274  const double t, double const dt, Eigen::VectorXd const& local_x,
275  Eigen::VectorXd const& local_xdot, std::vector<double>& local_b_data,
276  std::vector<double>& local_Jac_data);
277 
299  const double t, double const dt, Eigen::VectorXd const& local_x,
300  Eigen::VectorXd const& local_xdot, std::vector<double>& local_b_data,
301  std::vector<double>& local_Jac_data);
302 
303  std::size_t setSigma(double const* values);
304 
305  // TODO (naumov) This method is same as getIntPtSigma but for arguments and
306  // the ordering of the cache_mat.
307  // There should be only one.
308  std::vector<double> getSigma() const override;
309 
310  std::vector<double> const& getIntPtSigma(
311  const double t,
312  std::vector<GlobalVector*> const& x,
313  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
314  std::vector<double>& cache) const override;
315 
316  std::size_t setEpsilon(double const* values);
317 
318  std::vector<double> getEpsilon() const override;
319 
320  std::vector<double> const& getIntPtEpsilon(
321  const double t,
322  std::vector<GlobalVector*> const& x,
323  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
324  std::vector<double>& cache) const override;
325 
326  std::vector<double> const& getIntPtEpsilonMechanical(
327  const double t,
328  std::vector<GlobalVector*> const& x,
329  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
330  std::vector<double>& cache) const override;
331 
332  std::size_t setEpsilonMechanical(double const* values);
333 
334  std::vector<double> getEpsilonMechanical() const override;
335 
336  unsigned getNumberOfIntegrationPoints() const override;
337 
339  DisplacementDim>::MaterialStateVariables const&
340  getMaterialStateVariablesAt(unsigned integration_point) const override;
341 
342 private:
343 
345 
346  std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
347 
348  IntegrationMethod _integration_method;
352 
353  static const int temperature_index = 0;
354  static const int temperature_size = ShapeFunction::NPOINTS;
355  static const int displacement_index = ShapeFunction::NPOINTS;
356  static const int displacement_size =
357  ShapeFunction::NPOINTS * DisplacementDim;
358 };
359 
360 } // namespace ThermoMechanics
361 } // namespace ProcessLib
362 
363 #include "ThermoMechanicsFEM-impl.h"
#define OGS_FATAL(...)
Definition: Error.h:26
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
Definition: BMatrixPolicy.h:44
Local assembler of ThermoMechanics process.
typename BMatricesType::NodalForceVectorType NodalForceVectorType
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
typename ShapeMatricesType::template MatrixType< ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim, ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim > JacobianMatrix
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 assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) 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
ThermoMechanicsLocalAssembler(ThermoMechanicsLocalAssembler const &)=delete
void assembleWithJacobianForHeatConductionEquations(const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void postTimestepConcrete(Eigen::VectorXd const &, double const, double const) override
ThermoMechanicsLocalAssembler(ThermoMechanicsLocalAssembler &&)=delete
std::size_t setIPDataInitialConditions(std::string const &name, double const *values, int const integration_order) override
Returns number of read integration points.
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void assembleWithJacobianForDeformationEquations(const double t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
void assembleWithJacobianForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
ThermoMechanicsProcessData< DisplacementDim > & _process_data
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
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
typename ShapeMatricesType::template VectorType< ShapeFunction::NPOINTS+ShapeFunction::NPOINTS *DisplacementDim > RhsVector
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Definition: KelvinVector.h:170
MathLib::TemplatePoint< double, 3 > Point3d
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
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
Definition: SecondaryData.h:28
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