OGS
Process.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <string>
14 #include <tuple>
15 
20 #include "ParameterLib/Parameter.h"
26 #include "ProcessVariable.h"
27 #include "VectorMatrixAssembler.h"
28 
29 namespace MeshLib
30 {
31 class Mesh;
32 }
33 
34 namespace ProcessLib
35 {
36 struct CoupledSolutionsForStaggeredScheme;
37 
38 class Process
39  : public NumLib::ODESystem< // TODO: later on use a simpler ODE system
40  NumLib::ODESystemTag::FirstOrderImplicitQuasilinear,
41  NumLib::NonlinearSolverTag::Newton>
42 {
43 public:
46 
47  Process(std::string name_,
48  MeshLib::Mesh& mesh,
49  std::unique_ptr<AbstractJacobianAssembler>&& jacobian_assembler,
50  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const&
51  parameters,
52  unsigned const integration_order,
53  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
54  process_variables,
55  SecondaryVariableCollection&& secondary_variables,
56  const bool use_monolithic_scheme = true);
57 
59  void preTimestep(std::vector<GlobalVector*> const& x, const double t,
60  const double delta_t, const int process_id);
61 
63  void postTimestep(std::vector<GlobalVector*> const& x, const double t,
64  const double delta_t, int const process_id);
65 
68  void postNonLinearSolver(GlobalVector const& x, GlobalVector const& xdot,
69  const double t, double const dt,
70  int const process_id);
71 
72  void preIteration(const unsigned iter, GlobalVector const& x) final;
73 
75  void computeSecondaryVariable(double const t,
76  double const dt,
77  std::vector<GlobalVector*> const& x,
78  GlobalVector const& x_dot,
79  int const process_id);
80 
82 
83  void initialize();
84 
86  std::vector<GlobalVector*>& process_solutions,
87  std::vector<GlobalVector*> const& process_solutions_prev,
88  double const t,
89  int const process_id);
90 
92  const int process_id) const override;
93 
95  CoupledSolutionsForStaggeredScheme* const coupled_solutions)
96  {
97  _coupled_solutions = coupled_solutions;
98  }
99 
100  void updateDeactivatedSubdomains(double const time, const int process_id);
101 
104  int const /*process_id*/)
105  {
106  }
107 
109  const double /*t*/,
110  std::vector<GlobalVector*> const& /*integration_point_values_vectors*/,
111  std::vector<GlobalVector*>& /*nodal_values_vectors*/)
112  {
113  }
114 
115  void preAssemble(const double t, double const dt,
116  GlobalVector const& x) final;
117  void assemble(const double t, double const dt,
118  std::vector<GlobalVector*> const& x,
119  std::vector<GlobalVector*> const& xdot, int const process_id,
120  GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) final;
121 
122  void assembleWithJacobian(const double t, double const dt,
123  std::vector<GlobalVector*> const& x,
124  std::vector<GlobalVector*> const& xdot,
125  const double dxdot_dx, const double dx_dx,
126  int const process_id, GlobalMatrix& M,
127  GlobalMatrix& K, GlobalVector& b,
128  GlobalMatrix& Jac) final;
129 
130  std::vector<NumLib::IndexValueVector<GlobalIndexType>> const*
131  getKnownSolutions(double const t, GlobalVector const& x,
132  int const process_id) const final
133  {
134  return _boundary_conditions[process_id].getKnownSolutions(t, x);
135  }
136 
138  const int /*process_id*/) const
139  {
141  }
142 
143  MeshLib::Mesh& getMesh() const { return _mesh; }
144  std::vector<std::reference_wrapper<ProcessVariable>> const&
145  getProcessVariables(const int process_id) const
146  {
147  return _process_variables[process_id];
148  }
149 
151  {
152  return _secondary_variables;
153  }
154 
155  std::vector<std::unique_ptr<IntegrationPointWriter>> const*
157  {
158  if (mesh == _mesh)
159  {
161  }
162  return nullptr;
163  }
164 
165  // Used as a call back for CalculateSurfaceFlux process.
166 
167  virtual Eigen::Vector3d getFlux(
168  std::size_t /*element_id*/,
169  MathLib::Point3d const& /*p*/,
170  double const /*t*/,
171  std::vector<GlobalVector*> const& /*x*/) const
172  {
173  return Eigen::Vector3d{};
174  }
175 
176  virtual void solveReactionEquation(
177  std::vector<GlobalVector*>& /*x*/,
178  std::vector<GlobalVector*> const& /*x_prev*/, double const /*t*/,
179  double const /*dt*/, NumLib::EquationSystem& /*ode_sys*/,
180  int const /*process_id*/)
181  {
182  }
183 
184 protected:
186  {
188  }
189 
191  {
193  }
194 
201  const NumLib::LocalToGlobalIndexMap& dof_table, const int process_id);
202 
203 private:
206  NumLib::LocalToGlobalIndexMap const& dof_table,
207  MeshLib::Mesh const& mesh,
208  unsigned const integration_order) = 0;
209 
212  virtual void initializeBoundaryConditions();
213 
215  std::vector<GlobalVector*>& /*x*/,
216  double const /*t*/,
217  int const /*process_id*/)
218  {
219  }
220 
221  virtual void preAssembleConcreteProcess(const double /*t*/,
222  double const /*dt*/,
223  GlobalVector const& /*x*/)
224  {
225  }
226 
227  virtual void assembleConcreteProcess(const double t, double const dt,
228  std::vector<GlobalVector*> const& x,
229  std::vector<GlobalVector*> const& xdot,
230  int const process_id, GlobalMatrix& M,
231  GlobalMatrix& K, GlobalVector& b) = 0;
232 
234  const double t, double const dt, std::vector<GlobalVector*> const& x,
235  std::vector<GlobalVector*> const& xdot, const double dxdot_dx,
236  const double dx_dx, int const process_id, GlobalMatrix& M,
237  GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac) = 0;
238 
240  std::vector<GlobalVector*> const& /*x*/,
241  const double /*t*/,
242  const double /*dt*/,
243  const int /*process_id*/)
244  {
245  }
246 
248  std::vector<GlobalVector*> const& /*x*/,
249  const double /*t*/,
250  const double /*dt*/,
251  int const /*process_id*/)
252  {
253  }
254 
256  GlobalVector const& /*x*/, GlobalVector const& /*xdot*/,
257  const double /*t*/, double const /*dt*/, int const /*process_id*/)
258  {
259  }
260 
261  virtual void preIterationConcreteProcess(const unsigned /*iter*/,
262  GlobalVector const& /*x*/)
263  {
264  }
265 
267  double const /*t*/,
268  double const /*dt*/,
269  std::vector<GlobalVector*> const& /*x*/,
270  GlobalVector const& /*x_dot*/,
271  int const /*process_id*/)
272  {
273  }
274 
276  GlobalVector const& /*x*/)
277  {
279  }
280 
281 protected:
288  virtual void constructDofTable();
289 
295 
302  const int specified_prosess_id);
303 
312  virtual std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
314 
315 private:
316  void initializeExtrapolator();
317 
320  void computeSparsityPattern();
321 
322 public:
323  std::string const name;
324 
325 protected:
327  std::unique_ptr<MeshLib::MeshSubset const> _mesh_subset_all_nodes;
328 
329  std::unique_ptr<NumLib::LocalToGlobalIndexMap> _local_to_global_index_map;
330 
332 
334 
336 
340 
344  unsigned const _integration_order;
345 
349  std::vector<std::unique_ptr<IntegrationPointWriter>>
351 
353 
354 protected:
359  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
361 
365  std::vector<BoundaryConditionCollection> _boundary_conditions;
366 
367 private:
371  std::vector<SourceTermCollection> _source_term_collections;
372 
374 };
375 
376 } // namespace ProcessLib
MathLib::SparsityPattern< GlobalIndexType > GlobalSparsityPattern
Global vector based on Eigen vector.
Definition: EigenVector.h:26
NumLib::LocalToGlobalIndexMap const & getDOFTable() const
NumLib::Extrapolator & getExtrapolator() const
std::string const name
Definition: Process.h:323
virtual void initializeBoundaryConditions()
Definition: Process.cpp:81
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions(double const t, GlobalVector const &x, int const process_id) const final
Definition: Process.h:131
virtual NumLib::IterationResult postIterationConcreteProcess(GlobalVector const &)
Definition: Process.h:275
virtual Eigen::Vector3d getFlux(std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
Definition: Process.h:167
virtual void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int)
Definition: Process.h:239
std::vector< BoundaryConditionCollection > _boundary_conditions
Definition: Process.h:365
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_prosess_id)
Definition: Process.cpp:303
virtual void postTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, int const)
Definition: Process.h:247
void postNonLinearSolver(GlobalVector const &x, GlobalVector const &xdot, const double t, double const dt, int const process_id)
Definition: Process.cpp:402
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
virtual void preIterationConcreteProcess(const unsigned, GlobalVector const &)
Definition: Process.h:261
virtual void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)=0
void assembleWithJacobian(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
Definition: Process.cpp:230
void preAssemble(const double t, double const dt, GlobalVector const &x) final
Definition: Process.cpp:200
virtual void solveReactionEquation(std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
Definition: Process.h:176
virtual void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order)=0
Process specific initialization called by initialize().
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
Definition: Process.cpp:67
void computeSecondaryVariable(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
compute secondary variables for the coupled equations or for output.
Definition: Process.cpp:411
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable() const
Definition: Process.h:190
virtual void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &, GlobalVector const &, int const)
Definition: Process.h:266
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:327
MeshLib::Mesh & _mesh
Definition: Process.h:326
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition: Process.h:360
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const
Definition: Process.h:137
virtual void postNonLinearSolverConcreteProcess(GlobalVector const &, GlobalVector const &, const double, double const, int const)
Definition: Process.h:255
virtual void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &, double const, int const)
Definition: Process.h:214
NumLib::IterationResult postIteration(GlobalVector const &x) final
Definition: Process.cpp:430
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Definition: Process.h:350
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
CoupledSolutionsForStaggeredScheme * _coupled_solutions
Definition: Process.h:339
void setCoupledSolutionsForStaggeredScheme(CoupledSolutionsForStaggeredScheme *const coupled_solutions)
Definition: Process.h:94
void preTimestep(std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
Preprocessing before starting assembly for new timestep.
Definition: Process.cpp:382
void initializeExtrapolator()
Definition: Process.cpp:359
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
MeshLib::Mesh & getMesh() const
Definition: Process.h:143
virtual void constructDofTable()
Definition: Process.cpp:256
void constructMonolithicProcessDofTable()
Definition: Process.cpp:270
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
virtual void preAssembleConcreteProcess(const double, double const, GlobalVector const &)
Definition: Process.h:221
virtual void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)=0
unsigned const _integration_order
Definition: Process.h:344
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable >>> &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition: Process.cpp:22
std::vector< SourceTermCollection > _source_term_collections
Definition: Process.h:371
virtual void extrapolateIntegrationPointValuesToNodes(const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
Definition: Process.h:108
void updateDeactivatedSubdomains(double const time, const int process_id)
Definition: Process.cpp:190
std::vector< std::unique_ptr< IntegrationPointWriter > > const * getIntegrationPointWriter(MeshLib::Mesh const &mesh) const
Definition: Process.h:156
SecondaryVariableCollection const & getSecondaryVariables() const
Definition: Process.h:150
void computeSparsityPattern()
Definition: Process.cpp:376
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition: Process.cpp:182
void setInitialConditions(std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
Definition: Process.cpp:114
bool isMonolithicSchemeUsed() const
Definition: Process.h:102
void assemble(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
Definition: Process.cpp:206
void postTimestep(std::vector< GlobalVector * > const &x, const double t, const double delta_t, int const process_id)
Postprocessing after a complete timestep.
Definition: Process.cpp:392
ExtrapolatorData _extrapolator_data
Definition: Process.h:373
GlobalSparsityPattern _sparsity_pattern
Definition: Process.h:352
virtual void setCoupledTermForTheStaggeredSchemeToLocalAssemblers(int const)
Definition: Process.h:103
void preIteration(const unsigned iter, GlobalVector const &x) final
Definition: Process.cpp:424
const bool _use_monolithic_scheme
Definition: Process.h:335
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const
Definition: Process.cpp:335
Handles configuration of several secondary variables from the project file.
IterationResult
Status flags telling the NonlinearSolver if an iteration succeeded.