OGS
ThermoHydroMechanicsFEM.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <memory>
14 #include <vector>
15 
16 #include "IntegrationPointData.h"
20 #include "MathLib/KelvinVector.h"
25 #include "ParameterLib/Parameter.h"
32 
33 namespace ProcessLib
34 {
35 namespace ThermoHydroMechanics
36 {
37 namespace MPL = MaterialPropertyLib;
38 
41 template <typename ShapeMatrixType>
43 {
44  std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
45 };
46 
47 template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
48  typename IntegrationMethod, int DisplacementDim>
50 {
51 public:
54 
55  // Types for pressure.
58 
61 
64 
65  static int const KelvinVectorSize =
68 
69  using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
70 
72  ThermoHydroMechanicsLocalAssembler const&) = delete;
74  delete;
75 
77  MeshLib::Element const& e,
78  std::size_t const /*local_matrix_size*/,
79  bool const is_axially_symmetric,
80  unsigned const integration_order,
82 
84  std::size_t setIPDataInitialConditions(
85  std::string const& name,
86  double const* values,
87  int const integration_order) override;
88 
89  void assemble(double const /*t*/, double const /*dt*/,
90  std::vector<double> const& /*local_x*/,
91  std::vector<double> const& /*local_xdot*/,
92  std::vector<double>& /*local_M_data*/,
93  std::vector<double>& /*local_K_data*/,
94  std::vector<double>& /*local_rhs_data*/) override
95  {
96  OGS_FATAL(
97  "ThermoHydroMechanicsLocalAssembler: assembly without Jacobian is "
98  "not implemented.");
99  }
100 
101  void assembleWithJacobian(double const t, double const dt,
102  std::vector<double> const& local_x,
103  std::vector<double> const& local_xdot,
104  const double /*dxdot_dx*/, const double /*dx_dx*/,
105  std::vector<double>& /*local_M_data*/,
106  std::vector<double>& /*local_K_data*/,
107  std::vector<double>& local_rhs_data,
108  std::vector<double>& local_Jac_data) override;
109 
110  void initializeConcrete() override
111  {
112  unsigned const n_integration_points =
113  _integration_method.getNumberOfPoints();
114 
115  for (unsigned ip = 0; ip < n_integration_points; ip++)
116  {
117  auto& ip_data = _ip_data[ip];
118 
120  if (_process_data.initial_stress != nullptr)
121  {
122  ParameterLib::SpatialPosition const x_position{
123  std::nullopt, _element.getID(), ip,
125  ShapeFunctionDisplacement,
127  _element, ip_data.N_u))};
128 
129  ip_data.sigma_eff =
131  DisplacementDim>((*_process_data.initial_stress)(
132  std::numeric_limits<double>::quiet_NaN() /* time
133  independent
134  */
135  ,
136  x_position));
137  }
138 
139  ip_data.pushBackState();
140  }
141  }
142 
143  void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
144  double const /*t*/,
145  double const /*dt*/) override
146  {
147  unsigned const n_integration_points =
148  _integration_method.getNumberOfPoints();
149 
150  for (unsigned ip = 0; ip < n_integration_points; ip++)
151  {
152  _ip_data[ip].pushBackState();
153  }
154  }
155 
157  double const t, double const dt, Eigen::VectorXd const& local_x,
158  Eigen::VectorXd const& local_x_dot) override;
159 
160  void postNonLinearSolverConcrete(std::vector<double> const& local_x,
161  std::vector<double> const& local_xdot,
162  double const t, double const dt,
163  bool const use_monolithic_scheme,
164  int const process_id) override;
165 
166  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
167  const unsigned integration_point) const override
168  {
169  auto const& N_u = _secondary_data.N_u[integration_point];
170 
171  // assumes N is stored contiguously in memory
172  return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
173  }
174 
175  std::vector<double> const& getIntPtDarcyVelocity(
176  const double t,
177  std::vector<GlobalVector*> const& x,
178  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
179  std::vector<double>& cache) const override;
180 
181  // TODO (naumov) This method is same as getIntPtSigma but for arguments and
182  // the ordering of the cache_mat.
183  // There should be only one.
184  std::vector<double> getSigma() const override
185  {
186  constexpr int kelvin_vector_size =
188 
189  return transposeInPlace<kelvin_vector_size>(
190  [this](std::vector<double>& values)
191  { return getIntPtSigma(0, {}, {}, values); });
192  }
193 
194 private:
195  std::size_t setSigma(double const* values)
196  {
197  return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
198  values, _ip_data, &IpData::sigma_eff);
199  }
200 
201  std::vector<double> const& getIntPtSigma(
202  const double /*t*/,
203  std::vector<GlobalVector*> const& /*x*/,
204  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
205  std::vector<double>& cache) const override
206  {
207  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
208  _ip_data, &IpData::sigma_eff, cache);
209  }
210 
211  std::vector<double> getEpsilon() const override
212  {
213  constexpr int kelvin_vector_size =
215 
216  return transposeInPlace<kelvin_vector_size>(
217  [this](std::vector<double>& values)
218  { return getIntPtEpsilon(0, {}, {}, values); });
219  }
220 
221  virtual std::vector<double> const& getIntPtEpsilon(
222  const double /*t*/,
223  std::vector<GlobalVector*> const& /*x*/,
224  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
225  std::vector<double>& cache) const override
226  {
227  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
228  _ip_data, &IpData::eps, cache);
229  }
230 
231 private:
233 
236  using IpData =
238  ShapeMatricesTypePressure, DisplacementDim,
239  ShapeFunctionDisplacement::NPOINTS>;
240  std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
241 
242  IntegrationMethod _integration_method;
248 
249  // The shape function of pressure has the same form with the shape function
250  // of temperature
251  static const int temperature_index = 0;
252  static const int temperature_size = ShapeFunctionPressure::NPOINTS;
253  static const int pressure_index = ShapeFunctionPressure::NPOINTS;
254  static const int pressure_size = ShapeFunctionPressure::NPOINTS;
255  static const int displacement_index = ShapeFunctionPressure::NPOINTS * 2;
256  static const int displacement_size =
257  ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
258 };
259 
260 } // namespace ThermoHydroMechanics
261 } // namespace ProcessLib
262 
#define OGS_FATAL(...)
Definition: Error.h:26
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
ThermoHydroMechanicsLocalAssembler(ThermoHydroMechanicsLocalAssembler &&)=delete
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void postNonLinearSolverConcrete(std::vector< double > const &local_x, std::vector< double > const &local_xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id) override
void postTimestepConcrete(Eigen::VectorXd const &, double const, double const) override
ThermoHydroMechanicsLocalAssembler(ThermoHydroMechanicsLocalAssembler const &)=delete
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double, const double, std::vector< double > &, std::vector< double > &, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
virtual std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
ThermoHydroMechanicsProcessData< DisplacementDim > & _process_data
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_dot) override
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
std::size_t setIPDataInitialConditions(std::string const &name, double const *values, int const integration_order) override
Returns number of read integration points.
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)
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u