OGS
TH2MFEM.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <memory>
14 #include <vector>
15 
16 #include "ConstitutiveVariables.h"
17 #include "IntegrationPointData.h"
21 #include "MathLib/KelvinVector.h"
26 #include "ParameterLib/Parameter.h"
32 #include "TH2MProcessData.h"
33 
34 namespace ProcessLib
35 {
36 namespace TH2M
37 {
40 template <typename ShapeMatrixType>
42 {
43  std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
44 };
45 
46 template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
47  typename IntegrationMethod, int DisplacementDim>
49 {
50 public:
53 
56 
57  template <int N>
58  using VectorType =
59  typename ShapeMatricesTypePressure::template VectorType<N>;
60 
61  template <int M, int N>
62  using MatrixType =
63  typename ShapeMatricesTypePressure::template MatrixType<M, N>;
64 
67 
70 
71  static int const KelvinVectorSize =
73  using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
74 
76 
79 
81  std::size_t const /*local_matrix_size*/,
82  bool const is_axially_symmetric,
83  unsigned const integration_order,
84  TH2MProcessData<DisplacementDim>& process_data);
85 
86 private:
88  std::size_t setIPDataInitialConditions(
89  std::string const& name,
90  double const* values,
91  int const integration_order) override;
92 
93  void setInitialConditionsConcrete(std::vector<double> const& local_x,
94  double const t,
95  bool const /*use_monolithic_scheme*/,
96  int const /*process_id*/) override;
97 
98  void assemble(double const /*t*/, double const /*dt*/,
99  std::vector<double> const& /*local_x*/,
100  std::vector<double> const& /*local_xdot*/,
101  std::vector<double>& /*local_M_data*/,
102  std::vector<double>& /*local_K_data*/,
103  std::vector<double>& /*local_rhs_data*/) override;
104 
105  void assembleWithJacobian(double const t, double const dt,
106  std::vector<double> const& local_x,
107  std::vector<double> const& local_xdot,
108  const double /*dxdot_dx*/, const double /*dx_dx*/,
109  std::vector<double>& /*local_M_data*/,
110  std::vector<double>& /*local_K_data*/,
111  std::vector<double>& local_rhs_data,
112  std::vector<double>& local_Jac_data) override;
113 
114  void initializeConcrete() override
115  {
116  unsigned const n_integration_points =
117  _integration_method.getNumberOfPoints();
118 
119  for (unsigned ip = 0; ip < n_integration_points; ip++)
120  {
121  auto& ip_data = _ip_data[ip];
122 
124  if (_process_data.initial_stress != nullptr)
125  {
126  ParameterLib::SpatialPosition const x_position{
127  std::nullopt, _element.getID(), ip,
129  ShapeFunctionDisplacement,
131  _element, ip_data.N_u))};
132 
133  ip_data.sigma_eff =
135  DisplacementDim>((*_process_data.initial_stress)(
136  std::numeric_limits<
137  double>::quiet_NaN() /* time independent */,
138  x_position));
139  }
140 
141  ip_data.pushBackState();
142  }
143  }
144 
145  void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
146  double const /*t*/,
147  double const /*dt*/) override
148  {
149  unsigned const n_integration_points =
150  _integration_method.getNumberOfPoints();
151 
152  for (unsigned ip = 0; ip < n_integration_points; ip++)
153  {
154  _ip_data[ip].pushBackState();
155  }
156  }
157 
159  double const t, double const dt, Eigen::VectorXd const& local_x,
160  Eigen::VectorXd const& local_x_dot) override;
161 
162  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
163  const unsigned integration_point) const override
164  {
165  auto const& N_u = _secondary_data.N_u[integration_point];
166 
167  // assumes N is stored contiguously in memory
168  return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
169  }
170 
171  std::vector<double> const& getIntPtDarcyVelocityGas(
172  const double t,
173  std::vector<GlobalVector*> const& x,
174  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
175  std::vector<double>& cache) const override;
176  std::vector<double> const& getIntPtDarcyVelocityLiquid(
177  const double t,
178  std::vector<GlobalVector*> const& x,
179  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
180  std::vector<double>& cache) const override;
181 
182  std::vector<ConstitutiveVariables<DisplacementDim>>
183  updateConstitutiveVariables(Eigen::VectorXd const& local_x,
184  Eigen::VectorXd const& local_x_dot,
185  double const t, double const dt);
186 
187  std::size_t setSigma(double const* values)
188  {
189  return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
190  values, _ip_data, &IpData::sigma_eff);
191  }
192 
193  // TODO (naumov) This method is same as getIntPtSigma but for arguments and
194  // the ordering of the cache_mat.
195  // There should be only one.
196  std::vector<double> getSigma() const override
197  {
198  {
199  constexpr int kelvin_vector_size =
201  DisplacementDim);
202 
203  return transposeInPlace<kelvin_vector_size>(
204  [this](std::vector<double>& values)
205  { return getIntPtSigma(0, {}, {}, values); });
206  }
207  }
208 
209  std::vector<double> const& getIntPtSigma(
210  const double /*t*/,
211  std::vector<GlobalVector*> const& /*x*/,
212  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
213  std::vector<double>& cache) const override
214  {
215  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
216  _ip_data, &IpData::sigma_eff, cache);
217  }
218 
219  std::vector<double> getEpsilon() const override
220  {
221  constexpr int kelvin_vector_size =
223 
224  return transposeInPlace<kelvin_vector_size>(
225  [this](std::vector<double>& values)
226  { return getIntPtEpsilon(0, {}, {}, values); });
227  }
228 
229  virtual std::vector<double> const& getIntPtEpsilon(
230  const double /*t*/,
231  std::vector<GlobalVector*> const& /*x*/,
232  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
233  std::vector<double>& cache) const override
234  {
235  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
236  _ip_data, &IpData::eps, cache);
237  }
238 
239  // TODO: Here is some refactoring potential. All secondary variables could
240  // be stored in some container to avoid defining one method for each
241  // variable.
242 
243  virtual std::vector<double> const& getIntPtLiquidDensity(
244  const double /*t*/,
245  std::vector<GlobalVector*> const& /*x*/,
246  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
247  std::vector<double>& cache) const override
248  {
250  &IpData::rhoLR, cache);
251  }
252 
253  virtual std::vector<double> const& getIntPtGasDensity(
254  const double /*t*/,
255  std::vector<GlobalVector*> const& /*x*/,
256  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
257  std::vector<double>& cache) const override
258  {
260  &IpData::rhoGR, cache);
261  }
262 
263  virtual std::vector<double> const& getIntPtSolidDensity(
264  const double /*t*/,
265  std::vector<GlobalVector*> const& /*x*/,
266  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
267  std::vector<double>& cache) const override
268  {
270  &IpData::rhoSR, cache);
271  }
272 
273  virtual std::vector<double> const& getIntPtVapourPressure(
274  const double /*t*/,
275  std::vector<GlobalVector*> const& /*x*/,
276  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
277  std::vector<double>& cache) const override
278  {
280  &IpData::pWGR, cache);
281  }
282 
283  virtual std::vector<double> const& getIntPtPorosity(
284  const double /*t*/,
285  std::vector<GlobalVector*> const& /*x*/,
286  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
287  std::vector<double>& cache) const override
288  {
290  cache);
291  }
292 
293  std::vector<double> getSaturation() const override
294  {
295  std::vector<double> result;
296  getIntPtSaturation(0, {}, {}, result);
297  return result;
298  }
299 
300  virtual std::vector<double> const& getIntPtSaturation(
301  const double /*t*/,
302  std::vector<GlobalVector*> const& /*x*/,
303  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
304  std::vector<double>& cache) const override
305  {
307  cache);
308  }
309 
310  virtual std::vector<double> const& getIntPtMoleFractionGas(
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  {
317  &IpData::xnCG, cache);
318  }
319  virtual std::vector<double> const& getIntPtMassFractionGas(
320  const double /*t*/,
321  std::vector<GlobalVector*> const& /*x*/,
322  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
323  std::vector<double>& cache) const override
324  {
326  &IpData::xmCG, cache);
327  }
328  virtual std::vector<double> const& getIntPtMassFractionLiquid(
329  const double /*t*/,
330  std::vector<GlobalVector*> const& /*x*/,
331  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
332  std::vector<double>& cache) const override
333  {
335  &IpData::xmWL, cache);
336  }
337 
338  virtual std::vector<double> const& getIntPtRelativePermeabilityGas(
339  const double /*t*/,
340  std::vector<GlobalVector*> const& /*x*/,
341  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
342  std::vector<double>& cache) const override
343  {
345  _ip_data, &IpData::k_rel_G, cache);
346  }
347  virtual std::vector<double> const& getIntPtRelativePermeabilityLiquid(
348  const double /*t*/,
349  std::vector<GlobalVector*> const& /*x*/,
350  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
351  std::vector<double>& cache) const override
352  {
354  _ip_data, &IpData::k_rel_L, cache);
355  }
356 
357  virtual std::vector<double> const& getIntPtEnthalpyGas(
358  const double /*t*/,
359  std::vector<GlobalVector*> const& /*x*/,
360  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
361  std::vector<double>& cache) const override
362  {
364  cache);
365  }
366  virtual std::vector<double> const& getIntPtEnthalpyLiquid(
367  const double /*t*/,
368  std::vector<GlobalVector*> const& /*x*/,
369  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
370  std::vector<double>& cache) const override
371  {
373  cache);
374  }
375  virtual std::vector<double> const& getIntPtEnthalpySolid(
376  const double /*t*/,
377  std::vector<GlobalVector*> const& /*x*/,
378  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
379  std::vector<double>& cache) const override
380  {
382  cache);
383  }
384 
385 private:
387 
390  using IpData =
392  ShapeMatricesTypePressure, DisplacementDim,
393  ShapeFunctionDisplacement::NPOINTS>;
394  std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
395 
396  IntegrationMethod _integration_method;
402 
403  // The shape function of pressure has the same form with the shape function
404  // of temperature
405  static const int gas_pressure_index = 0;
406  static const int gas_pressure_size = ShapeFunctionPressure::NPOINTS;
407  static const int capillary_pressure_index = ShapeFunctionPressure::NPOINTS;
408  static const int capillary_pressure_size = ShapeFunctionPressure::NPOINTS;
409  static const int temperature_index = 2 * ShapeFunctionPressure::NPOINTS;
410  static const int temperature_size = ShapeFunctionPressure::NPOINTS;
411  static const int displacement_index = ShapeFunctionPressure::NPOINTS * 3;
412  static const int displacement_size =
413  ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
414 
415  static const int C_index = 0;
416  static const int C_size = ShapeFunctionPressure::NPOINTS;
417  static const int W_index = ShapeFunctionPressure::NPOINTS;
418  static const int W_size = ShapeFunctionPressure::NPOINTS;
419 };
420 
421 } // namespace TH2M
422 } // namespace ProcessLib
423 
424 #include "TH2MFEM-impl.h"
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
static const int temperature_index
Definition: TH2MFEM.h:409
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:209
virtual std::vector< double > const & getIntPtMassFractionGas(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:319
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
Definition: TH2MFEM.h:66
std::vector< double > const & getIntPtDarcyVelocityLiquid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Definition: TH2MFEM.h:401
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
Definition: TH2MFEM.h:55
std::size_t setSigma(double const *values)
Definition: TH2MFEM.h:187
std::vector< double > getSigma() const override
Definition: TH2MFEM.h:196
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_dot) override
std::vector< double > getSaturation() const override
Definition: TH2MFEM.h:293
typename ShapeMatricesTypePressure::template MatrixType< M, N > MatrixType
Definition: TH2MFEM.h:63
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
Definition: TH2MFEM.h:389
MeshLib::Element const & _element
Definition: TH2MFEM.h:397
virtual std::vector< double > const & getIntPtMoleFractionGas(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:310
virtual std::vector< double > const & getIntPtSaturation(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:300
static const int gas_pressure_index
Definition: TH2MFEM.h:405
virtual std::vector< double > const & getIntPtRelativePermeabilityGas(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:338
TH2MLocalAssembler(TH2MLocalAssembler &&)=delete
std::vector< double > const & getIntPtDarcyVelocityGas(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
Definition: TH2MFEM.h:52
typename ShapeMatricesTypePressure::template VectorType< N > VectorType
Definition: TH2MFEM.h:59
virtual std::vector< double > const & getIntPtRelativePermeabilityLiquid(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:347
static const int gas_pressure_size
Definition: TH2MFEM.h:406
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
Definition: TH2MFEM-impl.h:772
void setInitialConditionsConcrete(std::vector< double > const &local_x, double const t, bool const, int const) override
Definition: TH2MFEM-impl.h:744
std::vector< double > getEpsilon() const override
Definition: TH2MFEM.h:219
static const int displacement_size
Definition: TH2MFEM.h:412
virtual std::vector< double > const & getIntPtEnthalpySolid(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:375
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
Definition: TH2MFEM.h:69
virtual std::vector< double > const & getIntPtPorosity(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:283
TH2MProcessData< DisplacementDim > & _process_data
Definition: TH2MFEM.h:386
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
Definition: TH2MFEM.h:394
std::vector< ConstitutiveVariables< DisplacementDim > > updateConstitutiveVariables(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_dot, double const t, double const dt)
Definition: TH2MFEM-impl.h:100
std::size_t setIPDataInitialConditions(std::string const &name, double const *values, int const integration_order) override
Definition: TH2MFEM-impl.h:699
virtual std::vector< double > const & getIntPtSolidDensity(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:263
virtual std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:229
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
TH2MLocalAssembler(TH2MLocalAssembler const &)=delete
virtual std::vector< double > const & getIntPtVapourPressure(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:273
virtual std::vector< double > const & getIntPtMassFractionLiquid(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:328
void postTimestepConcrete(Eigen::VectorXd const &, double const, double const) override
Definition: TH2MFEM.h:145
static const int capillary_pressure_index
Definition: TH2MFEM.h:407
void initializeConcrete() override
Definition: TH2MFEM.h:114
virtual std::vector< double > const & getIntPtEnthalpyLiquid(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:366
virtual std::vector< double > const & getIntPtEnthalpyGas(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:357
static const int temperature_size
Definition: TH2MFEM.h:410
static const int displacement_index
Definition: TH2MFEM.h:411
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
Definition: TH2MFEM.h:73
virtual std::vector< double > const & getIntPtLiquidDensity(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:243
IntegrationMethod _integration_method
Definition: TH2MFEM.h:396
virtual std::vector< double > const & getIntPtGasDensity(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Definition: TH2MFEM.h:253
static const int capillary_pressure_size
Definition: TH2MFEM.h:408
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
Definition: TH2MFEM.h:162
static int const KelvinVectorSize
Definition: TH2MFEM.h:71
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)
std::vector< double > const & getIntegrationPointScalarData(std::vector< IntegrationPointData, Eigen::aligned_allocator< IntegrationPointData >> const &ip_data, MemberType member, std::vector< double > &cache)
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
BMatricesType::KelvinVectorType eps
BMatricesType::KelvinVectorType sigma_eff
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
Definition: TH2MFEM.h:43