Cantera  4.0.0a1
Loading...
Searching...
No Matches
SurfPhase.cpp
Go to the documentation of this file.
1/**
2 * @file SurfPhase.cpp
3 * Definitions for a simple thermodynamic model of a surface phase
4 * derived from ThermoPhase, assuming an ideal solution model
5 * (see @ref thermoprops and class
6 * @link Cantera::SurfPhase SurfPhase@endlink).
7 */
8
9// This file is part of Cantera. See License.txt in the top-level directory or
10// at https://cantera.org/license.txt for license and copyright information.
11
18
19namespace Cantera
20{
21
22SurfPhase::SurfPhase(const string& infile, const string& id_)
23{
24 setNDim(2);
25 initThermoFile(infile, id_);
26}
27
29{
30 if (m_n0 <= 0.0) {
31 return 0.0;
32 }
34 return mean_X(m_h0);
35}
36
38{
39 return enthalpy_mole();
40}
41
43{
45 double s = 0.0;
46 for (size_t k = 0; k < m_kk; k++) {
47 s += moleFraction(k) * (m_s0[k] -
48 GasConstant * log(std::max(concentration(k) * size(k)/m_n0, SmallNumber)));
49 }
50 return s;
51}
52
53double SurfPhase::cp_mole() const
54{
56 return mean_X(m_cp0);
57}
58
59double SurfPhase::cv_mole() const
60{
61 return cp_mole();
62}
63
64void SurfPhase::getPartialMolarEnthalpies(span<double> hbar) const
65{
66 getEnthalpy_RT(hbar);
67 for (size_t k = 0; k < m_kk; k++) {
68 hbar[k] *= RT();
69 }
70}
71
72void SurfPhase::getPartialMolarEntropies(span<double> sbar) const
73{
74 getEntropy_R(sbar);
76 for (size_t k = 0; k < m_kk; k++) {
77 sbar[k] -= log(std::max(m_work[k], SmallNumber)) - logStandardConc(k);
78 sbar[k] *= GasConstant;
79 }
80}
81
82void SurfPhase::getPartialMolarCp(span<double> cpbar) const
83{
84 getCp_R(cpbar);
85 for (size_t k = 0; k < m_kk; k++) {
86 cpbar[k] *= GasConstant;
87 }
88}
89
90// HKM 9/1/11 The partial molar volumes returned here are really partial molar areas.
91// Partial molar volumes for this phase should actually be equal to zero.
92void SurfPhase::getPartialMolarVolumes(span<double> vbar) const
93{
95}
96
97void SurfPhase::getStandardChemPotentials(span<double> mu0) const
98{
99 checkArraySize("SurfPhase::getStandardChemPotentials", mu0.size(), m_kk);
101 copy(m_mu0.begin(), m_mu0.end(), mu0.begin());
102}
103
104void SurfPhase::getChemPotentials(span<double> mu) const
105{
106 checkArraySize("SurfPhase::getChemPotentials", mu.size(), m_kk);
108 copy(m_mu0.begin(), m_mu0.end(), mu.begin());
110 for (size_t k = 0; k < m_kk; k++) {
111 mu[k] += RT() * (log(std::max(m_work[k], SmallNumber)) - logStandardConc(k));
112 }
113}
114
115void SurfPhase::getActivityConcentrations(span<double> c) const
116{
118}
119
121{
122 return m_n0/size(k);
123}
124
125double SurfPhase::logStandardConc(size_t k) const
126{
127 return m_logn0 - m_logsize[k];
128}
129
130void SurfPhase::getGibbs_RT(span<double> grt) const
131{
132 checkArraySize("SurfPhase::getGibbs_RT", grt.size(), m_kk);
134 scale(m_mu0.begin(), m_mu0.end(), grt.begin(), 1.0/RT());
135}
136
137void SurfPhase::getEnthalpy_RT(span<double> hrt) const
138{
139 checkArraySize("SurfPhase::getEnthalpy_RT", hrt.size(), m_kk);
141 scale(m_h0.begin(), m_h0.end(), hrt.begin(), 1.0/RT());
142}
143
144void SurfPhase::getEntropy_R(span<double> sr) const
145{
146 checkArraySize("SurfPhase::getEntropy_R", sr.size(), m_kk);
148 scale(m_s0.begin(), m_s0.end(), sr.begin(), 1.0/GasConstant);
149}
150
151void SurfPhase::getCp_R(span<double> cpr) const
152{
153 checkArraySize("SurfPhase::getCp_R", cpr.size(), m_kk);
155 scale(m_cp0.begin(), m_cp0.end(), cpr.begin(), 1.0/GasConstant);
156}
157
158void SurfPhase::getStandardVolumes(span<double> vol) const
159{
160 checkArraySize("SurfPhase::getStandardVolumes", vol.size(), m_kk);
161 for (size_t k = 0; k < m_kk; k++) {
162 vol[k] = 0.0;
163 }
164}
165
166void SurfPhase::getGibbs_RT_ref(span<double> grt) const
167{
168 getGibbs_RT(grt);
169}
170
171void SurfPhase::getEnthalpy_RT_ref(span<double> hrt) const
172{
173 getEnthalpy_RT(hrt);
174}
175
176void SurfPhase::getEntropy_R_ref(span<double> sr) const
177{
178 getEntropy_R(sr);
179}
180
181void SurfPhase::getCp_R_ref(span<double> cprt) const
182{
183 getCp_R(cprt);
184}
185
186bool SurfPhase::addSpecies(shared_ptr<Species> spec)
187{
188 bool added = ThermoPhase::addSpecies(spec);
189 if (added) {
190 m_h0.push_back(0.0);
191 m_s0.push_back(0.0);
192 m_cp0.push_back(0.0);
193 m_mu0.push_back(0.0);
194 m_work.push_back(0.0);
195 m_speciesSize.push_back(spec->size);
196 m_logsize.push_back(log(spec->size));
197 if (m_kk == 1) {
198 vector<double> cov{1.0};
199 setCoverages({cov});
200 }
201 }
202 return added;
203}
204
206{
207 if (n0 <= 0.0) {
208 throw CanteraError("SurfPhase::setSiteDensity",
209 "Site density must be positive. Got {}", n0);
210 }
211 m_n0 = n0;
212 m_logn0 = log(m_n0);
213 compositionChanged(); // trigger update of density
214}
215
216void SurfPhase::setCoverages(span<const double> theta)
217{
218 checkArraySize("SurfPhase::setCoverages", theta.size(), m_kk);
219 double sum = 0.0;
220 for (size_t k = 0; k < m_kk; k++) {
221 sum += theta[k] / size(k);
222 }
223 if (sum <= 0.0) {
224 throw CanteraError("SurfPhase::setCoverages",
225 "Sum of Coverage fractions is zero or negative");
226 }
227 for (size_t k = 0; k < m_kk; k++) {
228 m_work[k] = theta[k] / (sum * size(k));
229 }
231}
232
233void SurfPhase::setCoveragesNoNorm(span<const double> theta)
234{
235 checkArraySize("SurfPhase::setCoveragesNoNorm", theta.size(), m_kk);
236 double sum = 0.0;
237 double sum2 = 0.0;
238 for (size_t k = 0; k < m_kk; k++) {
239 sum += theta[k] / size(k);
240 sum2 += theta[k];
241 }
242 if (sum <= 0.0) {
243 throw CanteraError("SurfPhase::setCoverages",
244 "Sum of Coverage fractions is zero or negative");
245 }
246 for (size_t k = 0; k < m_kk; k++) {
247 m_work[k] = theta[k] * sum2 / (sum * size(k));
248 }
250}
251
252void SurfPhase::getCoverages(span<double> theta) const
253{
254 checkArraySize("SurfPhase::getCoverages", theta.size(), m_kk);
255 double sum_X = 0.0;
256 double sum_X_s = 0.0;
257 getMoleFractions(theta);
258 for (size_t k = 0; k < m_kk; k++) {
259 sum_X += theta[k];
260 sum_X_s += theta[k] * size(k);
261 }
262 for (size_t k = 0; k < m_kk; k++) {
263 theta[k] *= size(k) * sum_X / sum_X_s;
264 }
265}
266
267void SurfPhase::setCoveragesByName(const string& cov)
268{
270}
271
273{
274 vector<double> cv(m_kk, 0.0);
275 bool ifound = false;
276 for (size_t k = 0; k < m_kk; k++) {
277 double c = getValue(cov, speciesName(k), 0.0);
278 if (c > 0.0) {
279 ifound = true;
280 cv[k] = c;
281 }
282 }
283 if (!ifound) {
284 throw CanteraError("SurfPhase::setCoveragesByName",
285 "Input coverages are all zero or negative");
286 }
287 setCoverages(cv);
288}
289
290void SurfPhase::setState(const AnyMap& state) {
291 if (state.hasKey("coverages")) {
292 if (state["coverages"].is<string>()) {
293 setCoveragesByName(state["coverages"].asString());
294 } else {
295 setCoveragesByName(state["coverages"].asMap<double>());
296 }
297 }
299}
300
302{
305 double q = 0;
306 double sumX = 0;
307 for (size_t k = 0; k < m_kk; k++) {
308 q += m_work[k] * size(k);
309 // This term accounts for unnormalized coverages used in calculation of
310 // derivatives using finite differences.
311 sumX += m_work[k];
312 }
313 assignDensity(m_n0 * meanMolecularWeight() * sumX / q);
314}
315
316void SurfPhase::_updateThermo(bool force) const
317{
318 double tnow = temperature();
319 if (m_tlast != tnow || force) {
321 m_tlast = tnow;
322 for (size_t k = 0; k < m_kk; k++) {
323 m_h0[k] *= GasConstant * tnow;
324 m_s0[k] *= GasConstant;
325 m_cp0[k] *= GasConstant;
326 m_mu0[k] = m_h0[k] - tnow*m_s0[k];
327 }
328 m_tlast = tnow;
329 }
330}
331
333{
334 if (m_input.hasKey("site-density")) {
335 // Units are kmol/m^2 for surface phases or kmol/m for edge phases
336 setSiteDensity(m_input.convert("site-density",
337 Units(1.0, 0, -static_cast<double>(m_ndim), 0, 0, 0, 1)));
338 }
339}
340
341void SurfPhase::getParameters(AnyMap& phaseNode) const
342{
344 phaseNode["site-density"].setQuantity(
345 m_n0, Units(1.0, 0, -static_cast<double>(m_ndim), 0, 0, 0, 1));
346}
347
348EdgePhase::EdgePhase(const string& infile, const string& id_)
349{
350 setNDim(1);
351 initThermoFile(infile, id_);
352}
353
354}
Declarations for the EdgePhase ThermoPhase object, which models the interface between two surfaces (s...
Declaration for class Cantera::Species.
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase,...
Headers for the factory class that can create known ThermoPhase objects (see Thermodynamic Properties...
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
Definition AnyMap.cpp:1477
double convert(const string &key, const string &units) const
Convert the item stored by the given key to the units specified in units.
Definition AnyMap.cpp:1595
Base class for exceptions thrown by Cantera classes.
EdgePhase(const string &infile="", const string &id="")
Construct and initialize an EdgePhase directly from an input file.
virtual void update(double T, span< double > cp_R, span< double > h_RT, span< double > s_R) const
Compute the reference-state properties for all species.
void getMoleFractions(span< double > x) const
Get the species mole fraction vector.
Definition Phase.cpp:441
void assignDensity(const double density_)
Set the internally stored constant density (kg/m^3) of the phase.
Definition Phase.cpp:608
void setNDim(size_t ndim)
Set the number of spatial dimensions (1, 2, or 3).
Definition Phase.h:576
size_t m_kk
Number of species in the phase.
Definition Phase.h:875
size_t m_ndim
Dimensionality of the phase.
Definition Phase.h:879
double temperature() const
Temperature (K).
Definition Phase.h:585
double meanMolecularWeight() const
The mean molecular weight. Units: (kg/kmol)
Definition Phase.h:676
virtual double concentration(const size_t k) const
Concentration of species k.
Definition Phase.cpp:485
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:143
double mean_X(span< const double > Q) const
Evaluate the mole-fraction-weighted mean of an array Q.
Definition Phase.cpp:627
double moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition Phase.cpp:447
virtual void getConcentrations(span< double > c) const
Get the species concentrations (kmol/m^3).
Definition Phase.cpp:491
virtual void compositionChanged()
Apply changes to the state which are needed after the composition changes.
Definition Phase.cpp:925
virtual void setMoleFractions_NoNorm(span< const double > x)
Set the mole fractions to the specified values without normalizing.
Definition Phase.cpp:316
virtual void setMoleFractions(span< const double > x)
Set the mole fractions to the specified values.
Definition Phase.cpp:283
const vector< string > & speciesNames() const
Return a const reference to the vector of species names.
Definition Phase.cpp:151
SurfPhase(const string &infile="", const string &id="")
Construct and initialize a SurfPhase ThermoPhase object directly from an input file.
Definition SurfPhase.cpp:22
double enthalpy_mole() const override
Return the Molar Enthalpy. Units: J/kmol.
Definition SurfPhase.cpp:28
void setSiteDensity(double n0)
Set the site density of the surface phase (kmol m-2)
double logStandardConc(size_t k=0) const override
Natural logarithm of the standard concentration of the kth species.
void setState(const AnyMap &state) override
Set the state using an AnyMap containing any combination of properties supported by the thermodynamic...
void getCp_R(span< double > cpr) const override
Get the nondimensional Heat Capacities at constant pressure for the species standard states at the cu...
void getStandardChemPotentials(span< double > mu0) const override
Get the array of chemical potentials at unit activity for the species at their standard states at the...
Definition SurfPhase.cpp:97
void getEntropy_R_ref(span< double > er) const override
Returns the vector of nondimensional entropies of the reference state at the current temperature of t...
void getPartialMolarEnthalpies(span< double > hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
Definition SurfPhase.cpp:64
void setCoverages(span< const double > theta)
Set the surface site fractions to a specified state.
void getPartialMolarCp(span< double > cpbar) const override
Return an array of partial molar heat capacities for the species in the mixture.
Definition SurfPhase.cpp:82
vector< double > m_logsize
vector storing the log of the size of each species.
Definition SurfPhase.h:360
void getGibbs_RT(span< double > grt) const override
Get the nondimensional Gibbs functions for the species in their standard states at the current T and ...
void getCp_R_ref(span< double > cprt) const override
Returns the vector of nondimensional constant pressure heat capacities of the reference state at the ...
vector< double > m_work
Temporary work array.
Definition SurfPhase.h:353
void getParameters(AnyMap &phaseNode) const override
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
double m_n0
Surface site density (kmol m-2)
Definition SurfPhase.h:329
double size(size_t k) const
Returns the number of sites occupied by one molecule of species k.
Definition SurfPhase.h:237
vector< double > m_h0
Temporary storage for the reference state enthalpies.
Definition SurfPhase.h:341
void getEnthalpy_RT_ref(span< double > hrt) const override
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of ...
double cv_mole() const override
Molar heat capacity at constant volume and composition [J/kmol/K].
Definition SurfPhase.cpp:59
void getEnthalpy_RT(span< double > hrt) const override
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
vector< double > m_s0
Temporary storage for the reference state entropies.
Definition SurfPhase.h:344
vector< double > m_cp0
Temporary storage for the reference state heat capacities.
Definition SurfPhase.h:347
vector< double > m_speciesSize
Vector of species sizes (number of sites occupied). length m_kk.
Definition SurfPhase.h:332
void getEntropy_R(span< double > sr) const override
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
double intEnergy_mole() const override
Return the Molar Internal Energy. Units: J/kmol.
Definition SurfPhase.cpp:37
double entropy_mole() const override
Return the Molar Entropy. Units: J/kmol-K.
Definition SurfPhase.cpp:42
void _updateThermo(bool force=false) const
Update the species reference state thermodynamic functions.
void getGibbs_RT_ref(span< double > grt) const override
Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temper...
double cp_mole() const override
Molar heat capacity at constant pressure and composition [J/kmol/K].
Definition SurfPhase.cpp:53
void setCoveragesNoNorm(span< const double > theta)
Set the surface site fractions to a specified state.
void getPartialMolarVolumes(span< double > vbar) const override
Return an array of partial molar volumes for the species in the mixture.
Definition SurfPhase.cpp:92
void getStandardVolumes(span< double > vol) const override
Get the molar volumes of the species standard states at the current T and P of the solution.
void compositionChanged() override
Apply changes to the state which are needed after the composition changes.
void getCoverages(span< double > theta) const
Return a vector of surface coverages.
void getPartialMolarEntropies(span< double > sbar) const override
Returns an array of partial molar entropies of the species in the solution.
Definition SurfPhase.cpp:72
double standardConcentration(size_t k=0) const override
Return the standard concentration for the kth species.
vector< double > m_mu0
Temporary storage for the reference state Gibbs energies.
Definition SurfPhase.h:350
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
void setCoveragesByName(const string &cov)
Set the coverages from a string of colon-separated name:value pairs.
void getChemPotentials(span< double > mu) const override
Get the species chemical potentials. Units: J/kmol.
void getActivityConcentrations(span< double > c) const override
Return a vector of activity concentrations for each species.
double m_logn0
log of the surface site density
Definition SurfPhase.h:335
virtual void getParameters(AnyMap &phaseNode) const
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
virtual void setState(const AnyMap &state)
Set the state using an AnyMap containing any combination of properties supported by the thermodynamic...
double RT() const
Return the Gas Constant multiplied by the current temperature.
double m_tlast
last value of the temperature processed by reference state
void initThermoFile(const string &inputFile, const string &id)
Initialize a ThermoPhase object using an input file.
MultiSpeciesThermo m_spthermo
Pointer to the calculation manager for species reference-state thermodynamic properties.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
AnyMap m_input
Data supplied via setParameters.
A representation of the units associated with a dimensional quantity.
Definition Units.h:35
Composition parseCompString(const string &ss, const vector< string > &names)
Parse a composition string into a map consisting of individual key:composition pairs.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition utilities.h:118
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const double SmallNumber
smallest number to compare to zero.
Definition ct_defs.h:161
const U & getValue(const map< T, U > &m, const T &key, const U &default_val)
Const accessor for a value in a map.
Definition utilities.h:223
map< string, double > Composition
Map from string names to doubles.
Definition ct_defs.h:180
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...