Reaktoro
A unified framework for modeling chemically reactive systems
Phreeqc.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 // Reaktoro includes
21 #include <Reaktoro/Interfaces/Interface.hpp>
22 
23 // Forward declarations
24 class PHREEQC;
25 
26 namespace Reaktoro {
27 
28 // Forward declarations
29 class ReactionEquation;
30 
31 class Phreeqc : public Interface
32 {
33 public:
35  Phreeqc();
36 
39  Phreeqc(std::string database);
40 
42  virtual ~Phreeqc();
43 
45  virtual auto temperature() const -> double;
46 
48  virtual auto pressure() const -> double;
49 
51  virtual auto speciesAmounts() const -> Vector;
52 
54  virtual auto numElements() const -> unsigned;
55 
57  virtual auto numSpecies() const -> unsigned;
58 
60  virtual auto numPhases() const -> unsigned;
61 
63  virtual auto numSpeciesInPhase(Index iphase) const -> unsigned;
64 
66  virtual auto elementName(Index ielement) const -> std::string;
67 
69  virtual auto elementMolarMass(Index ielement) const -> double;
70 
72  virtual auto elementStoichiometry(Index ispecies, Index ielement) const -> double;
73 
75  virtual auto speciesName(Index ispecies) const -> std::string;
76 
78  virtual auto phaseName(Index iphase) const -> std::string;
79 
84  virtual auto properties(ThermoModelResult& res, double T, double P) -> void;
85 
90  virtual auto properties(ChemicalModelResult& res, double T, double P, VectorConstRef n) -> void;
91 
93  virtual auto clone() const -> std::shared_ptr<Interface>;
94 
101  auto set(double T, double P) -> void;
102 
111  auto set(double T, double P, VectorConstRef n) -> void;
112 
119  auto load(std::string database) -> void;
120 
128  auto execute(std::string input, std::string output) -> void;
129 
132  auto execute(std::string input) -> void;
133 
138  auto reset() -> void;
139 
141  auto reactions() const -> const std::vector<ReactionEquation>&;
142 
144  auto stoichiometricMatrix() const -> Matrix;
145 
147  auto standardMolarGibbsEnergies() const -> Vector;
148 
150  auto standardMolarEnthalpies() const -> Vector;
151 
153  auto standardMolarVolumes() const -> Vector;
154 
156  auto standardMolarHeatCapacitiesConstP() const -> Vector;
157 
159  auto standardMolarHeatCapacitiesConstV() const -> Vector;
160 
162  auto lnActivityCoefficients() const -> Vector;
163 
165  auto lnActivityConstants() const -> Vector;
166 
168  auto lnActivities() const -> Vector;
169 
171  auto lnEquilibriumConstants() const -> Vector;
172 
174  auto phaseMolarVolumes() const -> Vector;
175 
177  auto phreeqc() -> PHREEQC&;
178 
180  auto phreeqc() const -> const PHREEQC&;
181 
182 private:
183  struct Impl;
184 
185  std::shared_ptr<Impl> pimpl;
186 };
187 
188 } // namespace Reaktoro
auto stoichiometricMatrix() const -> Matrix
Return the stoichiometric matrix of the system of reactions.
Definition: Phreeqc.cpp:1246
virtual auto elementStoichiometry(Index ispecies, Index ielement) const -> double
Return the stoichiometry of an element in a species.
Definition: Phreeqc.cpp:1195
virtual auto numPhases() const -> unsigned
Return the number of phases.
Definition: Phreeqc.cpp:1160
virtual auto elementName(Index ielement) const -> std::string
Return the name of an element.
Definition: Phreeqc.cpp:1185
auto phaseMolarVolumes() const -> Vector
Return the molar volumes of the phases.
Definition: Phreeqc.cpp:1296
The result of a chemical model function that calculates the chemical properties of species.
Definition: ChemicalModel.hpp:30
virtual auto phaseName(Index iphase) const -> std::string
Return the name of a phase.
Definition: Phreeqc.cpp:1205
virtual auto elementMolarMass(Index ielement) const -> double
Return the molar mass of an element (in units of kg/mol)
Definition: Phreeqc.cpp:1190
auto database(std::string name) -> std::string
Return the contents of a built-in database as a string.
Definition: DatabaseUtils.cpp:65
auto execute(std::string input, std::string output) -> void
Execute a PHREEQC input script either provided as a file or input string.
Definition: Phreeqc.cpp:1226
virtual auto numElements() const -> unsigned
Return the number of elements.
Definition: Phreeqc.cpp:1150
virtual auto temperature() const -> double
Return the temperature (in units of K)
Definition: Phreeqc.cpp:1135
Phreeqc()
Construct a default Phreeqc instance.
Definition: Phreeqc.cpp:1124
The result of a thermodynamic model function that calculates standard thermodynamic properties of spe...
Definition: ThermoModel.hpp:30
auto standardMolarHeatCapacitiesConstV() const -> Vector
Return the standard molar isochoric heat capacities of the species (in units of J/(mol*K))
Definition: Phreeqc.cpp:1271
virtual auto numSpecies() const -> unsigned
Return the number of species.
Definition: Phreeqc.cpp:1155
auto lnEquilibriumConstants() const -> Vector
Return the ln equilibrium constants of the reactions.
Definition: Phreeqc.cpp:1291
The namespace containing all components of the Reaktoro library.
Definition: ChemicalScalar.hpp:24
auto reactions() const -> const std::vector< ReactionEquation > &
Return the system of reactions.
Definition: Phreeqc.cpp:1241
virtual auto properties(ThermoModelResult &res, double T, double P) -> void
Return the thermodynamic properties of the phases and its species.
Definition: Phreeqc.cpp:1301
auto load(std::string database) -> void
Load a PHREEQC database.
Definition: Phreeqc.cpp:1220
virtual auto pressure() const -> double
Return the pressure (in units of Pa)
Definition: Phreeqc.cpp:1140
auto standardMolarVolumes() const -> Vector
Return the standard molar volumes of the species (in units of m3/mol)
Definition: Phreeqc.cpp:1261
auto lnActivityConstants() const -> Vector
Return the ln activity contants of the species.
Definition: Phreeqc.cpp:1281
A class used to interface other codes with Reaktoro.
Definition: Interface.hpp:37
Definition: Phreeqc.hpp:32
virtual ~Phreeqc()
Destroy this Phreeqc instance.
Definition: Phreeqc.cpp:1132
virtual auto speciesName(Index ispecies) const -> std::string
Return the name of a species.
Definition: Phreeqc.cpp:1200
std::size_t Index
Define a type that represents an index.
Definition: Index.hpp:26
auto standardMolarGibbsEnergies() const -> Vector
Return the standard molar Gibbs free energies of the species (in units of J/mol)
Definition: Phreeqc.cpp:1251
auto standardMolarHeatCapacitiesConstP() const -> Vector
Return the standard molar isobaric heat capacities of the species (in units of J/(mol*K))
Definition: Phreeqc.cpp:1266
auto set(double T, double P) -> void
Set the temperature and pressure of the interfaced code.
Definition: Phreeqc.cpp:1210
auto lnActivityCoefficients() const -> Vector
Return the ln activity coefficients of the species.
Definition: Phreeqc.cpp:1276
auto reset() -> void
Reset this Phreeqc instance to a clean state.
Definition: Phreeqc.cpp:1236
auto standardMolarEnthalpies() const -> Vector
Return the standard molar enthalpies of the species (in units of J/mol)
Definition: Phreeqc.cpp:1256
Eigen::MatrixXd Matrix
Alias to Eigen type MatrixXd.
Definition: Matrix.hpp:42
auto phreeqc() -> PHREEQC &
Return a reference to the low-level Phreeqc instance.
Definition: Phreeqc.cpp:1348
virtual auto clone() const -> std::shared_ptr< Interface >
Return a clone of this Phreeqc instance.
Definition: Phreeqc.cpp:1343
auto lnActivities() const -> Vector
Return the ln activities of the species.
Definition: Phreeqc.cpp:1286
Eigen::Ref< const Eigen::VectorXd > VectorConstRef
< Alias to Eigen type Ref<VectorXd>.
Definition: Matrix.hpp:31
virtual auto numSpeciesInPhase(Index iphase) const -> unsigned
Return the number of species in a phase.
Definition: Phreeqc.cpp:1165
virtual auto speciesAmounts() const -> Vector
Return the amounts of the species (in units of mol)
Definition: Phreeqc.cpp:1145