OGS
RichardsMechanicsFEM.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"
30 
31 namespace ProcessLib
32 {
33 namespace RichardsMechanics
34 {
35 namespace MPL = MaterialPropertyLib;
36 
39 template <typename ShapeMatrixType>
41 {
42  std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
43 };
44 
45 template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
46  typename IntegrationMethod, int DisplacementDim>
48  : public LocalAssemblerInterface<DisplacementDim>
49 {
50 public:
55 
58 
59  using BMatricesType =
62 
63  using IpData =
65  ShapeMatricesTypePressure, DisplacementDim,
66  ShapeFunctionDisplacement::NPOINTS>;
67 
68  static int const KelvinVectorSize =
71 
72  using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
73 
75  delete;
77 
79  MeshLib::Element const& e,
80  std::size_t const /*local_matrix_size*/,
81  bool const is_axially_symmetric,
82  unsigned const integration_order,
84 
86  std::size_t setIPDataInitialConditions(
87  std::string const& name,
88  double const* values,
89  int const integration_order) override;
90 
91  void setInitialConditionsConcrete(std::vector<double> const& local_x,
92  double const t,
93  bool const use_monolithic_scheme,
94  int const process_id) override;
95 
96  void assemble(double const t, double const dt,
97  std::vector<double> const& local_x,
98  std::vector<double> const& local_xdot,
99  std::vector<double>& local_M_data,
100  std::vector<double>& local_K_data,
101  std::vector<double>& local_rhs_data) override;
102 
103  void assembleWithJacobian(double const t, double const dt,
104  std::vector<double> const& local_x,
105  std::vector<double> const& local_xdot,
106  const double /*dxdot_dx*/, const double /*dx_dx*/,
107  std::vector<double>& /*local_M_data*/,
108  std::vector<double>& /*local_K_data*/,
109  std::vector<double>& local_rhs_data,
110  std::vector<double>& local_Jac_data) override;
111 
113  double const t, double const dt, Eigen::VectorXd const& local_x,
114  Eigen::VectorXd const& local_xdot, const double dxdot_dx,
115  const double dx_dx, int const process_id,
116  std::vector<double>& local_M_data, std::vector<double>& local_K_data,
117  std::vector<double>& local_b_data,
118  std::vector<double>& local_Jac_data) override;
119 
120  void initializeConcrete() override
121  {
122  unsigned const n_integration_points =
123  _integration_method.getNumberOfPoints();
124 
125  for (unsigned ip = 0; ip < n_integration_points; ip++)
126  {
127  auto& ip_data = _ip_data[ip];
128 
130  if (_process_data.initial_stress != nullptr)
131  {
132  ParameterLib::SpatialPosition const x_position{
133  std::nullopt, _element.getID(), ip,
135  ShapeFunctionDisplacement,
137  _element, ip_data.N_u))};
138 
139  ip_data.sigma_eff =
141  DisplacementDim>((*_process_data.initial_stress)(
142  std::numeric_limits<
143  double>::quiet_NaN() /* time independent */,
144  x_position));
145  }
146 
147  ip_data.pushBackState();
148  }
149  }
150 
151  void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
152  double const /*t*/,
153  double const /*dt*/) override
154  {
155  unsigned const n_integration_points =
156  _integration_method.getNumberOfPoints();
157 
158  for (unsigned ip = 0; ip < n_integration_points; ip++)
159  {
160  _ip_data[ip].pushBackState();
161  }
162  }
163 
165  double const t, double const dt, Eigen::VectorXd const& local_x,
166  Eigen::VectorXd const& local_x_dot) override;
167 
168  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
169  const unsigned integration_point) const override
170  {
171  auto const& N_u = _secondary_data.N_u[integration_point];
172 
173  // assumes N is stored contiguously in memory
174  return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
175  }
176 
177  std::vector<double> getSigma() const override;
178 
179  std::vector<double> const& getIntPtDarcyVelocity(
180  const double t,
181  std::vector<GlobalVector*> const& x,
182  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
183  std::vector<double>& cache) const override;
184 
185  std::vector<double> getSaturation() const override;
186  std::vector<double> const& getIntPtSaturation(
187  const double t,
188  std::vector<GlobalVector*> const& x,
189  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
190  std::vector<double>& cache) const override;
191 
192  std::vector<double> getMicroSaturation() const override;
193  std::vector<double> const& getIntPtMicroSaturation(
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> getMicroPressure() const override;
200  std::vector<double> const& getIntPtMicroPressure(
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> getPorosity() const override;
207  std::vector<double> const& getIntPtPorosity(
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> getTransportPorosity() const override;
214  std::vector<double> const& getIntPtTransportPorosity(
215  const double t,
216  std::vector<GlobalVector*> const& x,
217  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
218  std::vector<double>& cache) const override;
219 
220  std::vector<double> const& getIntPtSigma(
221  const double t,
222  std::vector<GlobalVector*> const& x,
223  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
224  std::vector<double>& cache) const override;
225 
226  std::vector<double> getSwellingStress() const override;
227  std::vector<double> const& getIntPtSwellingStress(
228  const double t,
229  std::vector<GlobalVector*> const& x,
230  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
231  std::vector<double>& cache) const override;
232 
233  std::vector<double> getEpsilon() const override;
234  std::vector<double> const& getIntPtEpsilon(
235  const double t,
236  std::vector<GlobalVector*> const& x,
237  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
238  std::vector<double>& cache) const override;
239 
240  std::vector<double> getMaterialStateVariableInternalState(
241  std::function<BaseLib::DynamicSpan<double>(
243  MaterialStateVariables&)> const& get_values_span,
244  int const& n_components) const override;
245 
246  std::vector<double> const& getIntPtDryDensitySolid(
247  const double t,
248  std::vector<GlobalVector*> const& x,
249  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
250  std::vector<double>& cache) const override;
251 
252 private:
279  double const t, double const dt, Eigen::VectorXd const& local_x,
280  Eigen::VectorXd const& local_xdot, const double dxdot_dx,
281  const double dx_dx, std::vector<double>& local_M_data,
282  std::vector<double>& local_K_data, std::vector<double>& local_b_data,
283  std::vector<double>& local_Jac_data);
284 
314  double const t, double const dt, Eigen::VectorXd const& local_x,
315  Eigen::VectorXd const& local_xdot, const double dxdot_dx,
316  const double dx_dx, std::vector<double>& local_M_data,
317  std::vector<double>& local_K_data, std::vector<double>& local_b_data,
318  std::vector<double>& local_Jac_data);
319 
320  unsigned getNumberOfIntegrationPoints() const override;
321 
323  DisplacementDim>::MaterialStateVariables const&
324  getMaterialStateVariablesAt(unsigned integration_point) const override;
325 
326 private:
328 
329  std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
330 
331  IntegrationMethod _integration_method;
337 
338  static const int pressure_index = 0;
339  static const int pressure_size = ShapeFunctionPressure::NPOINTS;
340  static const int displacement_index = ShapeFunctionPressure::NPOINTS;
341  static const int displacement_size =
342  ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
343 };
344 
345 } // namespace RichardsMechanics
346 } // namespace ProcessLib
347 
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
void postTimestepConcrete(Eigen::VectorXd const &, double const, double const) override
RichardsMechanicsProcessData< DisplacementDim > & _process_data
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
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
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
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
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
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
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
void assembleWithJacobianForDeformationEquations(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd 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_b_data, std::vector< double > &local_Jac_data)
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_dot) override
std::vector< double > getMaterialStateVariableInternalState(std::function< BaseLib::DynamicSpan< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const override
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 & getIntPtSaturation(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 & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
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 & getIntPtMicroPressure(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
RichardsMechanicsLocalAssembler(RichardsMechanicsLocalAssembler &&)=delete
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
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
std::vector< double > const & getIntPtMicroSaturation(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 t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data) override
void assembleWithJacobianForPressureEquations(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd 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_b_data, std::vector< double > &local_Jac_data)
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
std::size_t setIPDataInitialConditions(std::string const &name, double const *values, int const integration_order) override
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
RichardsMechanicsLocalAssembler(RichardsMechanicsLocalAssembler const &)=delete
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
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
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
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u