Reaktoro  v2.11.0
A unified framework for modeling chemically reactive systems
Phreeqc.hpp
1 // Reaktoro is a unified framework for modeling chemically reactive systems.
2 //
3 // Copyright © 2014-2024 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 // Reaktoro includes
21 #include <Reaktoro/Common/Matrix.hpp>
22 #include <Reaktoro/Common/Types.hpp>
23 
24 // Forward declarations
25 class PHREEQC;
26 
27 namespace Reaktoro {
28 
29 // Forward declarations
30 class ReactionEquation;
31 
32 // Used to execute PHREEQC from C++/Python and inspect its internal state.
33 class Phreeqc
34 {
35 public:
38 
41  Phreeqc(String database);
42 
44  virtual ~Phreeqc();
45 
47  virtual auto temperature() const -> double;
48 
50  virtual auto pressure() const -> double;
51 
53  virtual auto speciesAmounts() const -> ArrayXd;
54 
56  virtual auto numElements() const -> unsigned;
57 
59  virtual auto numSpecies() const -> unsigned;
60 
62  virtual auto numPhases() const -> unsigned;
63 
65  virtual auto numSpeciesInPhase(Index iphase) const -> unsigned;
66 
68  virtual auto elementName(Index ielement) const -> String;
69 
71  virtual auto elementMolarMass(Index ielement) const -> double;
72 
74  virtual auto elementStoichiometry(Index ispecies, Index ielement) const -> double;
75 
77  virtual auto speciesName(Index ispecies) const -> String;
78 
80  virtual auto phaseName(Index iphase) const -> String;
81 
82  // /// Return the thermodynamic properties of the phases and its species.
83  // /// @param iphase The index of the phase
84  // /// @param T The temperature (in units of K)
85  // /// @param P The pressure (in units of Pa)
86  // virtual auto properties(ThermoModelResult& res, double T, double P) -> void;
87 
88  // /// Return the chemical properties of the phases and its species.
89  // /// @param T The temperature (in units of K)
90  // /// @param P The pressure (in units of Pa)
91  // /// @param n The amounts of the species (in units of mol)
92  // virtual auto properties(ChemicalModelResult& res, double T, double P, ArrayXdConstRef n) -> void;
93 
100  auto set(double T, double P) -> void;
101 
110  auto set(double T, double P, ArrayXdConstRef n) -> void;
111 
118  auto load(String database) -> void;
119 
127  auto execute(String input, String output) -> void;
128 
131  auto execute(String input) -> void;
132 
137  auto reset() -> void;
138 
140  auto reactions() const -> const Vec<Pairs<String, double>>&;
141 
144 
147 
150 
152  auto standardMolarVolumes() const -> ArrayXd;
153 
156 
159 
162 
164  auto lnActivityConstants() const -> ArrayXd;
165 
167  auto lnActivities() const -> ArrayXd;
168 
171 
173  auto phaseMolarVolumes() const -> ArrayXd;
174 
176  auto phreeqc() -> PHREEQC&;
177 
179  auto phreeqc() const -> const PHREEQC&;
180 
181 private:
182  struct Impl;
183 
184  std::shared_ptr<Impl> pimpl;
185 };
186 
187 } // namespace Reaktoro
Definition: Phreeqc.hpp:34
auto execute(String input, String output) -> void
Execute a PHREEQC input script either provided as a file or input string.
auto stoichiometricMatrix() const -> MatrixXd
Return the stoichiometric matrix of the system of reactions.
Phreeqc(String database)
Construct a Phreeqc instance with given PHREEQC database.
auto standardMolarGibbsEnergies() const -> ArrayXd
Return the standard molar Gibbs free energies of the species (in units of J/mol)
auto standardMolarVolumes() const -> ArrayXd
Return the standard molar volumes of the species (in units of m3/mol)
auto lnActivityCoefficients() const -> ArrayXd
Return the ln activity coefficients of the species.
auto load(String database) -> void
Load a PHREEQC database.
virtual auto numElements() const -> unsigned
Return the number of elements.
auto phreeqc() -> PHREEQC &
Return a reference to the low-level Phreeqc instance.
auto standardMolarHeatCapacitiesConstV() const -> ArrayXd
Return the standard molar isochoric heat capacities of the species (in units of J/(mol*K))
auto phaseMolarVolumes() const -> ArrayXd
Return the molar volumes of the phases.
virtual auto speciesName(Index ispecies) const -> String
Return the name of a species.
virtual auto numPhases() const -> unsigned
Return the number of phases.
virtual auto numSpeciesInPhase(Index iphase) const -> unsigned
Return the number of species in a phase.
auto lnEquilibriumConstants() const -> ArrayXd
Return the ln equilibrium constants of the reactions.
Phreeqc()
Construct a default Phreeqc instance.
auto set(double T, double P) -> void
Set the temperature and pressure of the interfaced code.
virtual auto elementMolarMass(Index ielement) const -> double
Return the molar mass of an element (in units of kg/mol)
virtual auto pressure() const -> double
Return the pressure (in units of Pa)
auto reset() -> void
Reset this Phreeqc instance to a clean state.
auto standardMolarEnthalpies() const -> ArrayXd
Return the standard molar enthalpies of the species (in units of J/mol)
virtual auto numSpecies() const -> unsigned
Return the number of species.
auto reactions() const -> const Vec< Pairs< String, double >> &
Return the system of reactions.
auto lnActivityConstants() const -> ArrayXd
Return the ln activity contants of the species.
virtual auto elementName(Index ielement) const -> String
Return the name of an element.
virtual auto phaseName(Index iphase) const -> String
Return the name of a phase.
auto standardMolarHeatCapacitiesConstP() const -> ArrayXd
Return the standard molar isobaric heat capacities of the species (in units of J/(mol*K))
virtual auto elementStoichiometry(Index ispecies, Index ielement) const -> double
Return the stoichiometry of an element in a species.
auto lnActivities() const -> ArrayXd
Return the ln activities of the species.
virtual auto speciesAmounts() const -> ArrayXd
Return the amounts of the species (in units of mol)
virtual ~Phreeqc()
Destroy this Phreeqc instance.
virtual auto temperature() const -> double
Return the temperature (in units of K)
The namespace containing all components of the Reaktoro library.
Definition: Algorithms.hpp:29
Vec< Pair< T, U > > Pairs
Convenient alias for std::vector<std::pair<T, U>>.
Definition: Types.hpp:90
std::vector< T > Vec
Convenient alias for std::vector<T>.
Definition: Types.hpp:66
Eigen::MatrixXd MatrixXd
Convenient alias to Eigen type.
Definition: Matrix.hpp:137
std::string String
Convenient alias for std::string.
Definition: Types.hpp:52
std::size_t Index
Define a type that represents an index.
Definition: Index.hpp:26
Eigen::Ref< const ArrayXd > ArrayXdConstRef
Convenient alias to Eigen type.
Definition: Matrix.hpp:105
Eigen::ArrayXd ArrayXd
Convenient alias to Eigen type.
Definition: Matrix.hpp:103