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 auto const* const bulk_element_ids = bulkElementIDs(_mesh);
40 if (bulk_element_ids == nullptr)
41 {
43 "The 'bulk_element_ids' property does not exist on the mesh "
44 "{:s}.",
45 _mesh.getName());
46 }
47 active_element_ids_.assign(bulk_element_ids->begin(),
48 bulk_element_ids->end());
49 }
50
51 std::vector<std::size_t> const& activeElementIDs() const
52 {
54 }
55
56 virtual void initialize() {}
57
59 std::vector<double> const& /*concentrations*/,
60 GlobalIndexType const& /*chemical_system_id*/,
61 MaterialPropertyLib::Medium const& /*medium*/,
62 ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/)
63 {
64 }
65
67 std::vector<double> const& /*concentrations*/,
68 GlobalIndexType const& /*chemical_system_id*/,
69 MaterialPropertyLib::Medium const* /*medium*/,
71 ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/,
72 double const /*dt*/)
73 {
74 }
75
77
78 virtual void executeSpeciationCalculation([[maybe_unused]] double const dt)
79 {
80 }
81
82 virtual double getConcentration(
83 int const /*component_id*/,
84 GlobalIndexType const /*chemical_system_id*/) const
85 {
86 return std::numeric_limits<double>::quiet_NaN();
87 }
88
89 virtual std::vector<std::string> const getComponentList() const
90 {
91 return {};
92 }
93
94 virtual Eigen::SparseMatrix<double> const* getStoichiometricMatrix() const
95 {
96 return nullptr;
97 }
98
99 virtual double getKineticPrefactor(
100 [[maybe_unused]] std::size_t reaction_id) const
101 {
102 return std::numeric_limits<double>::quiet_NaN();
103 }
104
106 GlobalIndexType const& /*chemical_system_id*/,
107 MaterialPropertyLib::Medium const& /*medium*/,
108 ParameterLib::SpatialPosition const& /*pos*/, double const /*porosity*/,
109 double const /*t*/, double const /*dt*/)
110 {
111 }
112
114 GlobalIndexType const& /*chemical_system_id*/,
115 MaterialPropertyLib::Medium const& /*medium*/,
116 double& /*porosity*/)
117 {
118 }
119
121 std::size_t const /*ele_id*/,
122 std::vector<GlobalIndexType> const& /*chemical_system_indices*/)
123 {
124 }
125
126 virtual ~ChemicalSolverInterface() = default;
127
128public:
129 std::vector<GlobalIndexType> chemical_system_index_map;
134
135private:
136 std::vector<std::size_t> active_element_ids_;
137};
138} // namespace ChemistryLib
#define OGS_FATAL(...)
Definition Error.h:26
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
std::vector< std::size_t > const & activeElementIDs() const
virtual void initializeChemicalSystemConcrete(std::vector< double > const &, GlobalIndexType const &, MaterialPropertyLib::Medium const &, ParameterLib::SpatialPosition const &, double const)
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:105