OGS
ThermoRichardsMechanicsFEM.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <memory>
14 #include <vector>
15 
16 #include "IntegrationPointData.h"
19 #include "MathLib/KelvinVector.h"
24 #include "ParameterLib/Parameter.h"
29 
30 namespace ProcessLib
31 {
32 namespace ThermoRichardsMechanics
33 {
34 namespace MPL = MaterialPropertyLib;
35 
38 template <typename ShapeMatrixType>
40 {
41  std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
42 };
43 
44 template <typename ShapeFunctionDisplacement, typename ShapeFunction,
45  typename IntegrationMethod, int DisplacementDim>
47  : public LocalAssemblerInterface<DisplacementDim>
48 {
49 public:
52  // Note: temperature variable uses the same shape functions as that are used
53  // by pressure variable.
56 
59 
60  using BMatricesType =
63 
64  using IpData =
66  ShapeMatricesType, DisplacementDim,
67  ShapeFunctionDisplacement::NPOINTS>;
68 
69  static int const KelvinVectorSize =
72 
73  using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
74 
79 
81  MeshLib::Element const& e,
82  std::size_t const /*local_matrix_size*/,
83  bool const is_axially_symmetric,
84  unsigned const integration_order,
86 
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 assembleWithJacobian(double const t, double const dt,
99  std::vector<double> const& local_x,
100  std::vector<double> const& local_xdot,
101  const double /*dxdot_dx*/, const double /*dx_dx*/,
102  std::vector<double>& /*local_M_data*/,
103  std::vector<double>& /*local_K_data*/,
104  std::vector<double>& local_rhs_data,
105  std::vector<double>& local_Jac_data) override;
106 
107  void initializeConcrete() override
108  {
109  unsigned const n_integration_points =
110  integration_method_.getNumberOfPoints();
111 
112  for (unsigned ip = 0; ip < n_integration_points; ip++)
113  {
114  auto& ip_data = ip_data_[ip];
115 
117  if (process_data_.initial_stress != nullptr)
118  {
119  ParameterLib::SpatialPosition const x_position{
120  std::nullopt, element_.getID(), ip,
122  ShapeFunctionDisplacement,
124  element_, ip_data.N_u))};
125 
126  ip_data.sigma_eff =
128  DisplacementDim>((*process_data_.initial_stress)(
129  std::numeric_limits<
130  double>::quiet_NaN() /* time independent */,
131  x_position));
132  }
133 
134  ip_data.pushBackState();
135  }
136  }
137 
138  void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
139  double const /*t*/,
140  double const /*dt*/) override
141  {
142  unsigned const n_integration_points =
143  integration_method_.getNumberOfPoints();
144 
145  for (unsigned ip = 0; ip < n_integration_points; ip++)
146  {
147  ip_data_[ip].pushBackState();
148  }
149  }
150 
152  double const t, double const dt, Eigen::VectorXd const& local_x,
153  Eigen::VectorXd const& local_x_dot) override;
154 
155  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
156  const unsigned integration_point) const override
157  {
158  auto const& N_u = secondary_data_.N_u[integration_point];
159 
160  // assumes N is stored contiguously in memory
161  return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
162  }
163 
164  std::vector<double> getSigma() const override;
165 
166  std::vector<double> const& getIntPtDarcyVelocity(
167  const double t,
168  std::vector<GlobalVector*> const& x,
169  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
170  std::vector<double>& cache) const override;
171 
172  std::vector<double> getSaturation() const override;
173  std::vector<double> const& getIntPtSaturation(
174  const double t,
175  std::vector<GlobalVector*> const& x,
176  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
177  std::vector<double>& cache) const override;
178 
179  std::vector<double> getPorosity() const override;
180  std::vector<double> const& getIntPtPorosity(
181  const double t,
182  std::vector<GlobalVector*> const& x,
183  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
184  std::vector<double>& cache) const override;
185 
186  std::vector<double> getTransportPorosity() const override;
187  std::vector<double> const& getIntPtTransportPorosity(
188  const double t,
189  std::vector<GlobalVector*> const& x,
190  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
191  std::vector<double>& cache) const override;
192 
193  std::vector<double> const& getIntPtSigma(
194  const double t,
195  std::vector<GlobalVector*> const& x,
196  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
197  std::vector<double>& cache) const override;
198 
199  std::vector<double> getSwellingStress() const override;
200  std::vector<double> const& getIntPtSwellingStress(
201  const double t,
202  std::vector<GlobalVector*> const& x,
203  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
204  std::vector<double>& cache) const override;
205 
206  std::vector<double> getEpsilon() const override;
207  std::vector<double> const& getIntPtEpsilon(
208  const double t,
209  std::vector<GlobalVector*> const& x,
210  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
211  std::vector<double>& cache) const override;
212 
213  std::vector<double> const& getIntPtDryDensitySolid(
214  const double t,
215  std::vector<GlobalVector*> const& x,
216  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
217  std::vector<double>& cache) const override;
218 
219 private:
220  unsigned getNumberOfIntegrationPoints() const override;
221 
223  DisplacementDim>::MaterialStateVariables const&
224  getMaterialStateVariablesAt(unsigned integration_point) const override;
225 
226 private:
228 
229  std::vector<IpData, Eigen::aligned_allocator<IpData>> ip_data_;
230 
231  IntegrationMethod integration_method_;
237 
238  static const int temperature_index = 0;
239  static const int temperature_size = ShapeFunction::NPOINTS;
240  static const int pressure_index = temperature_size;
241  static const int pressure_size = ShapeFunction::NPOINTS;
242  static const int displacement_index = 2 * ShapeFunction::NPOINTS;
243  static const int displacement_size =
244  ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
245 };
246 
247 } // namespace ThermoRichardsMechanics
248 } // namespace ProcessLib
249 
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
VectorType< _kelvin_vector_size > KelvinVectorType
Definition: BMatrixPolicy.h:48
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
std::vector< double > const & getIntPtDryDensitySolid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
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
ThermoRichardsMechanicsProcessData< DisplacementDim > & process_data_
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_dot) override
std::vector< double > const & getIntPtPorosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
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
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > secondary_data_
void setInitialConditionsConcrete(std::vector< double > const &local_x, double const t, bool const use_monolithic_scheme, 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
std::vector< IpData, Eigen::aligned_allocator< IpData > > ip_data_
ThermoRichardsMechanicsLocalAssembler(ThermoRichardsMechanicsLocalAssembler const &)=delete
std::size_t setIPDataInitialConditions(std::string const &name, double const *values, int const integration_order) override
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
ThermoRichardsMechanicsLocalAssembler(ThermoRichardsMechanicsLocalAssembler &&)=delete
std::vector< double > const & getIntPtTransportPorosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtSwellingStress(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 &, double const, double const) override
std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
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