OGS
LocalAssemblerInterface.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <cassert>
7
10
11namespace ProcessLib
12{
14 double const /*t*/,
15 double const /*dt*/,
16 std::vector<double> const& /*local_x*/,
17 std::vector<double> const& /*local_x_prev*/,
18 std::vector<double>& /*local_M_data*/,
19 std::vector<double>& /*local_K_data*/,
20 std::vector<double>& /*local_b_data*/)
21{
23 "The assemble() function is not implemented in the local assembler.");
24}
25
27 double const /*t*/, double const /*dt*/, Eigen::VectorXd const& /*local_x*/,
28 Eigen::VectorXd const& /*local_x_prev*/, int const /*process_id*/,
29 std::vector<double>& /*local_M_data*/,
30 std::vector<double>& /*local_K_data*/,
31 std::vector<double>& /*local_b_data*/)
32{
34 "The assembleForStaggeredScheme() function is not implemented in the "
35 "local assembler.");
36}
37
39 double const /*t*/, double const /*dt*/,
40 std::vector<double> const& /*local_x*/,
41 std::vector<double> const& /*local_x_prev*/,
42 std::vector<double>& /*local_b_data*/,
43 std::vector<double>& /*local_Jac_data*/)
44{
46 "The assembleWithJacobian() function is not implemented in the local "
47 "assembler.");
48}
49
51 double const /*t*/, double const /*dt*/, Eigen::VectorXd const& /*local_x*/,
52 Eigen::VectorXd const& /*local_x_prev*/, int const /*process_id*/,
53 std::vector<double>& /*local_b_data*/,
54 std::vector<double>& /*local_Jac_data*/)
55{
57 "The assembleWithJacobianForStaggeredScheme() function is not "
58 "implemented in the local assembler.");
59}
60
62 std::size_t const mesh_item_id,
63 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
64 double const t, double const dt, std::vector<GlobalVector*> const& x,
65 GlobalVector const& x_prev, int const process_id)
66{
67 auto const local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
68
69 // Todo: A more decent way is to directly pass x_prevs as done for x
70 auto const indices =
71 NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
72 auto const local_x_prev_vec = x_prev.get(indices);
73 auto const local_x_prev = MathLib::toVector(local_x_prev_vec);
74
75 computeSecondaryVariableConcrete(t, dt, local_x, local_x_prev);
76}
77
79 std::size_t const mesh_item_id,
80 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
81 std::vector<GlobalVector*> const& x, double const t, int const process_id)
82{
83 auto local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
84
85 setInitialConditionsConcrete(local_x, t, process_id);
86}
87
89 std::size_t const /*mesh_item_id*/,
90 NumLib::LocalToGlobalIndexMap const& /*dof_table*/)
91{
93}
94
96 std::size_t const mesh_item_id,
97 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
98 double const t, double const delta_t)
99{
100 auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
101 auto const local_x = x.get(indices);
102
103 preTimestepConcrete(local_x, t, delta_t);
104}
105
107 std::size_t const mesh_item_id,
108 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
109 std::vector<GlobalVector*> const& x,
110 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
111 int const process_id)
112{
113 auto const local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
114 auto const local_x_prev =
115 NumLib::getLocalX(mesh_item_id, dof_tables, x_prev);
116
117 postTimestepConcrete(local_x, local_x_prev, t, dt, process_id);
118}
119
121 std::size_t const mesh_item_id,
122 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
123 std::vector<GlobalVector*> const& x,
124 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
125 int const process_id)
126{
127 auto const local_x = NumLib::getLocalX(mesh_item_id, dof_tables, x);
128 auto const local_x_prev =
129 NumLib::getLocalX(mesh_item_id, dof_tables, x_prev);
130
131 postNonLinearSolverConcrete(local_x, local_x_prev, t, dt, process_id);
132}
133
134} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenVector GlobalVector
double get(IndexType rowId) const
get entry
Definition EigenVector.h:52
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 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 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)
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)