OGS
StokesFlowProcess.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include "StokesFlowFEM.h"
14 #include "StokesFlowProcessData.h"
15 
16 #include "ProcessLib/Process.h"
17 
18 namespace ProcessLib
19 {
20 namespace StokesFlow
21 {
25 template <int GlobalDim>
26 class StokesFlowProcess final : public Process
27 {
28 public:
30  std::string name,
31  MeshLib::Mesh& mesh,
32  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&&
33  jacobian_assembler,
34  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const&
35  parameters,
36  unsigned const integration_order,
37  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
38  process_variables,
39  StokesFlowProcessData&& process_data,
40  SecondaryVariableCollection&& secondary_variables,
41  bool const use_monolithic_scheme);
42 
43  bool isLinear() const override { return false; }
44 
46  const int /*process_id*/) const override;
47 
49  const int /*process_id*/) const override;
50 
51  void computeSecondaryVariableConcrete(double const /*t*/,
52  double const /*dt*/,
53  std::vector<GlobalVector*> const& x,
54  GlobalVector const& /*x_dot*/,
55  int const /*process_id*/) override;
56 
57  void postTimestepConcreteProcess(std::vector<GlobalVector*> const& x,
58  const double t,
59  const double dt,
60  int const process_id) override;
61 
62 private:
63  void constructDofTable() override;
64 
66  NumLib::LocalToGlobalIndexMap const& dof_table,
67  MeshLib::Mesh const& mesh,
68  unsigned const integration_order) override;
69 
70  void initializeBoundaryConditions() override;
71 
72  void assembleConcreteProcess(const double t, double const dt,
73  std::vector<GlobalVector*> const& x,
74  std::vector<GlobalVector*> const& xdot,
75  int const process_id, GlobalMatrix& M,
76  GlobalMatrix& K, GlobalVector& b) override;
77 
79  const double t, double const dt, std::vector<GlobalVector*> const& x,
80  std::vector<GlobalVector*> const& xdot, int const process_id,
82  GlobalMatrix& Jac) override;
83 
84  std::vector<MeshLib::Node*> _base_nodes;
85  std::unique_ptr<MeshLib::MeshSubset const> _mesh_subset_base_nodes;
86 
88 
89  std::vector<std::unique_ptr<StokesFlowLocalAssemblerInterface>>
91 
92  std::unique_ptr<NumLib::LocalToGlobalIndexMap>
94 
97  std::unique_ptr<NumLib::LocalToGlobalIndexMap>
99 };
100 
101 extern template class StokesFlowProcess<2>;
102 } // namespace StokesFlow
103 } // namespace ProcessLib
Global vector based on Eigen vector.
Definition: EigenVector.h:28
std::string const name
Definition: Process.h:328
Handles configuration of several secondary variables from the project file.
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
StokesFlowProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::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, StokesFlowProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
std::vector< MeshLib::Node * > _base_nodes
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &x, GlobalVector const &, int const) override
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
void assembleWithJacobianConcreteProcess(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, GlobalMatrix &Jac) override
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) override
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const override
std::vector< std::unique_ptr< StokesFlowLocalAssemblerInterface > > _local_assemblers
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, const double t, const double dt, int const process_id) override
MathLib::MatrixSpecifications getMatrixSpecifications(const int) const override
static const double t