20 SurfPhase::SurfPhase(doublereal n0):
47 "Couldn't find phase name in file:" + id_);
51 string model = th[
"model"];
52 if (model !=
"Surface" && model !=
"Edge") {
54 "thermo model attribute must be Surface or Edge");
67 string model = th[
"model"];
68 if (model !=
"Surface" && model !=
"Edge") {
70 "thermo model attribute must be Surface or Edge");
77 m_logn0(right.m_logn0),
78 m_press(right.m_press),
79 m_tlast(right.m_tlast)
126 for (
size_t k = 0; k <
m_kk; k++) {
134 for (
size_t k = 0; k <
m_kk; k++) {
142 for (
size_t k = 0; k <
m_kk; k++) {
165 for (
size_t k = 0; k <
m_kk; k++) {
191 "Bad value for number of parameter");
228 for (
size_t k = 0; k <
m_kk; k++) {
257 "Number of species is equal to zero");
268 for (
size_t k = 0; k <
m_kk; k++) {
277 "Bad value for parameter");
286 for (
size_t k = 0; k <
m_kk; k++) {
290 for (
size_t k = 0; k <
m_kk; k++) {
291 cout <<
"theta(" << k <<
") = " << theta[k] << endl;
294 "Sum of Coverage fractions is zero or negative");
296 for (
size_t k = 0; k <
m_kk; k++) {
308 for (
size_t k = 0; k <
m_kk; k++) {
321 for (
size_t k = 0; k <
m_kk; k++) {
333 for (
size_t k = 0; k < kk; k++) {
342 "Input coverages are all zero or negative");
350 if (
m_tlast != tnow || force) {
355 for (
size_t k = 0; k <
m_kk; k++) {
367 eosdata.
_require(
"model",
"Surface");
368 doublereal n =
getFloat(eosdata,
"site_density",
"toSI");
371 "missing or negative site density");
404 if (&right !=
this) {
419 doublereal n =
getFloat(eosdata,
"site_density",
"toSI");
422 "missing or negative site density");
std::map< std::string, doublereal > compositionMap
Map connecting a string name with a double.
void setCoveragesNoNorm(const doublereal *theta)
Set the surface site fractions to a specified state.
ThermoPhase * duplMyselfAsThermoPhase() const
Duplicator from a ThermoPhase object.
void _require(const std::string &a, const std::string &v) const
Require that the current xml node have an attribute named by the first argument, a, and that this attribute have the the string value listed in the second argument, v.
XML_Node * get_XML_File(const std::string &file, int debug)
Return a pointer to the XML tree for a Cantera input file.
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase, assuming an ideal solution model (see Thermodynamic Properties and class SurfPhase).
virtual void getCp_R_ref(doublereal *cprt) const
Returns the vector of nondimensional constant pressure heat capacities of the reference state at the ...
const doublereal OneAtm
One atmosphere [Pa].
virtual void getGibbs_RT(doublereal *grt) const
Get the nondimensional Gibbs functions for the species in their standard states at the current T and ...
SurfPhase(doublereal n0=0.0)
Constructor.
ThermoPhase & operator=(const ThermoPhase &right)
Assignment operator.
void setCoveragesByName(const std::string &cov)
Set the coverages from a string of colon-separated name:value pairs.
Headers for the factory class that can create known ThermoPhase objects (see Thermodynamic Properties...
doublereal m_press
Current value of the pressure (Pa)
void setCoverages(const doublereal *theta)
Set the surface site fractions to a specified state.
virtual void getEntropy_R_ref(doublereal *er) const
Returns the vector of nondimensional entropies of the reference state at the current temperature of t...
virtual doublereal logStandardConc(size_t k=0) const
Return the log of the standard concentration for the kth species.
doublereal m_tlast
Current value of the temperature (Kelvin)
Class XML_Node is a tree-based representation of the contents of an XML file.
bool getOptionalFloat(const Cantera::XML_Node &parent, const std::string &name, doublereal &fltRtn, const std::string &type)
Get an optional floating-point value from a child element.
doublereal size(size_t k) const
This routine returns the size of species k.
virtual void getStandardVolumes(doublereal *vol) const
Get the molar volumes of the species standard states at the current T and P of the solution...
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
doublereal getFloat(const Cantera::XML_Node &parent, const std::string &name, const std::string &type)
Get a floating-point value from a child element.
virtual void getGibbs_RT_ref(doublereal *grt) const
Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temper...
vector_fp m_work
Temporary work array.
void getConcentrations(doublereal *const c) const
Get the species concentrations (kmol/m^3).
void getCoverages(doublereal *theta) const
Return a vector of surface coverages.
virtual void getActivityConcentrations(doublereal *c) const
Return a vector of activity concentrations for each species.
XML_Node & child(const size_t n) const
Return a changeable reference to the n'th child of the current node.
Base class for a phase with thermodynamic properties.
EdgePhase & operator=(const EdgePhase &right)
Assignment Operator.
doublereal mean_X(const doublereal *const Q) const
Evaluate the mole-fraction-weighted mean of an array Q.
doublereal m_n0
Surface site density (kmol m-2)
A thermodynamic phase representing a one dimensional edge between two surfaces.
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
bool importPhase(XML_Node &phase, ThermoPhase *th, SpeciesThermoFactory *spfactory)
Import a phase information into an empty thermophase object.
SurfPhase & operator=(const SurfPhase &right)
Assignment operator.
doublereal m_logn0
log of the surface site density
virtual void initThermo()
Initialize the SurfPhase object after all species have been set up.
virtual void setStateFromXML(const XML_Node &state)
Set the initial state of the Surface Phase from an XML_Node.
virtual void getCp_R(doublereal *cpr) const
Get the nondimensional Heat Capacities at constant pressure for the species standard states at the cu...
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
Base class for exceptions thrown by Cantera classes.
vector_fp m_logsize
vector storing the log of the size of each species.
virtual void getEnthalpy_RT(doublereal *hrt) const
Get the nondimensional Enthalpy functions for the species standard states at their standard states at...
compositionMap parseCompString(const std::string &ss, const std::vector< std::string > &names)
Parse a composition string into a map consisting of individual key:composition pairs.
virtual void getEnthalpy_RT_ref(doublereal *hrt) const
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of ...
ThermoPhase * duplMyselfAsThermoPhase() const
Duplicator from the ThermoPhase parent class.
virtual void setConcentrations(const doublereal *const conc)
Set the concentrations to the specified values within the phase.
virtual doublereal intEnergy_mole() const
Return the Molar Internal Energy. Units: J/kmol.
bool hasChild(const std::string &ch) const
Tests whether the current node has a child node with a particular name.
void _updateThermo(bool force=false) const
Update the species reference state thermodynamic functions.
vector_fp m_s0
Temporary storage for the reference state entropies.
virtual void getPartialMolarVolumes(doublereal *vbar) const
Return an array of partial molar volumes for the species in the mixture.
virtual void setParametersFromXML(const XML_Node &thermoData)
Set the Equation-of-State parameters by reading an XML Node Input.
size_t nSpecies() const
Returns the number of species in the phase.
void setNDim(size_t ndim)
Set the number of spatial dimensions (1, 2, or 3).
const std::vector< std::string > & speciesNames() const
Return a const reference to the vector of species names.
doublereal temperature() const
Temperature (K).
vector_fp m_h0
Temporary storage for the reference state enthalpies.
void setSiteDensity(doublereal n0)
Set the site density of the surface phase (kmol m-2)
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
virtual void setTemperature(const doublereal temp)
Set the internally stored temperature of the phase (K).
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
vector_fp m_mu0
Temporary storage for the reference state gibbs energies.
const doublereal GasConstant
Universal Gas Constant. [J/kmol/K].
Contains declarations for string manipulation functions within Cantera.
virtual doublereal standardConcentration(size_t k=0) const
Return the standard concentration for the kth species.
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
virtual void setParametersFromXML(const XML_Node &thermoData)
Set the Equation-of-State parameters by reading an XML Node Input.
EdgePhase(doublereal n0=0.0)
Constructor.
virtual void getEntropy_R(doublereal *sr) const
Get the array of nondimensional Entropy functions for the species standard states at the current T an...
size_t m_kk
Number of species in the phase.
virtual doublereal enthalpy_mole() const
Return the Molar Enthalpy. Units: J/kmol.
virtual void update(doublereal T, doublereal *cp_R, doublereal *h_RT, doublereal *s_R) const =0
Compute the reference-state properties for all species.
std::string getChildValue(const Cantera::XML_Node &parent, const std::string &nameString)
This function reads a child node with the name, nameString, and returns its xml value as the return s...
virtual void getPartialMolarEntropies(doublereal *sbar) const
Returns an array of partial molar entropies of the species in the solution.
virtual void getPartialMolarCp(doublereal *cpbar) const
Return an array of partial molar heat capacities for the species in the mixture.
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
SpeciesThermo * m_spthermo
Pointer to the calculation manager for species reference-state thermodynamic properties.
std::string speciesName(size_t k) const
Name of the species with index k.
XML_Node * get_XML_NameID(const std::string &nameTarget, const std::string &file_ID, XML_Node *root)
This routine will locate an XML node in either the input XML tree or in another input file specified ...
Declarations for the EdgePhase ThermoPhase object, which models the interface between two surfaces (s...
vector_fp m_cp0
Temporary storage for the reference state heat capacities.
virtual void getStandardChemPotentials(doublereal *mu0) const
Get the array of chemical potentials at unit activity for the standard state species at the current T...
virtual void setParameters(int n, doublereal *const c)
Set the equation of state parameters from the argument list.