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_b_data,
66 std::vector<double>& local_Jac_data);
67
69 double const t, double const dt, Eigen::VectorXd const& local_x,
70 Eigen::VectorXd const& local_x_prev, int const process_id,
71 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
72
73 virtual void computeSecondaryVariable(
74 std::size_t const mesh_item_id,
75 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
76 double const t, double const dt, std::vector<GlobalVector*> const& x,
77 GlobalVector const& x_prev, int const process_id);
78
79 virtual void preTimestep(std::size_t const mesh_item_id,
80 NumLib::LocalToGlobalIndexMap const& dof_table,
81 GlobalVector const& x, double const t,
82 double const delta_t);
83
84 virtual void postTimestep(
85 std::size_t const mesh_item_id,
86 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
87 std::vector<GlobalVector*> const& x,
88 std::vector<GlobalVector*> const& x_prev, double const t,
89 double const dt, int const process_id);
90
92 std::size_t const mesh_item_id,
93 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
94 std::vector<GlobalVector*> const& x,
95 std::vector<GlobalVector*> const& x_prev, double const t,
96 double const dt, int const process_id);
97
101 virtual Eigen::Vector3d getFlux(
102 MathLib::Point3d const& /*p_local_coords*/,
103 double const /*t*/,
104 std::vector<double> const& /*local_x*/) const
105 {
106 return Eigen::Vector3d{};
107 }
108
110 virtual Eigen::Vector3d getFlux(
111 MathLib::Point3d const& /*p_local_coords*/,
112 double const /*t*/,
113 std::vector<std::vector<double>> const& /*local_xs*/) const
114 {
115 return Eigen::Vector3d{};
116 }
117
118private:
119 virtual void setInitialConditionsConcrete(Eigen::VectorXd const /*local_x*/,
120 double const /*t*/,
121 int const /*process_id*/)
122 {
123 }
124
125 virtual void initializeConcrete() {}
126 virtual void preTimestepConcrete(std::vector<double> const& /*local_x*/,
127 double const /*t*/, double const /*dt*/)
128 {
129 }
130
131 virtual void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
132 Eigen::VectorXd const& /*local_x_prev*/,
133 double const /*t*/, double const /*dt*/,
134 int const /*process_id*/)
135 {
136 }
137
139 Eigen::VectorXd const& /*local_x*/,
140 Eigen::VectorXd const& /*local_x_prev*/, double const /*t*/,
141 double const /*dt*/, int const /*process_id*/)
142 {
143 }
144
146 double const /*t*/,
147 double const /*dt*/,
148 Eigen::VectorXd const& /*local_x*/,
149 Eigen::VectorXd const& /*local_x_prev*/)
150 {
151 }
152};
153
154} // namespace ProcessLib
Definition of the Point3d class.
Global vector based on Eigen vector.
Definition EigenVector.h:25
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 assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, 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 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_b_data, std::vector< double > &local_Jac_data)
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)