OGS
LocalAssemblerInterface.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14
15#include "MathLib/Point3d.h"
17
18namespace NumLib
19{
20class LocalToGlobalIndexMap;
21} // namespace NumLib
22
23namespace ProcessLib
24{
31{
32public:
33 virtual ~LocalAssemblerInterface() = default;
34
35 virtual void setInitialConditions(
36 std::size_t const mesh_item_id,
37 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
38 std::vector<GlobalVector*> const& x, double const t,
39 int const process_id);
40
41 virtual void initialize(std::size_t const mesh_item_id,
42 NumLib::LocalToGlobalIndexMap const& dof_table);
43
44 virtual void preAssemble(double const /*t*/, double const /*dt*/,
45 std::vector<double> const& /*local_x*/){};
46
47 virtual void assemble(double const t, double const dt,
48 std::vector<double> const& local_x,
49 std::vector<double> const& local_x_prev,
50 std::vector<double>& local_M_data,
51 std::vector<double>& local_K_data,
52 std::vector<double>& local_b_data);
53
54 virtual void assembleForStaggeredScheme(double const t, double const dt,
55 Eigen::VectorXd const& local_x,
56 Eigen::VectorXd const& local_x_prev,
57 int const process_id,
58 std::vector<double>& local_M_data,
59 std::vector<double>& local_K_data,
60 std::vector<double>& local_b_data);
61
62 virtual void assembleWithJacobian(double const t, double const dt,
63 std::vector<double> const& local_x,
64 std::vector<double> const& local_x_prev,
65 std::vector<double>& local_M_data,
66 std::vector<double>& local_K_data,
67 std::vector<double>& local_b_data,
68 std::vector<double>& local_Jac_data);
69
71 double const t, double const dt, Eigen::VectorXd const& local_x,
72 Eigen::VectorXd const& local_x_prev, int const process_id,
73 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
74 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
75
76 virtual void computeSecondaryVariable(
77 std::size_t const mesh_item_id,
78 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
79 double const t, double const dt, std::vector<GlobalVector*> const& x,
80 GlobalVector const& x_prev, int const process_id);
81
82 virtual void preTimestep(std::size_t const mesh_item_id,
83 NumLib::LocalToGlobalIndexMap const& dof_table,
84 GlobalVector const& x, double const t,
85 double const delta_t);
86
87 virtual void postTimestep(
88 std::size_t const mesh_item_id,
89 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
90 std::vector<GlobalVector*> const& x,
91 std::vector<GlobalVector*> const& x_prev, double const t,
92 double const dt, int const process_id);
93
95 std::size_t const mesh_item_id,
96 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
97 std::vector<GlobalVector*> const& x,
98 std::vector<GlobalVector*> const& x_prev, double const t,
99 double const dt, int const process_id);
100
104 virtual Eigen::Vector3d getFlux(
105 MathLib::Point3d const& /*p_local_coords*/,
106 double const /*t*/,
107 std::vector<double> const& /*local_x*/) const
108 {
109 return Eigen::Vector3d{};
110 }
111
113 virtual Eigen::Vector3d getFlux(
114 MathLib::Point3d const& /*p_local_coords*/,
115 double const /*t*/,
116 std::vector<std::vector<double>> const& /*local_xs*/) const
117 {
118 return Eigen::Vector3d{};
119 }
120
121private:
122 virtual void setInitialConditionsConcrete(Eigen::VectorXd const /*local_x*/,
123 double const /*t*/,
124 int const /*process_id*/)
125 {
126 }
127
128 virtual void initializeConcrete() {}
129 virtual void preTimestepConcrete(std::vector<double> const& /*local_x*/,
130 double const /*t*/, double const /*dt*/)
131 {
132 }
133
134 virtual void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
135 Eigen::VectorXd const& /*local_x_prev*/,
136 double const /*t*/, double const /*dt*/,
137 int const /*process_id*/)
138 {
139 }
140
142 Eigen::VectorXd const& /*local_x*/,
143 Eigen::VectorXd const& /*local_x_prev*/, double const /*t*/,
144 double const /*dt*/, int const /*process_id*/)
145 {
146 }
147
149 double const /*t*/,
150 double const /*dt*/,
151 Eigen::VectorXd const& /*local_x*/,
152 Eigen::VectorXd const& /*local_x_prev*/)
153 {
154 }
155};
156
157} // namespace ProcessLib
Definition of the Point3d class.
Global vector based on Eigen vector.
Definition EigenVector.h:25
virtual void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, 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 preTimestepConcrete(std::vector< double > const &, double const, double const)
virtual ~LocalAssemblerInterface()=default
virtual void setInitialConditionsConcrete(Eigen::VectorXd const, double const, int const)
virtual void postNonLinearSolverConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const)
void postNonLinearSolver(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
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 Eigen::Vector3d getFlux(MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
Fits to staggered scheme.
virtual void assembleWithJacobianForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, 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 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_prev, int const process_id)
virtual void preAssemble(double const, double const, std::vector< double > const &)
virtual void assembleForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const)
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)
virtual void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void setInitialConditions(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)