OGS
LocalAssemblerInterface.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <Eigen/Dense>
14 
15 #include "NumLib/NumericsConfig.h"
16 #include "MathLib/Point3d.h"
17 
18 namespace NumLib
19 {
20 class LocalToGlobalIndexMap;
21 } // NumLib
22 
23 namespace ProcessLib
24 {
25 struct CoupledSolutionsForStaggeredScheme;
26 struct LocalCoupledSolutions;
27 
34 {
35 public:
36  virtual ~LocalAssemblerInterface() = default;
37 
38  void setInitialConditions(std::size_t const mesh_item_id,
39  NumLib::LocalToGlobalIndexMap const& dof_table,
40  GlobalVector const& x, double const t,
41  bool const use_monolithic_scheme,
42  int const process_id);
43 
44  virtual void initialize(std::size_t const mesh_item_id,
45  NumLib::LocalToGlobalIndexMap const& dof_table);
46 
47  virtual void preAssemble(double const /*t*/, double const /*dt*/,
48  std::vector<double> const& /*local_x*/){};
49 
50  virtual void assemble(double const t, double const dt,
51  std::vector<double> const& local_x,
52  std::vector<double> const& local_xdot,
53  std::vector<double>& local_M_data,
54  std::vector<double>& local_K_data,
55  std::vector<double>& local_b_data);
56 
57  virtual void assembleForStaggeredScheme(double const t, double const dt,
58  Eigen::VectorXd const& local_x,
59  Eigen::VectorXd const& local_xdot,
60  int const process_id,
61  std::vector<double>& local_M_data,
62  std::vector<double>& local_K_data,
63  std::vector<double>& local_b_data);
64 
65  virtual void assembleWithJacobian(double const t, double const dt,
66  std::vector<double> const& local_x,
67  std::vector<double> const& local_xdot,
68  const double dxdot_dx, const double dx_dx,
69  std::vector<double>& local_M_data,
70  std::vector<double>& local_K_data,
71  std::vector<double>& local_b_data,
72  std::vector<double>& local_Jac_data);
73 
75  double const t, double const dt, Eigen::VectorXd const& local_x,
76  Eigen::VectorXd const& local_xdot, const double dxdot_dx,
77  const double dx_dx, int const process_id,
78  std::vector<double>& local_M_data, std::vector<double>& local_K_data,
79  std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
80 
81  virtual void computeSecondaryVariable(
82  std::size_t const mesh_item_id,
83  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
84  double const t, double const dt, std::vector<GlobalVector*> const& x,
85  GlobalVector const& x_dot, int const process_id);
86 
87  virtual void preTimestep(std::size_t const mesh_item_id,
88  NumLib::LocalToGlobalIndexMap const& dof_table,
89  GlobalVector const& x, double const t,
90  double const delta_t);
91 
92  virtual void postTimestep(
93  std::size_t const mesh_item_id,
94  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
95  std::vector<GlobalVector*> const& x, double const t, double const dt);
96 
97  void postNonLinearSolver(std::size_t const mesh_item_id,
98  NumLib::LocalToGlobalIndexMap const& dof_table,
99  GlobalVector const& x, GlobalVector const& xdot,
100  double const t, double const dt,
101  bool const use_monolithic_scheme,
102  int const process_id);
103 
104  virtual std::vector<double> interpolateNodalValuesToIntegrationPoints(
105  std::vector<double> const& /*local_x*/)
106  {
107  return {};
108  }
109 
113  virtual Eigen::Vector3d getFlux(
114  MathLib::Point3d const& /*p_local_coords*/,
115  double const /*t*/,
116  std::vector<double> const& /*local_x*/) const
117  {
118  return Eigen::Vector3d{};
119  }
120 
122  virtual Eigen::Vector3d getFlux(
123  MathLib::Point3d const& /*p_local_coords*/,
124  double const /*t*/,
125  std::vector<std::vector<double>> const& /*local_xs*/) const
126  {
127  return Eigen::Vector3d{};
128  }
129 
130 private:
132  std::vector<double> const& /*local_x*/, double const /*t*/,
133  bool const /*use_monolithic_scheme*/, int const /*process_id*/)
134  {
135  }
136 
137  virtual void initializeConcrete() {}
138 
139  virtual void preTimestepConcrete(std::vector<double> const& /*local_x*/,
140  double const /*t*/, double const /*dt*/)
141  {
142  }
143 
144  virtual void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
145  double const /*t*/, double const /*dt*/)
146  {
147  }
148 
150  std::vector<double> const& /*local_x*/,
151  std::vector<double> const& /*local_xdot*/, double const /*t*/,
152  double const /*dt*/, bool const /*use_monolithic_scheme*/,
153  int const /*process_id*/)
154  {
155  }
156 
158  double const /*t*/,
159  double const /*dt*/,
160  Eigen::VectorXd const& /*local_x*/,
161  Eigen::VectorXd const& /*local_x_dot*/)
162  {
163  }
164 
166  double const /*t*/, std::vector<std::vector<double>> const&
167  /*coupled_local_solutions*/)
168  {
169  }
170 };
171 
172 } // namespace ProcessLib
Definition of the Point3d class.
Global vector based on Eigen vector.
Definition: EigenVector.h:26
virtual void preTimestepConcrete(std::vector< double > const &, double const, double const)
virtual ~LocalAssemblerInterface()=default
virtual Eigen::Vector3d getFlux(MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const
Fits to staggered scheme.
virtual void assembleForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void setInitialConditionsConcrete(std::vector< double > const &, double const, bool const, int const)
virtual void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > 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)
virtual void computeSecondaryVariableConcrete(double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &)
virtual Eigen::Vector3d getFlux(MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual void preAssemble(double const, double const, std::vector< double > const &)
virtual void postNonLinearSolverConcrete(std::vector< double > const &, std::vector< double > const &, double const, double const, bool const, int const)
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
void setInitialConditions(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id)
virtual 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_b_data)
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
void postNonLinearSolver(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
virtual void postTimestepConcrete(Eigen::VectorXd const &, double const, double const)
virtual 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)
virtual void computeSecondaryVariableWithCoupledProcessConcrete(double const, std::vector< std::vector< double >> const &)
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints(std::vector< double > const &)