OGS
LocalAssemblerInterface.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14#include <optional>
15
16#include "MathLib/Point3d.h"
18
19namespace NumLib
20{
21class LocalToGlobalIndexMap;
22} // namespace NumLib
23
24namespace ProcessLib
25{
27{
29 int size;
30};
37{
38public:
39 virtual ~LocalAssemblerInterface() = default;
40
41 virtual void setInitialConditions(
42 std::size_t const mesh_item_id,
43 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
44 std::vector<GlobalVector*> const& x, double const t,
45 int const process_id);
46
47 virtual void initialize(std::size_t const mesh_item_id,
48 NumLib::LocalToGlobalIndexMap const& dof_table);
49
50 virtual void preAssemble(double const /*t*/, double const /*dt*/,
51 std::vector<double> const& /*local_x*/){};
52
53 virtual void assemble(double const t, double const dt,
54 std::vector<double> const& local_x,
55 std::vector<double> const& local_x_prev,
56 std::vector<double>& local_M_data,
57 std::vector<double>& local_K_data,
58 std::vector<double>& local_b_data);
59
60 virtual void assembleForStaggeredScheme(double const t, double const dt,
61 Eigen::VectorXd const& local_x,
62 Eigen::VectorXd const& local_x_prev,
63 int const process_id,
64 std::vector<double>& local_M_data,
65 std::vector<double>& local_K_data,
66 std::vector<double>& local_b_data);
67
68 virtual void assembleWithJacobian(double const t, double const dt,
69 std::vector<double> const& local_x,
70 std::vector<double> const& local_x_prev,
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_x_prev, int const process_id,
77 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
78
79 virtual void computeSecondaryVariable(
80 std::size_t const mesh_item_id,
81 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
82 double const t, double const dt, std::vector<GlobalVector*> const& x,
83 GlobalVector const& x_prev, int const process_id);
84
85 virtual void preTimestep(std::size_t const mesh_item_id,
86 NumLib::LocalToGlobalIndexMap const& dof_table,
87 GlobalVector const& x, double const t,
88 double const delta_t);
89
90 virtual void postTimestep(
91 std::size_t const mesh_item_id,
92 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
93 std::vector<GlobalVector*> const& x,
94 std::vector<GlobalVector*> const& x_prev, double const t,
95 double const dt, int const process_id);
96
98 std::size_t const mesh_item_id,
99 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
100 std::vector<GlobalVector*> const& x,
101 std::vector<GlobalVector*> const& x_prev, double const t,
102 double const dt, int const process_id);
103
107 virtual Eigen::Vector3d getFlux(
108 MathLib::Point3d const& /*p_local_coords*/,
109 double const /*t*/,
110 std::vector<double> const& /*local_x*/) const
111 {
112 return Eigen::Vector3d{};
113 }
114
116 virtual Eigen::Vector3d getFlux(
117 MathLib::Point3d const& /*p_local_coords*/,
118 double const /*t*/,
119 std::vector<std::vector<double>> const& /*local_xs*/) const
120 {
121 return Eigen::Vector3d{};
122 }
123
124 virtual std::optional<VectorSegment> getVectorDeformationSegment() const
125 {
126 return std::nullopt;
127 }
128
129private:
130 virtual void setInitialConditionsConcrete(Eigen::VectorXd const /*local_x*/,
131 double const /*t*/,
132 int const /*process_id*/)
133 {
134 }
135
136 virtual void initializeConcrete() {}
137 virtual void preTimestepConcrete(std::vector<double> const& /*local_x*/,
138 double const /*t*/, double const /*dt*/)
139 {
140 }
141
142 virtual void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
143 Eigen::VectorXd const& /*local_x_prev*/,
144 double const /*t*/, double const /*dt*/,
145 int const /*process_id*/)
146 {
147 }
148
150 Eigen::VectorXd const& /*local_x*/,
151 Eigen::VectorXd const& /*local_x_prev*/, double const /*t*/,
152 double const /*dt*/, int const /*process_id*/)
153 {
154 }
155
157 double const /*t*/,
158 double const /*dt*/,
159 Eigen::VectorXd const& /*local_x*/,
160 Eigen::VectorXd const& /*local_x_prev*/)
161 {
162 }
163};
164
165} // 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 std::optional< VectorSegment > getVectorDeformationSegment() const
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)