Reaktoro
A unified framework for modeling chemically reactive systems
Interface.hpp
1 // Reaktoro is a unified framework for modeling chemically reactive systems.
2 //
3 // Copyright (C) 2014-2018 Allan Leal
4 //
5 // This library is free software; you can redistribute it and/or
6 // modify it under the terms of the GNU Lesser General Public
7 // License as published by the Free Software Foundation; either
8 // version 2.1 of the License, or (at your option) any later version.
9 //
10 // This library is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 // Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with this library. If not, see <http://www.gnu.org/licenses/>.
17 
18 #pragma once
19 
20 // C++ includes
21 #include <string>
22 #include <memory>
23 
24 // Reaktoro includes
25 #include <Reaktoro/Math/Matrix.hpp>
26 #include <Reaktoro/Thermodynamics/Models/ChemicalModel.hpp>
27 #include <Reaktoro/Thermodynamics/Models/ThermoModel.hpp>
28 
29 namespace Reaktoro {
30 
31 // Forward declarations
32 class ChemicalState;
33 class ChemicalSystem;
34 
36 class Interface
37 {
38 public:
40  virtual ~Interface() = 0;
41 
43  virtual auto temperature() const -> double = 0;
44 
46  virtual auto pressure() const -> double = 0;
47 
49  virtual auto speciesAmounts() const -> Vector = 0;
50 
52  virtual auto numElements() const -> unsigned = 0;
53 
55  virtual auto numSpecies() const -> unsigned = 0;
56 
58  virtual auto numPhases() const -> unsigned = 0;
59 
61  virtual auto numSpeciesInPhase(Index iphase) const -> unsigned = 0;
62 
64  virtual auto elementName(Index ielement) const -> std::string = 0;
65 
67  virtual auto elementMolarMass(Index ielement) const -> double = 0;
68 
70  virtual auto elementStoichiometry(Index ispecies, Index ielement) const -> double = 0;
71 
73  virtual auto speciesName(Index ispecies) const -> std::string = 0;
74 
76  virtual auto phaseName(Index iphase) const -> std::string = 0;
77 
82  virtual auto properties(ThermoModelResult& res, double T, double P) -> void = 0;
83 
88  virtual auto properties(ChemicalModelResult& res, double T, double P, VectorConstRef n) -> void = 0;
89 
91  virtual auto clone() const -> std::shared_ptr<Interface> = 0;
92 
94  auto formulaMatrix() const -> Matrix;
95 
97  auto indexElement(std::string element) const -> Index;
98 
100  auto indexSpecies(std::string species) const -> Index;
101 
103  auto indexPhase(std::string phase) const -> Index;
104 
106  auto indexPhaseWithSpecies(Index ispecies) const -> Index;
107 
109  auto indexFirstSpeciesInPhase(Index iphase) const -> Index;
110 
112  auto system() const -> ChemicalSystem;
113 
116  auto state(const ChemicalSystem& system) const -> ChemicalState;
117 
119  operator ChemicalSystem() const;
120 };
121 
122 } // namespace Reaktoro
virtual auto phaseName(Index iphase) const -> std::string=0
Return the name of a phase.
virtual auto elementStoichiometry(Index ispecies, Index ielement) const -> double=0
Return the stoichiometry of an element in a species.
virtual auto properties(ThermoModelResult &res, double T, double P) -> void=0
Return the thermodynamic properties of the phases and its species.
A class to represent a system and its attributes and properties.
Definition: ChemicalSystem.hpp:38
The result of a chemical model function that calculates the chemical properties of species.
Definition: ChemicalModel.hpp:30
Provides a computational representation of the state of a multiphase chemical system.
Definition: ChemicalState.hpp:61
auto indexSpecies(std::string species) const -> Index
Return the index of a species.
Definition: Interface.cpp:80
virtual auto speciesName(Index ispecies) const -> std::string=0
Return the name of a species.
virtual auto numElements() const -> unsigned=0
Return the number of elements.
virtual auto temperature() const -> double=0
Return the temperature (in units of K)
virtual auto pressure() const -> double=0
Return the pressure (in units of Pa)
virtual auto clone() const -> std::shared_ptr< Interface >=0
Return a clone of this Interface instance.
The result of a thermodynamic model function that calculates standard thermodynamic properties of spe...
Definition: ThermoModel.hpp:30
The namespace containing all components of the Reaktoro library.
Definition: ChemicalScalar.hpp:24
auto state(const ChemicalSystem &system) const -> ChemicalState
Return a ChemicalState instance created from an instance of a class derived from Interface.
Definition: Interface.cpp:191
A class used to interface other codes with Reaktoro.
Definition: Interface.hpp:37
virtual auto numSpecies() const -> unsigned=0
Return the number of species.
virtual auto elementMolarMass(Index ielement) const -> double=0
Return the molar mass of an element (in units of kg/mol)
auto formulaMatrix() const -> Matrix
Return the formula matrix of the species.
Definition: Interface.cpp:60
virtual auto numSpeciesInPhase(Index iphase) const -> unsigned=0
Return the number of species in a phase.
virtual auto numPhases() const -> unsigned=0
Return the number of phases.
auto indexFirstSpeciesInPhase(Index iphase) const -> Index
Return the index of the first species in a phase.
Definition: Interface.cpp:112
auto indexPhaseWithSpecies(Index ispecies) const -> Index
Return the index of the phase with a species.
Definition: Interface.cpp:98
auto indexPhase(std::string phase) const -> Index
Return the index of a phase.
Definition: Interface.cpp:89
std::size_t Index
Define a type that represents an index.
Definition: Index.hpp:26
auto system() const -> ChemicalSystem
Return a ChemicalSystem instance created from an instance of a class derived from Interface.
Definition: Interface.cpp:122
virtual auto speciesAmounts() const -> Vector=0
Return the amounts of the species (in units of mol)
virtual auto elementName(Index ielement) const -> std::string=0
Return the name of an element.
Eigen::MatrixXd Matrix
Alias to Eigen type MatrixXd.
Definition: Matrix.hpp:42
auto indexElement(std::string element) const -> Index
Return the index of an element.
Definition: Interface.cpp:71
virtual ~Interface()=0
Virtual destructor.
Definition: Interface.cpp:57
Eigen::Ref< const Eigen::VectorXd > VectorConstRef
< Alias to Eigen type Ref<VectorXd>.
Definition: Matrix.hpp:31