OGS
LocalAssemblerInterface.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <Eigen/Core>
7
8#include "MathLib/Point3d.h"
10
11namespace NumLib
12{
14} // namespace NumLib
15
16namespace ProcessLib
17{
24{
25public:
26 virtual ~LocalAssemblerInterface() = default;
27
28 virtual void setInitialConditions(
29 std::size_t const mesh_item_id,
30 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
31 std::vector<GlobalVector*> const& x, double const t,
32 int const process_id);
33
34 virtual void initialize(std::size_t const mesh_item_id,
35 NumLib::LocalToGlobalIndexMap const& dof_table);
36
37 virtual void preAssemble(double const /*t*/, double const /*dt*/,
38 std::vector<double> const& /*local_x*/) {};
39
40 virtual void assemble(double const t, double const dt,
41 std::vector<double> const& local_x,
42 std::vector<double> const& local_x_prev,
43 std::vector<double>& local_M_data,
44 std::vector<double>& local_K_data,
45 std::vector<double>& local_b_data);
46
47 virtual void assembleForStaggeredScheme(double const t, double const dt,
48 Eigen::VectorXd const& local_x,
49 Eigen::VectorXd const& local_x_prev,
50 int const process_id,
51 std::vector<double>& local_M_data,
52 std::vector<double>& local_K_data,
53 std::vector<double>& local_b_data);
54
55 virtual void assembleWithJacobian(double const t, double const dt,
56 std::vector<double> const& local_x,
57 std::vector<double> const& local_x_prev,
58 std::vector<double>& local_b_data,
59 std::vector<double>& local_Jac_data);
60
62 double const t, double const dt, Eigen::VectorXd const& local_x,
63 Eigen::VectorXd const& local_x_prev, int const process_id,
64 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
65
66 virtual void computeSecondaryVariable(
67 std::size_t const mesh_item_id,
68 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
69 double const t, double const dt, std::vector<GlobalVector*> const& x,
70 GlobalVector const& x_prev, int const process_id);
71
72 virtual void preTimestep(std::size_t const mesh_item_id,
73 NumLib::LocalToGlobalIndexMap const& dof_table,
74 GlobalVector const& x, double const t,
75 double const delta_t);
76
77 virtual void postTimestep(
78 std::size_t const mesh_item_id,
79 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
80 std::vector<GlobalVector*> const& x,
81 std::vector<GlobalVector*> const& x_prev, double const t,
82 double const dt, int const process_id);
83
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
94 virtual Eigen::Vector3d getFlux(
95 MathLib::Point3d const& /*p_local_coords*/,
96 double const /*t*/,
97 std::vector<double> const& /*local_x*/) const
98 {
99 return Eigen::Vector3d{};
100 }
101
103 virtual Eigen::Vector3d getFlux(
104 MathLib::Point3d const& /*p_local_coords*/,
105 double const /*t*/,
106 std::vector<std::vector<double>> const& /*local_xs*/) const
107 {
108 return Eigen::Vector3d{};
109 }
110
111 virtual int getNumberOfVectorElementsForDeformation() const { return 0; }
112
113private:
114 virtual void setInitialConditionsConcrete(Eigen::VectorXd const /*local_x*/,
115 double const /*t*/,
116 int const /*process_id*/)
117 {
118 }
119
120 virtual void initializeConcrete() {}
121 virtual void preTimestepConcrete(std::vector<double> const& /*local_x*/,
122 double const /*t*/, double const /*dt*/)
123 {
124 }
125
126 virtual void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
127 Eigen::VectorXd const& /*local_x_prev*/,
128 double const /*t*/, double const /*dt*/,
129 int const /*process_id*/)
130 {
131 }
132
134 Eigen::VectorXd const& /*local_x*/,
135 Eigen::VectorXd const& /*local_x_prev*/, double const /*t*/,
136 double const /*dt*/, int const /*process_id*/)
137 {
138 }
139
141 double const /*t*/,
142 double const /*dt*/,
143 Eigen::VectorXd const& /*local_x*/,
144 Eigen::VectorXd const& /*local_x_prev*/)
145 {
146 }
147};
148
149} // namespace ProcessLib
MathLib::EigenVector GlobalVector
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 int getNumberOfVectorElementsForDeformation() const
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)