OGS
LocalAssemblerInterface.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14
17
18namespace ProcessLib
19{
21 double const /*t*/,
22 double const /*dt*/,
23 std::vector<double> const& /*local_x*/,
24 std::vector<double> const& /*local_x_prev*/,
25 std::vector<double>& /*local_M_data*/,
26 std::vector<double>& /*local_K_data*/,
27 std::vector<double>& /*local_b_data*/)
28{
30 "The assemble() function is not implemented in the local assembler.");
31}
32
34 double const /*t*/, double const /*dt*/, Eigen::VectorXd const& /*local_x*/,
35 Eigen::VectorXd const& /*local_x_prev*/, int const /*process_id*/,
36 std::vector<double>& /*local_M_data*/,
37 std::vector<double>& /*local_K_data*/,
38 std::vector<double>& /*local_b_data*/)
39{
41 "The assembleForStaggeredScheme() function is not implemented in the "
42 "local assembler.");
43}
44
46 double const /*t*/, double const /*dt*/,
47 std::vector<double> const& /*local_x*/,
48 std::vector<double> const& /*local_x_prev*/,
49 std::vector<double>& /*local_M_data*/,
50 std::vector<double>& /*local_K_data*/,
51 std::vector<double>& /*local_b_data*/,
52 std::vector<double>& /*local_Jac_data*/)
53{
55 "The assembleWithJacobian() function is not implemented in the local "
56 "assembler.");
57}
58
60 double const /*t*/, double const /*dt*/, Eigen::VectorXd const& /*local_x*/,
61 Eigen::VectorXd const& /*local_x_prev*/, int const /*process_id*/,
62 std::vector<double>& /*local_M_data*/,
63 std::vector<double>& /*local_K_data*/,
64 std::vector<double>& /*local_b_data*/,
65 std::vector<double>& /*local_Jac_data*/)
66{
68 "The assembleWithJacobianForStaggeredScheme() function is not "
69 "implemented in the local assembler.");
70}
71
73 std::size_t const mesh_item_id,
74 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
75 double const t, double const dt, std::vector<GlobalVector*> const& x,
76 GlobalVector const& x_prev, int const process_id)
77{
78 auto const local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
79
80 // Todo: A more decent way is to directly pass x_prevs as done for x
81 auto const indices =
82 NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
83 auto const local_x_prev_vec = x_prev.get(indices);
84 auto const local_x_prev = MathLib::toVector(local_x_prev_vec);
85
86 computeSecondaryVariableConcrete(t, dt, local_x, local_x_prev);
87}
88
90 std::size_t const mesh_item_id,
91 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
92 std::vector<GlobalVector*> const& x, double const t, int const process_id)
93{
94 auto local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
95
96 setInitialConditionsConcrete(local_x, t, process_id);
97}
98
100 std::size_t const /*mesh_item_id*/,
101 NumLib::LocalToGlobalIndexMap const& /*dof_table*/)
102{
104}
105
107 std::size_t const mesh_item_id,
108 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
109 double const t, double const delta_t)
110{
111 auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
112 auto const local_x = x.get(indices);
113
114 preTimestepConcrete(local_x, t, delta_t);
115}
116
118 std::size_t const mesh_item_id,
119 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
120 std::vector<GlobalVector*> const& x,
121 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
122 int const process_id)
123{
124 auto const local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
125 auto const local_x_prev =
126 NumLib::getLocalX(mesh_item_id, dof_tables, x_prev);
127
128 postTimestepConcrete(local_x, local_x_prev, t, dt, process_id);
129}
130
132 std::size_t const mesh_item_id,
133 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
134 std::vector<GlobalVector*> const& x,
135 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
136 int const process_id)
137{
138 auto const local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
139 auto const local_x_prev =
140 NumLib::getLocalX(mesh_item_id, dof_tables, x_prev);
141
142 postNonLinearSolverConcrete(local_x, local_x_prev, t, dt, process_id);
143}
144
145} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
Global vector based on Eigen vector.
Definition EigenVector.h:25
double get(IndexType rowId) const
get entry
Definition EigenVector.h:58
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 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 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 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)
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
Eigen::VectorXd getLocalX(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x)