OGS
ChemicalSolverInterface.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Sparse>
14
18#include "MeshLib/Mesh.h"
19
21{
22class Medium;
23}
24
25namespace ParameterLib
26{
27class SpatialPosition;
28}
29
30namespace ChemistryLib
31{
33{
34public:
36 GlobalLinearSolver& linear_solver_)
37 : _mesh(mesh), linear_solver(linear_solver_)
38 {
39 }
40
41 std::vector<std::size_t> const& getElementIDs()
42 {
43 return *_mesh.getProperties().template getPropertyVector<std::size_t>(
44 "bulk_element_ids", MeshLib::MeshItemType::Cell, 1);
45 }
46
47 virtual void initialize() {}
48
50 std::vector<double> const& /*concentrations*/,
51 GlobalIndexType const& /*chemical_system_id*/,
52 MaterialPropertyLib::Medium const& /*medium*/,
53 ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/)
54 {
55 }
56
58 std::vector<double> const& /*concentrations*/,
59 GlobalIndexType const& /*chemical_system_id*/,
60 MaterialPropertyLib::Medium const* /*medium*/,
62 ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/,
63 double const /*dt*/)
64 {
65 }
66
68
69 virtual void executeSpeciationCalculation([[maybe_unused]] double const dt)
70 {
71 }
72
73 virtual double getConcentration(
74 int const /*component_id*/,
75 GlobalIndexType const /*chemical_system_id*/) const
76 {
77 return std::numeric_limits<double>::quiet_NaN();
78 }
79
80 virtual std::vector<std::string> const getComponentList() const
81 {
82 return {};
83 }
84
85 virtual Eigen::SparseMatrix<double> const* getStoichiometricMatrix() const
86 {
87 return nullptr;
88 }
89
90 virtual double getKineticPrefactor(
91 [[maybe_unused]] std::size_t reaction_id) const
92 {
93 return std::numeric_limits<double>::quiet_NaN();
94 }
95
97 GlobalIndexType const& /*chemical_system_id*/,
98 MaterialPropertyLib::Medium const& /*medium*/,
99 ParameterLib::SpatialPosition const& /*pos*/, double const /*porosity*/,
100 double const /*t*/, double const /*dt*/)
101 {
102 }
103
105 GlobalIndexType const& /*chemical_system_id*/,
106 MaterialPropertyLib::Medium const& /*medium*/,
107 double& /*porosity*/)
108 {
109 }
110
112 std::size_t const /*ele_id*/,
113 std::vector<GlobalIndexType> const& /*chemical_system_indices*/)
114 {
115 }
116
117 virtual ~ChemicalSolverInterface() = default;
118
119public:
120 std::vector<GlobalIndexType> chemical_system_index_map;
125};
126} // namespace ChemistryLib
GlobalMatrix::IndexType GlobalIndexType
Definition of the Mesh class.
virtual void updatePorosityPostReaction(GlobalIndexType const &, MaterialPropertyLib::Medium const &, double &)
virtual double getKineticPrefactor(std::size_t reaction_id) const
virtual void computeSecondaryVariable(std::size_t const, std::vector< GlobalIndexType > const &)
virtual void setChemicalSystemConcrete(std::vector< double > const &, GlobalIndexType const &, MaterialPropertyLib::Medium const *, MaterialPropertyLib::VariableArray const &, ParameterLib::SpatialPosition const &, double const, double const)
virtual void updateVolumeFractionPostReaction(GlobalIndexType const &, MaterialPropertyLib::Medium const &, ParameterLib::SpatialPosition const &, double const, double const, double const)
std::vector< GlobalIndexType > chemical_system_index_map
virtual ~ChemicalSolverInterface()=default
virtual void executeSpeciationCalculation(double const dt)
virtual std::vector< std::string > const getComponentList() const
virtual Eigen::SparseMatrix< double > const * getStoichiometricMatrix() const
ChemicalSolverInterface(MeshLib::Mesh const &mesh, GlobalLinearSolver &linear_solver_)
virtual double getConcentration(int const, GlobalIndexType const) const
virtual void initializeChemicalSystemConcrete(std::vector< double > const &, GlobalIndexType const &, MaterialPropertyLib::Medium const &, ParameterLib::SpatialPosition const &, double const)
std::vector< std::size_t > const & getElementIDs()
Properties & getProperties()
Definition Mesh.h:134