OGS
HydroMechanicsLocalAssemblerMatrix.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <vector>
14 
17 
19 
22 
23 namespace ProcessLib
24 {
25 namespace LIE
26 {
27 namespace HydroMechanics
28 {
29 namespace MPL = MaterialPropertyLib;
30 
31 template <typename ShapeFunctionDisplacement,
32  typename ShapeFunctionPressure,
33  typename IntegrationMethod,
34  int GlobalDim>
37 {
38 public:
40  HydroMechanicsLocalAssemblerMatrix const&) = delete;
42  delete;
43 
45  MeshLib::Element const& e,
46  std::size_t const n_variables,
47  std::size_t const local_matrix_size,
48  std::vector<unsigned> const& dofIndex_to_localIndex,
49  bool const is_axially_symmetric,
50  unsigned const integration_order,
52 
53  void preTimestepConcrete(std::vector<double> const& /*local_x*/,
54  double const /*t*/,
55  double const /*delta_t*/) override
56  {
57  for (auto& data : _ip_data)
58  {
59  data.pushBackState();
60  }
61  }
62 
63 protected:
64  void assembleWithJacobianConcrete(double const t, double const dt,
65  Eigen::VectorXd const& local_x,
66  Eigen::VectorXd const& local_x_dot,
67  Eigen::VectorXd& local_rhs,
68  Eigen::MatrixXd& local_Jac) override;
69 
71  double const t, double const dt,
72  Eigen::Ref<const Eigen::VectorXd> const& p,
73  Eigen::Ref<const Eigen::VectorXd> const& p_dot,
74  Eigen::Ref<const Eigen::VectorXd> const& u,
75  Eigen::Ref<const Eigen::VectorXd> const& u_dot,
76  Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_u,
77  Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pu,
78  Eigen::Ref<Eigen::MatrixXd> J_uu, Eigen::Ref<Eigen::MatrixXd> J_up);
79 
81  double const t, double const dt,
82  Eigen::VectorXd const& local_x) override;
83 
85  double const t, double const dt,
86  Eigen::Ref<const Eigen::VectorXd> const& p,
87  Eigen::Ref<const Eigen::VectorXd> const& u);
88 
90  double const t, Eigen::Ref<Eigen::VectorXd> p);
91  void setPressureDotOfInactiveNodes(Eigen::Ref<Eigen::VectorXd> p_dot);
92 
93  // Types for displacement.
96  using BMatricesType =
98 
99  // Types for pressure.
102 
107  GlobalDim,
108  ShapeFunctionDisplacement::NPOINTS>;
109 
110  using GlobalDimVector = Eigen::Matrix<double, GlobalDim, 1>;
111 
113 
114  std::vector<IntegrationPointDataType,
115  Eigen::aligned_allocator<IntegrationPointDataType>>
117 
118  static const int pressure_index = 0;
119  static const int pressure_size = ShapeFunctionPressure::NPOINTS;
120  static const int displacement_index = ShapeFunctionPressure::NPOINTS;
121  static const int displacement_size =
122  ShapeFunctionDisplacement::NPOINTS * GlobalDim;
123  static const int kelvin_vector_size =
125 };
126 
127 } // namespace HydroMechanics
128 } // namespace LIE
129 } // namespace ProcessLib
130 
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assembleWithJacobianConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_dot, Eigen::VectorXd &local_rhs, Eigen::MatrixXd &local_Jac) override
BMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > BMatricesType
ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim > ShapeMatricesTypePressure
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix &&)=delete
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix const &)=delete
void postTimestepConcreteWithVector(double const t, double const dt, Eigen::VectorXd const &local_x) override
void postTimestepConcreteWithBlockVectors(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)
ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > ShapeMatricesTypeDisplacement
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
IntegrationPointDataMatrix< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, GlobalDim, ShapeFunctionDisplacement::NPOINTS > IntegrationPointDataType
void assembleBlockMatricesWithJacobian(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_dot, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_dot, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_u, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pu, Eigen::Ref< Eigen::MatrixXd > J_uu, Eigen::Ref< Eigen::MatrixXd > J_up)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23