MargulesVPSSTP is a derived class of GibbsExcessVPSSTP that employs the Margules approximation for the excess Gibbs free energy. More...
#include <MargulesVPSSTP.h>
MargulesVPSSTP is a derived class of GibbsExcessVPSSTP that employs the Margules approximation for the excess Gibbs free energy.
MargulesVPSSTP derives from class GibbsExcessVPSSTP which is derived from VPStandardStateTP, and overloads the virtual methods defined there with ones that use expressions appropriate for the Margules Excess Gibbs free energy approximation.
The independent unknowns are pressure, temperature, and mass fraction.
All species are defined to have standard states that depend upon both the temperature and the pressure. The Margules approximation assumes symmetric standard states, where all of the standard state assume that the species are in pure component states at the temperature and pressure of the solution. I don't think it prevents, however, some species from being dilute in the solution.
The molar excess Gibbs free energy is given by the following formula which is a sum over interactions i. Each of the interactions are binary interactions involving two of the species in the phase, denoted, Ai and Bi. This is the generalization of the Margules formulation for a phase that has more than 2 species.
\[ G^E = \sum_i \left( H_{Ei} - T S_{Ei} \right) \]
\[ H^E_i = n X_{Ai} X_{Bi} \left( h_{o,i} + h_{1,i} X_{Bi} \right) \]
\[ S^E_i = n X_{Ai} X_{Bi} \left( s_{o,i} + s_{1,i} X_{Bi} \right) \]
where n is the total moles in the solution.
The activity of a species defined in the phase is given by an excess Gibbs free energy formulation.
\[ a_k = \gamma_k X_k \]
where
\[ R T \ln( \gamma_k )= \frac{d(n G^E)}{d(n_k)}\Bigg|_{n_i} \]
Taking the derivatives results in the following expression
\[ R T \ln( \gamma_k )= \sum_i \left( \left( \delta_{Ai,k} X_{Bi} + \delta_{Bi,k} X_{Ai} - X_{Ai} X_{Bi} \right) \left( g^E_{o,i} + g^E_{1,i} X_{Bi} \right) + \left( \delta_{Bi,k} - X_{Bi} \right) X_{Ai} X_{Bi} g^E_{1,i} \right) \]
where \( g^E_{o,i} = h_{o,i} - T s_{o,i} \) and \( g^E_{1,i} = h_{1,i} - T s_{1,i} \) and where \( X_k \) is the mole fraction of species k.
This object inherits from the class VPStandardStateTP. Therefore, the specification and calculation of all standard state and reference state values are handled at that level. Various functional forms for the standard state are permissible. The chemical potential for species k is equal to
\[ \mu_k(T,P) = \mu^o_k(T, P) + R T \ln(\gamma_k X_k) \]
The partial molar entropy for species k is given by
\[ \tilde{s}_k(T,P) = s^o_k(T,P) - R \ln( \gamma_k X_k ) - R T \frac{d \ln(\gamma_k) }{dT} \]
The partial molar enthalpy for species k is given by
\[ \tilde{h}_k(T,P) = h^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} \]
The partial molar volume for species k is
\[ \tilde V_k(T,P) = V^o_k(T,P) + R T \frac{d \ln(\gamma_k) }{dP} \]
The partial molar Heat Capacity for species k is
\[ \tilde{C}_{p,k}(T,P) = C^o_{p,k}(T,P) - 2 R T \frac{d \ln( \gamma_k )}{dT} - R T^2 \frac{d^2 \ln(\gamma_k) }{{dT}^2} \]
\( C^a_k \) are defined such that \( a_k = C^a_k / C^s_k, \) where \( C^s_k \) is a standard concentration defined below and \( a_k \) are activities used in the thermodynamic functions. These activity (or generalized) concentrations are used by kinetics manager classes to compute the forward and reverse rates of elementary reactions. The activity concentration, \( C^a_k \),is given by the following expression.
\[ C^a_k = C^s_k X_k = \frac{P}{R T} X_k \]
The standard concentration for species k is independent of k and equal to
\[ C^s_k = C^s = \frac{P}{R T} \]
For example, a bulk-phase binary gas reaction between species j and k, producing a new gas species l would have the following equation for its rate of progress variable, \( R^1 \), which has units of kmol m-3 s-1.
\[ R^1 = k^1 C_j^a C_k^a = k^1 (C^s a_j) (C^s a_k) \]
where
\[ C_j^a = C^s a_j \mbox{\quad and \quad} C_k^a = C^s a_k \]
\( C_j^a \) is the activity concentration of species j, and \( C_k^a \) is the activity concentration of species k. \( C^s \) is the standard concentration. \( a_j \) is the activity of species j which is equal to the mole fraction of j.
The reverse rate constant can then be obtained from the law of microscopic reversibility and the equilibrium expression for the system.
\[ \frac{a_j a_k}{ a_l} = K_a^{o,1} = \exp(\frac{\mu^o_l - \mu^o_j - \mu^o_k}{R T} ) \]
\( K_a^{o,1} \) is the dimensionless form of the equilibrium constant, associated with the pressure dependent standard states \( \mu^o_l(T,P) \) and their associated activities, \( a_l \), repeated here:
\[ \mu_l(T,P) = \mu^o_l(T, P) + R T \ln a_l \]
We can switch over to expressing the equilibrium constant in terms of the reference state chemical potentials
\[ K_a^{o,1} = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) * \frac{P_{ref}}{P} \]
The concentration equilibrium constant, \( K_c \), may be obtained by changing over to activity concentrations. When this is done:
\[ \frac{C^a_j C^a_k}{ C^a_l} = C^o K_a^{o,1} = K_c^1 = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) * \frac{P_{ref}}{RT} \]
Kinetics managers will calculate the concentration equilibrium constant, \( K_c \), using the second and third part of the above expression as a definition for the concentration equilibrium constant.
For completeness, the pressure equilibrium constant may be obtained as well
\[ \frac{P_j P_k}{ P_l P_{ref}} = K_p^1 = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) \]
\( K_p \) is the simplest form of the equilibrium constant for ideal gases. However, it isn't necessarily the simplest form of the equilibrium constant for other types of phases; \( K_c \) is used instead because it is completely general.
The reverse rate of progress may be written down as
\[ R^{-1} = k^{-1} C_l^a = k^{-1} (C^o a_l) \]
where we can use the concept of microscopic reversibility to write the reverse rate constant in terms of the forward rate constant and the concentration equilibrium constant, \( K_c \).
\[ k^{-1} = k^1 K^1_c \]
\( k^{-1} \) has units of s-1.
Definition at line 214 of file MargulesVPSSTP.h.
Public Member Functions | |
MargulesVPSSTP (const string &inputFile="", const string &id="") | |
Construct a MargulesVPSSTP object from an input file. | |
string | type () const override |
String indicating the thermodynamic model implemented. | |
Molar Thermodynamic Properties | |
double | enthalpy_mole () const override |
Molar enthalpy. Units: J/kmol. | |
double | entropy_mole () const override |
Molar entropy. Units: J/kmol/K. | |
double | cp_mole () const override |
Molar heat capacity at constant pressure. Units: J/kmol/K. | |
double | cv_mole () const override |
Molar heat capacity at constant volume. Units: J/kmol/K. | |
Activities, Standard States, and Activity Concentrations | |
The activity \( a_k \) of a species in solution is related to the chemical potential by \[ \mu_k = \mu_k^0(T) + \hat R T \ln a_k. \] The quantity \( \mu_k^0(T,P) \) is the chemical potential at unit activity, which depends only on temperature and pressure. | |
void | getLnActivityCoefficients (double *lnac) const override |
Get the array of non-dimensional molar-based ln activity coefficients at the current solution temperature, pressure, and solution concentration. | |
Partial Molar Properties of the Solution | |
void | getChemPotentials (double *mu) const override |
Get the species chemical potentials. Units: J/kmol. | |
void | getPartialMolarEnthalpies (double *hbar) const override |
Returns an array of partial molar enthalpies for the species in the mixture. | |
void | getPartialMolarEntropies (double *sbar) const override |
Returns an array of partial molar entropies for the species in the mixture. | |
void | getPartialMolarCp (double *cpbar) const override |
Returns an array of partial molar entropies for the species in the mixture. | |
void | getPartialMolarVolumes (double *vbar) const override |
Return an array of partial molar volumes for the species in the mixture. | |
void | getd2lnActCoeffdT2 (double *d2lnActCoeffdT2) const |
Get the array of temperature second derivatives of the log activity coefficients. | |
void | getdlnActCoeffdT (double *dlnActCoeffdT) const override |
Get the array of temperature derivatives of the log activity coefficients. | |
Initialization | |
The following methods are used in the process of constructing the phase and setting its parameters from a specification in an input file. They are not normally used in application programs. To see how they are used, see importPhase() | |
void | initThermo () override |
Initialize the ThermoPhase object after all species have been set up. | |
void | getParameters (AnyMap &phaseNode) const override |
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function. | |
void | addBinaryInteraction (const string &speciesA, const string &speciesB, double h0, double h1, double s0, double s1, double vh0, double vh1, double vs0, double vs1) |
Add a binary species interaction with the specified parameters. | |
Derivatives of Thermodynamic Variables needed for Applications | |
void | getdlnActCoeffds (const double dTds, const double *const dXds, double *dlnActCoeffds) const override |
Get the change in activity coefficients wrt changes in state (temp, mole fraction, etc) along a line in parameter space or along a line in physical space. | |
void | getdlnActCoeffdlnX_diag (double *dlnActCoeffdlnX_diag) const override |
Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only. | |
void | getdlnActCoeffdlnN_diag (double *dlnActCoeffdlnN_diag) const override |
Get the array of log species mole number derivatives of the log activity coefficients. | |
void | getdlnActCoeffdlnN (const size_t ld, double *const dlnActCoeffdlnN) override |
Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers. | |
Public Member Functions inherited from GibbsExcessVPSSTP | |
bool | addSpecies (shared_ptr< Species > spec) override |
Add a Species to this Phase. | |
Units | standardConcentrationUnits () const override |
Returns the units of the "standard concentration" for this phase. | |
void | getActivityConcentrations (double *c) const override |
This method returns an array of generalized concentrations. | |
double | standardConcentration (size_t k=0) const override |
The standard concentration \( C^0_k \) used to normalize the generalized concentration. | |
double | logStandardConc (size_t k=0) const override |
Natural logarithm of the standard concentration of the kth species. | |
void | getActivities (double *ac) const override |
Get the array of non-dimensional activities (molality based for this class and classes that derive from it) at the current solution temperature, pressure, and solution concentration. | |
void | getActivityCoefficients (double *ac) const override |
Get the array of non-dimensional molar-based activity coefficients at the current solution temperature, pressure, and solution concentration. | |
virtual void | getdlnActCoeffdlnX (double *dlnActCoeffdlnX) const |
Get the array of log concentration-like derivatives of the log activity coefficients. | |
Public Member Functions inherited from VPStandardStateTP | |
void | setTemperature (const double temp) override |
Set the temperature of the phase. | |
void | setPressure (double p) override |
Set the internally stored pressure (Pa) at constant temperature and composition. | |
void | setState_TP (double T, double pres) override |
Set the temperature and pressure at the same time. | |
double | pressure () const override |
Returns the current pressure of the phase. | |
virtual void | updateStandardStateThermo () const |
Updates the standard state thermodynamic functions at the current T and P of the solution. | |
double | minTemp (size_t k=npos) const override |
Minimum temperature for which the thermodynamic data for the species or phase are valid. | |
double | maxTemp (size_t k=npos) const override |
Maximum temperature for which the thermodynamic data for the species are valid. | |
PDSS * | providePDSS (size_t k) |
const PDSS * | providePDSS (size_t k) const |
VPStandardStateTP () | |
Constructor. | |
bool | isCompressible () const override |
Return whether phase represents a compressible substance. | |
int | standardStateConvention () const override |
This method returns the convention used in specification of the standard state, of which there are currently two, temperature based, and variable pressure based. | |
void | getStandardChemPotentials (double *mu) const override |
Get the array of chemical potentials at unit activity for the species at their standard states at the current T and P of the solution. | |
void | getEnthalpy_RT (double *hrt) const override |
Get the nondimensional Enthalpy functions for the species at their standard states at the current T and P of the solution. | |
void | getEntropy_R (double *sr) const override |
Get the array of nondimensional Entropy functions for the standard state species at the current T and P of the solution. | |
void | getGibbs_RT (double *grt) const override |
Get the nondimensional Gibbs functions for the species in their standard states at the current T and P of the solution. | |
void | getPureGibbs (double *gpure) const override |
Get the Gibbs functions for the standard state of the species at the current T and P of the solution. | |
void | getIntEnergy_RT (double *urt) const override |
Returns the vector of nondimensional Internal Energies of the standard state species at the current T and P of the solution. | |
void | getCp_R (double *cpr) const override |
Get the nondimensional Heat Capacities at constant pressure for the species standard states at the current T and P of the solution. | |
void | getStandardVolumes (double *vol) const override |
Get the molar volumes of the species standard states at the current T and P of the solution. | |
virtual const vector< double > & | getStandardVolumes () const |
void | initThermo () override |
Initialize the ThermoPhase object after all species have been set up. | |
void | getSpeciesParameters (const string &name, AnyMap &speciesNode) const override |
Get phase-specific parameters of a Species object such that an identical one could be reconstructed and added to this phase. | |
bool | addSpecies (shared_ptr< Species > spec) override |
Add a Species to this Phase. | |
void | installPDSS (size_t k, unique_ptr< PDSS > &&pdss) |
Install a PDSS object for species k | |
virtual bool | addSpecies (shared_ptr< Species > spec) |
Add a Species to this Phase. | |
void | getEnthalpy_RT_ref (double *hrt) const override |
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of the solution and the reference pressure for the species. | |
void | getGibbs_RT_ref (double *grt) const override |
Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temperature of the solution and the reference pressure for the species. | |
void | getGibbs_ref (double *g) const override |
Returns the vector of the Gibbs function of the reference state at the current temperature of the solution and the reference pressure for the species. | |
void | getEntropy_R_ref (double *er) const override |
Returns the vector of nondimensional entropies of the reference state at the current temperature of the solution and the reference pressure for each species. | |
void | getCp_R_ref (double *cprt) const override |
Returns the vector of nondimensional constant pressure heat capacities of the reference state at the current temperature of the solution and reference pressure for each species. | |
void | getStandardVolumes_ref (double *vol) const override |
Get the molar volumes of the species reference states at the current T and P_ref of the solution. | |
Public Member Functions inherited from ThermoPhase | |
ThermoPhase ()=default | |
Constructor. | |
double | RT () const |
Return the Gas Constant multiplied by the current temperature. | |
double | equivalenceRatio () const |
Compute the equivalence ratio for the current mixture from available oxygen and required oxygen. | |
string | type () const override |
String indicating the thermodynamic model implemented. | |
virtual bool | isIdeal () const |
Boolean indicating whether phase is ideal. | |
virtual string | phaseOfMatter () const |
String indicating the mechanical phase of the matter in this Phase. | |
virtual double | refPressure () const |
Returns the reference pressure in Pa. | |
double | Hf298SS (const size_t k) const |
Report the 298 K Heat of Formation of the standard state of one species (J kmol-1) | |
virtual void | modifyOneHf298SS (const size_t k, const double Hf298New) |
Modify the value of the 298 K Heat of Formation of one species in the phase (J kmol-1) | |
virtual void | resetHf298 (const size_t k=npos) |
Restore the original heat of formation of one or more species. | |
bool | chargeNeutralityNecessary () const |
Returns the chargeNeutralityNecessity boolean. | |
virtual double | intEnergy_mole () const |
Molar internal energy. Units: J/kmol. | |
virtual double | gibbs_mole () const |
Molar Gibbs function. Units: J/kmol. | |
virtual double | isothermalCompressibility () const |
Returns the isothermal compressibility. Units: 1/Pa. | |
virtual double | thermalExpansionCoeff () const |
Return the volumetric thermal expansion coefficient. Units: 1/K. | |
virtual double | soundSpeed () const |
Return the speed of sound. Units: m/s. | |
void | setElectricPotential (double v) |
Set the electric potential of this phase (V). | |
double | electricPotential () const |
Returns the electric potential of this phase (V). | |
virtual int | activityConvention () const |
This method returns the convention used in specification of the activities, of which there are currently two, molar- and molality-based conventions. | |
void | getElectrochemPotentials (double *mu) const |
Get the species electrochemical potentials. | |
virtual void | getPartialMolarIntEnergies (double *ubar) const |
Return an array of partial molar internal energies for the species in the mixture. | |
virtual void | getIntEnergy_RT_ref (double *urt) const |
Returns the vector of nondimensional internal Energies of the reference state at the current temperature of the solution and the reference pressure for each species. | |
double | enthalpy_mass () const |
Specific enthalpy. Units: J/kg. | |
double | intEnergy_mass () const |
Specific internal energy. Units: J/kg. | |
double | entropy_mass () const |
Specific entropy. Units: J/kg/K. | |
double | gibbs_mass () const |
Specific Gibbs function. Units: J/kg. | |
double | cp_mass () const |
Specific heat at constant pressure. Units: J/kg/K. | |
double | cv_mass () const |
Specific heat at constant volume. Units: J/kg/K. | |
virtual void | setState_TPX (double t, double p, const double *x) |
Set the temperature (K), pressure (Pa), and mole fractions. | |
virtual void | setState_TPX (double t, double p, const Composition &x) |
Set the temperature (K), pressure (Pa), and mole fractions. | |
virtual void | setState_TPX (double t, double p, const string &x) |
Set the temperature (K), pressure (Pa), and mole fractions. | |
virtual void | setState_TPY (double t, double p, const double *y) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. | |
virtual void | setState_TPY (double t, double p, const Composition &y) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. | |
virtual void | setState_TPY (double t, double p, const string &y) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. | |
virtual void | setState_HP (double h, double p, double tol=1e-9) |
Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase. | |
virtual void | setState_UV (double u, double v, double tol=1e-9) |
Set the specific internal energy (J/kg) and specific volume (m^3/kg). | |
virtual void | setState_SP (double s, double p, double tol=1e-9) |
Set the specific entropy (J/kg/K) and pressure (Pa). | |
virtual void | setState_SV (double s, double v, double tol=1e-9) |
Set the specific entropy (J/kg/K) and specific volume (m^3/kg). | |
virtual void | setState_ST (double s, double t, double tol=1e-9) |
Set the specific entropy (J/kg/K) and temperature (K). | |
virtual void | setState_TV (double t, double v, double tol=1e-9) |
Set the temperature (K) and specific volume (m^3/kg). | |
virtual void | setState_PV (double p, double v, double tol=1e-9) |
Set the pressure (Pa) and specific volume (m^3/kg). | |
virtual void | setState_UP (double u, double p, double tol=1e-9) |
Set the specific internal energy (J/kg) and pressure (Pa). | |
virtual void | setState_VH (double v, double h, double tol=1e-9) |
Set the specific volume (m^3/kg) and the specific enthalpy (J/kg) | |
virtual void | setState_TH (double t, double h, double tol=1e-9) |
Set the temperature (K) and the specific enthalpy (J/kg) | |
virtual void | setState_SH (double s, double h, double tol=1e-9) |
Set the specific entropy (J/kg/K) and the specific enthalpy (J/kg) | |
virtual void | setState_DP (double rho, double p) |
Set the density (kg/m**3) and pressure (Pa) at constant composition. | |
virtual void | setState (const AnyMap &state) |
Set the state using an AnyMap containing any combination of properties supported by the thermodynamic model. | |
void | setMixtureFraction (double mixFrac, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) | |
void | setMixtureFraction (double mixFrac, const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) | |
void | setMixtureFraction (double mixFrac, const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) | |
double | mixtureFraction (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const |
Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. | |
double | mixtureFraction (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const |
Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. | |
double | mixtureFraction (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const |
Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. | |
void | setEquivalenceRatio (double phi, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. | |
void | setEquivalenceRatio (double phi, const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. | |
void | setEquivalenceRatio (double phi, const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. | |
double | equivalenceRatio (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. | |
double | equivalenceRatio (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. | |
double | equivalenceRatio (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. | |
double | stoichAirFuelRatio (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. | |
double | stoichAirFuelRatio (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. | |
double | stoichAirFuelRatio (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. | |
void | equilibrate (const string &XY, const string &solver="auto", double rtol=1e-9, int max_steps=50000, int max_iter=100, int estimate_equil=0, int log_level=0) |
Equilibrate a ThermoPhase object. | |
virtual void | setToEquilState (const double *mu_RT) |
This method is used by the ChemEquil equilibrium solver. | |
virtual bool | compatibleWithMultiPhase () const |
Indicates whether this phase type can be used with class MultiPhase for equilibrium calculations. | |
virtual double | critTemperature () const |
Critical temperature (K). | |
virtual double | critPressure () const |
Critical pressure (Pa). | |
virtual double | critVolume () const |
Critical volume (m3/kmol). | |
virtual double | critCompressibility () const |
Critical compressibility (unitless). | |
virtual double | critDensity () const |
Critical density (kg/m3). | |
virtual double | satTemperature (double p) const |
Return the saturation temperature given the pressure. | |
virtual double | satPressure (double t) |
Return the saturation pressure given the temperature. | |
virtual double | vaporFraction () const |
Return the fraction of vapor at the current conditions. | |
virtual void | setState_Tsat (double t, double x) |
Set the state to a saturated system at a particular temperature. | |
virtual void | setState_Psat (double p, double x) |
Set the state to a saturated system at a particular pressure. | |
void | setState_TPQ (double T, double P, double Q) |
Set the temperature, pressure, and vapor fraction (quality). | |
void | modifySpecies (size_t k, shared_ptr< Species > spec) override |
Modify the thermodynamic data associated with a species. | |
virtual MultiSpeciesThermo & | speciesThermo (int k=-1) |
Return a changeable reference to the calculation manager for species reference-state thermodynamic properties. | |
virtual const MultiSpeciesThermo & | speciesThermo (int k=-1) const |
void | initThermoFile (const string &inputFile, const string &id) |
Initialize a ThermoPhase object using an input file. | |
virtual void | setParameters (const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap()) |
Set equation of state parameters from an AnyMap phase description. | |
AnyMap | parameters (bool withInput=true) const |
Returns the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function. | |
const AnyMap & | input () const |
Access input data associated with the phase description. | |
AnyMap & | input () |
virtual void | getdlnActCoeffdlnN_numderiv (const size_t ld, double *const dlnActCoeffdlnN) |
virtual string | report (bool show_thermo=true, double threshold=-1e-14) const |
returns a summary of the state of the phase as a string | |
Public Member Functions inherited from Phase | |
Phase ()=default | |
Default constructor. | |
Phase (const Phase &)=delete | |
Phase & | operator= (const Phase &)=delete |
virtual bool | isPure () const |
Return whether phase represents a pure (single species) substance. | |
virtual bool | hasPhaseTransition () const |
Return whether phase represents a substance with phase transitions. | |
virtual bool | isCompressible () const |
Return whether phase represents a compressible substance. | |
virtual map< string, size_t > | nativeState () const |
Return a map of properties defining the native state of a substance. | |
string | nativeMode () const |
Return string acronym representing the native state of a Phase. | |
virtual vector< string > | fullStates () const |
Return a vector containing full states defining a phase. | |
virtual vector< string > | partialStates () const |
Return a vector of settable partial property sets within a phase. | |
virtual size_t | stateSize () const |
Return size of vector defining internal state of the phase. | |
void | saveState (vector< double > &state) const |
Save the current internal state of the phase. | |
virtual void | saveState (size_t lenstate, double *state) const |
Write to array 'state' the current internal state. | |
void | restoreState (const vector< double > &state) |
Restore a state saved on a previous call to saveState. | |
virtual void | restoreState (size_t lenstate, const double *state) |
Restore the state of the phase from a previously saved state vector. | |
double | molecularWeight (size_t k) const |
Molecular weight of species k . | |
void | getMolecularWeights (double *weights) const |
Copy the vector of molecular weights into array weights. | |
const vector< double > & | molecularWeights () const |
Return a const reference to the internal vector of molecular weights. | |
const vector< double > & | inverseMolecularWeights () const |
Return a const reference to the internal vector of molecular weights. | |
void | getCharges (double *charges) const |
Copy the vector of species charges into array charges. | |
virtual void | setMolesNoTruncate (const double *const N) |
Set the state of the object with moles in [kmol]. | |
double | elementalMassFraction (const size_t m) const |
Elemental mass fraction of element m. | |
double | elementalMoleFraction (const size_t m) const |
Elemental mole fraction of element m. | |
double | charge (size_t k) const |
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the magnitude of the electron charge. | |
double | chargeDensity () const |
Charge density [C/m^3]. | |
size_t | nDim () const |
Returns the number of spatial dimensions (1, 2, or 3) | |
void | setNDim (size_t ndim) |
Set the number of spatial dimensions (1, 2, or 3). | |
virtual bool | ready () const |
Returns a bool indicating whether the object is ready for use. | |
int | stateMFNumber () const |
Return the State Mole Fraction Number. | |
virtual void | invalidateCache () |
Invalidate any cached values which are normally updated only when a change in state is detected. | |
bool | caseSensitiveSpecies () const |
Returns true if case sensitive species names are enforced. | |
void | setCaseSensitiveSpecies (bool cflag=true) |
Set flag that determines whether case sensitive species are enforced in look-up operations, for example speciesIndex. | |
vector< double > | getCompositionFromMap (const Composition &comp) const |
Converts a Composition to a vector with entries for each species Species that are not specified are set to zero in the vector. | |
void | massFractionsToMoleFractions (const double *Y, double *X) const |
Converts a mixture composition from mole fractions to mass fractions. | |
void | moleFractionsToMassFractions (const double *X, double *Y) const |
Converts a mixture composition from mass fractions to mole fractions. | |
string | name () const |
Return the name of the phase. | |
void | setName (const string &nm) |
Sets the string name for the phase. | |
string | elementName (size_t m) const |
Name of the element with index m. | |
size_t | elementIndex (const string &name) const |
Return the index of element named 'name'. | |
const vector< string > & | elementNames () const |
Return a read-only reference to the vector of element names. | |
double | atomicWeight (size_t m) const |
Atomic weight of element m. | |
double | entropyElement298 (size_t m) const |
Entropy of the element in its standard state at 298 K and 1 bar. | |
int | atomicNumber (size_t m) const |
Atomic number of element m. | |
int | elementType (size_t m) const |
Return the element constraint type Possible types include: | |
int | changeElementType (int m, int elem_type) |
Change the element type of the mth constraint Reassigns an element type. | |
const vector< double > & | atomicWeights () const |
Return a read-only reference to the vector of atomic weights. | |
size_t | nElements () const |
Number of elements. | |
void | checkElementIndex (size_t m) const |
Check that the specified element index is in range. | |
void | checkElementArraySize (size_t mm) const |
Check that an array size is at least nElements(). | |
double | nAtoms (size_t k, size_t m) const |
Number of atoms of element m in species k . | |
size_t | speciesIndex (const string &name) const |
Returns the index of a species named 'name' within the Phase object. | |
string | speciesName (size_t k) const |
Name of the species with index k. | |
const vector< string > & | speciesNames () const |
Return a const reference to the vector of species names. | |
size_t | nSpecies () const |
Returns the number of species in the phase. | |
void | checkSpeciesIndex (size_t k) const |
Check that the specified species index is in range. | |
void | checkSpeciesArraySize (size_t kk) const |
Check that an array size is at least nSpecies(). | |
void | setMoleFractionsByName (const Composition &xMap) |
Set the species mole fractions by name. | |
void | setMoleFractionsByName (const string &x) |
Set the mole fractions of a group of species by name. | |
void | setMassFractionsByName (const Composition &yMap) |
Set the species mass fractions by name. | |
void | setMassFractionsByName (const string &x) |
Set the species mass fractions by name. | |
void | setState_TD (double t, double rho) |
Set the internally stored temperature (K) and density (kg/m^3) | |
Composition | getMoleFractionsByName (double threshold=0.0) const |
Get the mole fractions by name. | |
double | moleFraction (size_t k) const |
Return the mole fraction of a single species. | |
double | moleFraction (const string &name) const |
Return the mole fraction of a single species. | |
Composition | getMassFractionsByName (double threshold=0.0) const |
Get the mass fractions by name. | |
double | massFraction (size_t k) const |
Return the mass fraction of a single species. | |
double | massFraction (const string &name) const |
Return the mass fraction of a single species. | |
void | getMoleFractions (double *const x) const |
Get the species mole fraction vector. | |
virtual void | setMoleFractions (const double *const x) |
Set the mole fractions to the specified values. | |
virtual void | setMoleFractions_NoNorm (const double *const x) |
Set the mole fractions to the specified values without normalizing. | |
void | getMassFractions (double *const y) const |
Get the species mass fractions. | |
const double * | massFractions () const |
Return a const pointer to the mass fraction array. | |
virtual void | setMassFractions (const double *const y) |
Set the mass fractions to the specified values and normalize them. | |
virtual void | setMassFractions_NoNorm (const double *const y) |
Set the mass fractions to the specified values without normalizing. | |
virtual void | getConcentrations (double *const c) const |
Get the species concentrations (kmol/m^3). | |
virtual double | concentration (const size_t k) const |
Concentration of species k. | |
virtual void | setConcentrations (const double *const conc) |
Set the concentrations to the specified values within the phase. | |
virtual void | setConcentrationsNoNorm (const double *const conc) |
Set the concentrations without ignoring negative concentrations. | |
double | temperature () const |
Temperature (K). | |
virtual double | electronTemperature () const |
Electron Temperature (K) | |
virtual double | density () const |
Density (kg/m^3). | |
virtual double | molarDensity () const |
Molar density (kmol/m^3). | |
virtual double | molarVolume () const |
Molar volume (m^3/kmol). | |
virtual void | setDensity (const double density_) |
Set the internally stored density (kg/m^3) of the phase. | |
virtual void | setElectronTemperature (double etemp) |
Set the internally stored electron temperature of the phase (K). | |
double | mean_X (const double *const Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. | |
double | mean_X (const vector< double > &Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. | |
double | meanMolecularWeight () const |
The mean molecular weight. Units: (kg/kmol) | |
double | sum_xlogx () const |
Evaluate \( \sum_k X_k \ln X_k \). | |
size_t | addElement (const string &symbol, double weight=-12345.0, int atomicNumber=0, double entropy298=ENTROPY298_UNKNOWN, int elem_type=CT_ELEM_TYPE_ABSPOS) |
Add an element. | |
void | addSpeciesAlias (const string &name, const string &alias) |
Add a species alias (that is, a user-defined alternative species name). | |
virtual vector< string > | findIsomers (const Composition &compMap) const |
Return a vector with isomers names matching a given composition map. | |
virtual vector< string > | findIsomers (const string &comp) const |
Return a vector with isomers names matching a given composition string. | |
shared_ptr< Species > | species (const string &name) const |
Return the Species object for the named species. | |
shared_ptr< Species > | species (size_t k) const |
Return the Species object for species whose index is k. | |
void | ignoreUndefinedElements () |
Set behavior when adding a species containing undefined elements to just skip the species. | |
void | addUndefinedElements () |
Set behavior when adding a species containing undefined elements to add those elements to the phase. | |
void | throwUndefinedElements () |
Set the behavior when adding a species containing undefined elements to throw an exception. | |
Protected Attributes | |
size_t | numBinaryInteractions_ = 0 |
number of binary interaction expressions | |
vector< double > | m_HE_b_ij |
Enthalpy term for the binary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< double > | m_HE_c_ij |
Enthalpy term for the ternary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< double > | m_SE_b_ij |
Entropy term for the binary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< double > | m_SE_c_ij |
Entropy term for the ternary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< double > | m_VHE_b_ij |
Enthalpy term for the binary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< double > | m_VHE_c_ij |
Enthalpy term for the ternary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< double > | m_VSE_b_ij |
Entropy term for the binary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< double > | m_VSE_c_ij |
Entropy term for the ternary mole fraction interaction of the excess Gibbs free energy expression. | |
vector< size_t > | m_pSpecies_A_ij |
vector of species indices representing species A in the interaction | |
vector< size_t > | m_pSpecies_B_ij |
vector of species indices representing species B in the interaction | |
int | formMargules_ = 0 |
form of the Margules interaction expression | |
int | formTempModel_ = 0 |
form of the temperature dependence of the Margules interaction expression | |
Protected Attributes inherited from GibbsExcessVPSSTP | |
vector< double > | moleFractions_ |
Storage for the current values of the mole fractions of the species. | |
vector< double > | lnActCoeff_Scaled_ |
Storage for the current values of the activity coefficients of the species. | |
vector< double > | dlnActCoeffdT_Scaled_ |
Storage for the current derivative values of the gradients with respect to temperature of the log of the activity coefficients of the species. | |
vector< double > | d2lnActCoeffdT2_Scaled_ |
Storage for the current derivative values of the gradients with respect to temperature of the log of the activity coefficients of the species. | |
vector< double > | dlnActCoeffdlnN_diag_ |
Storage for the current derivative values of the gradients with respect to logarithm of the mole fraction of the log of the activity coefficients of the species. | |
vector< double > | dlnActCoeffdlnX_diag_ |
Storage for the current derivative values of the gradients with respect to logarithm of the mole fraction of the log of the activity coefficients of the species. | |
Array2D | dlnActCoeffdlnN_ |
Storage for the current derivative values of the gradients with respect to logarithm of the species mole number of the log of the activity coefficients of the species. | |
Protected Attributes inherited from VPStandardStateTP | |
double | m_Pcurrent = OneAtm |
Current value of the pressure - state variable. | |
double | m_minTemp = 0.0 |
The minimum temperature at which data for all species is valid. | |
double | m_maxTemp = BigNumber |
The maximum temperature at which data for all species is valid. | |
double | m_Tlast_ss = -1.0 |
The last temperature at which the standard state thermodynamic properties were calculated at. | |
double | m_Plast_ss = -1.0 |
The last pressure at which the Standard State thermodynamic properties were calculated at. | |
vector< unique_ptr< PDSS > > | m_PDSS_storage |
Storage for the PDSS objects for the species. | |
vector< double > | m_h0_RT |
Vector containing the species reference enthalpies at T = m_tlast and P = p_ref. | |
vector< double > | m_cp0_R |
Vector containing the species reference constant pressure heat capacities at T = m_tlast and P = p_ref. | |
vector< double > | m_g0_RT |
Vector containing the species reference Gibbs functions at T = m_tlast and P = p_ref. | |
vector< double > | m_s0_R |
Vector containing the species reference entropies at T = m_tlast and P = p_ref. | |
vector< double > | m_V0 |
Vector containing the species reference molar volumes. | |
vector< double > | m_hss_RT |
Vector containing the species Standard State enthalpies at T = m_tlast and P = m_plast. | |
vector< double > | m_cpss_R |
Vector containing the species Standard State constant pressure heat capacities at T = m_tlast and P = m_plast. | |
vector< double > | m_gss_RT |
Vector containing the species Standard State Gibbs functions at T = m_tlast and P = m_plast. | |
vector< double > | m_sss_R |
Vector containing the species Standard State entropies at T = m_tlast and P = m_plast. | |
vector< double > | m_Vss |
Vector containing the species standard state volumes at T = m_tlast and P = m_plast. | |
Protected Attributes inherited from ThermoPhase | |
MultiSpeciesThermo | m_spthermo |
Pointer to the calculation manager for species reference-state thermodynamic properties. | |
AnyMap | m_input |
Data supplied via setParameters. | |
double | m_phi = 0.0 |
Stored value of the electric potential for this phase. Units are Volts. | |
bool | m_chargeNeutralityNecessary = false |
Boolean indicating whether a charge neutrality condition is a necessity. | |
int | m_ssConvention = cSS_CONVENTION_TEMPERATURE |
Contains the standard state convention. | |
double | m_tlast = 0.0 |
last value of the temperature processed by reference state | |
Protected Attributes inherited from Phase | |
ValueCache | m_cache |
Cached for saved calculations within each ThermoPhase. | |
size_t | m_kk = 0 |
Number of species in the phase. | |
size_t | m_ndim = 3 |
Dimensionality of the phase. | |
vector< double > | m_speciesComp |
Atomic composition of the species. | |
vector< double > | m_speciesCharge |
Vector of species charges. length m_kk. | |
map< string, shared_ptr< Species > > | m_species |
UndefElement::behavior | m_undefinedElementBehavior = UndefElement::add |
Flag determining behavior when adding species with an undefined element. | |
bool | m_caseSensitiveSpecies = false |
Flag determining whether case sensitive species names are enforced. | |
Private Member Functions | |
void | initLengths () |
Initialize lengths of local variables after all species have been identified. | |
void | s_update_lnActCoeff () const |
Update the activity coefficients. | |
void | s_update_dlnActCoeff_dT () const |
Update the derivative of the log of the activity coefficients wrt T. | |
void | s_update_dlnActCoeff_dlnX_diag () const |
Update the derivative of the log of the activity coefficients wrt log(mole fraction) | |
void | s_update_dlnActCoeff_dlnN_diag () const |
Update the derivative of the log of the activity coefficients wrt log(moles) - diagonal only. | |
void | s_update_dlnActCoeff_dlnN () const |
Update the derivative of the log of the activity coefficients wrt log(moles_m) | |
Additional Inherited Members | |
Protected Member Functions inherited from GibbsExcessVPSSTP | |
void | compositionChanged () override |
Apply changes to the state which are needed after the composition changes. | |
void | calcDensity () override |
Calculate the density of the mixture using the partial molar volumes and mole fractions as input. | |
Protected Member Functions inherited from VPStandardStateTP | |
virtual void | calcDensity () |
Calculate the density of the mixture using the partial molar volumes and mole fractions as input. | |
virtual void | _updateStandardStateThermo () const |
Updates the standard state thermodynamic functions at the current T and P of the solution. | |
void | invalidateCache () override |
Invalidate any cached values which are normally updated only when a change in state is detected. | |
const vector< double > & | Gibbs_RT_ref () const |
virtual void | getParameters (AnyMap &phaseNode) const |
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function. | |
Protected Member Functions inherited from Phase | |
void | assertCompressible (const string &setter) const |
Ensure that phase is compressible. | |
void | assignDensity (const double density_) |
Set the internally stored constant density (kg/m^3) of the phase. | |
void | setMolecularWeight (const int k, const double mw) |
Set the molecular weight of a single species to a given value. | |
virtual void | compositionChanged () |
Apply changes to the state which are needed after the composition changes. | |
|
explicit |
Construct a MargulesVPSSTP object from an input file.
inputFile | Name of the input file containing the phase definition. If blank, an empty phase will be created. |
id | name (ID) of the phase in the input file. If empty, the first phase definition in the input file will be used. |
Definition at line 18 of file MargulesVPSSTP.cpp.
|
inlineoverridevirtual |
String indicating the thermodynamic model implemented.
Usually corresponds to the name of the derived class, less any suffixes such as "Phase", TP", "VPSS", etc.
Reimplemented from Phase.
Definition at line 226 of file MargulesVPSSTP.h.
|
overridevirtual |
Molar enthalpy. Units: J/kmol.
Reimplemented from ThermoPhase.
Definition at line 52 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Molar entropy. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 64 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Molar heat capacity at constant pressure. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 76 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Molar heat capacity at constant volume. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 88 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Get the array of non-dimensional molar-based ln activity coefficients at the current solution temperature, pressure, and solution concentration.
lnac | Output vector of ln activity coefficients. Length: m_kk. |
Reimplemented from ThermoPhase.
Definition at line 25 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Get the species chemical potentials. Units: J/kmol.
This function returns a vector of chemical potentials of the species in solution at the current temperature, pressure and mole fraction of the solution.
mu | Output vector of species chemical potentials. Length: m_kk. Units: J/kmol |
Reimplemented from ThermoPhase.
Definition at line 38 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Returns an array of partial molar enthalpies for the species in the mixture.
Units (J/kmol)
For this phase, the partial molar enthalpies are equal to the standard state enthalpies modified by the derivative of the molality-based activity coefficient wrt temperature
\[ \bar h_k(T,P) = h^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} \]
hbar | Vector of returned partial molar enthalpies (length m_kk, units = J/kmol) |
Reimplemented from ThermoPhase.
Definition at line 93 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Returns an array of partial molar entropies for the species in the mixture.
Units (J/kmol)
For this phase, the partial molar enthalpies are equal to the standard state enthalpies modified by the derivative of the activity coefficient wrt temperature
\[ \bar s_k(T,P) = s^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} - R \ln( \gamma_k X_k) - R T \frac{d \ln(\gamma_k) }{dT} \]
sbar | Vector of returned partial molar entropies (length m_kk, units = J/kmol/K) |
Reimplemented from ThermoPhase.
Definition at line 132 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Returns an array of partial molar entropies for the species in the mixture.
Units (J/kmol)
For this phase, the partial molar enthalpies are equal to the standard state enthalpies modified by the derivative of the activity coefficient wrt temperature
\[ ??????????????? \bar s_k(T,P) = s^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} - R \ln( \gamma_k X_k) - R T \frac{d \ln(\gamma_k) }{dT} ??????????????? \]
cpbar | Vector of returned partial molar heat capacities (length m_kk, units = J/kmol/K) |
Reimplemented from ThermoPhase.
Definition at line 112 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Return an array of partial molar volumes for the species in the mixture.
Units: m^3/kmol.
Frequently, for this class of thermodynamics representations, the excess Volume due to mixing is zero. Here, we set it as a default. It may be overridden in derived classes.
vbar | Output vector of species partial molar volumes. Length = m_kk. units are m^3/kmol. |
Reimplemented from GibbsExcessVPSSTP.
Definition at line 154 of file MargulesVPSSTP.cpp.
void getd2lnActCoeffdT2 | ( | double * | d2lnActCoeffdT2 | ) | const |
Get the array of temperature second derivatives of the log activity coefficients.
units = 1/Kelvin
d2lnActCoeffdT2 | Output vector of temperature 2nd derivatives of the log Activity Coefficients. length = m_kk |
Definition at line 324 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Get the array of temperature derivatives of the log activity coefficients.
units = 1/Kelvin
dlnActCoeffdT | Output vector of temperature derivatives of the log Activity Coefficients. length = m_kk |
Reimplemented from GibbsExcessVPSSTP.
Definition at line 316 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Initialize the ThermoPhase object after all species have been set up.
This method is provided to allow subclasses to perform any initialization required after all species have been added. For example, it might be used to resize internal work arrays that must have an entry for each species. The base class implementation does nothing, and subclasses that do not require initialization do not need to overload this method. Derived classes which do override this function should call their parent class's implementation of this function as their last action.
When importing from an AnyMap phase description (or from a YAML file), setupPhase() adds all the species, stores the input data in m_input, and then calls this method to set model parameters from the data stored in m_input.
Reimplemented from ThermoPhase.
Definition at line 179 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function.
This does not include user-defined fields available in input().
Reimplemented from ThermoPhase.
Definition at line 205 of file MargulesVPSSTP.cpp.
void addBinaryInteraction | ( | const string & | speciesA, |
const string & | speciesB, | ||
double | h0, | ||
double | h1, | ||
double | s0, | ||
double | s1, | ||
double | vh0, | ||
double | vh1, | ||
double | vs0, | ||
double | vs1 | ||
) |
Add a binary species interaction with the specified parameters.
speciesA | name of the first species |
speciesB | name of the second species |
h0 | first excess enthalpy coefficient [J/kmol] |
h1 | second excess enthalpy coefficient [J/kmol] |
s0 | first excess entropy coefficient [J/kmol/K] |
s1 | second excess entropy coefficient [J/kmol/K] |
vh0 | first enthalpy coefficient for excess volume [m^3/kmol] |
vh1 | second enthalpy coefficient for excess volume [m^3/kmol] |
vs0 | first entropy coefficient for excess volume [m^3/kmol/K] |
vs1 | second entropy coefficient for excess volume [m^3/kmol/K] |
Definition at line 239 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Get the change in activity coefficients wrt changes in state (temp, mole fraction, etc) along a line in parameter space or along a line in physical space.
dTds | Input of temperature change along the path |
dXds | Input vector of changes in mole fraction along the path. length = m_kk Along the path length it must be the case that the mole fractions sum to one. |
dlnActCoeffds | Output vector of the directional derivatives of the log Activity Coefficients along the path. length = m_kk units are 1/units(s). if s is a physical coordinate then the units are 1/m. |
Reimplemented from ThermoPhase.
Definition at line 332 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only.
For ideal mixtures (unity activity coefficients), this can return zero. Implementations should take the derivative of the logarithm of the activity coefficient with respect to the logarithm of the mole fraction variable that represents the standard state. This quantity is to be used in conjunction with derivatives of that mole fraction variable when the derivative of the chemical potential is taken.
units = dimensionless
dlnActCoeffdlnX_diag | Output vector of derivatives of the log Activity Coefficients wrt the mole fractions. length = m_kk |
Reimplemented from ThermoPhase.
Definition at line 460 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Get the array of log species mole number derivatives of the log activity coefficients.
For ideal mixtures (unity activity coefficients), this can return zero. Implementations should take the derivative of the logarithm of the activity coefficient with respect to the logarithm of the concentration- like variable (for example, moles) that represents the standard state. This quantity is to be used in conjunction with derivatives of that species mole number variable when the derivative of the chemical potential is taken.
units = dimensionless
dlnActCoeffdlnN_diag | Output vector of derivatives of the log Activity Coefficients. length = m_kk |
Reimplemented from ThermoPhase.
Definition at line 452 of file MargulesVPSSTP.cpp.
|
overridevirtual |
Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers.
Implementations should take the derivative of the logarithm of the activity coefficient with respect to a species log mole number (with all other species mole numbers held constant). The default treatment in the ThermoPhase object is to set this vector to zero.
units = 1 / kmol
dlnActCoeffdlnN[ ld * k + m] will contain the derivative of log act_coeff for the m-th species with respect to the number of moles of the k-th species.
\[ \frac{d \ln(\gamma_m) }{d \ln( n_k ) }\Bigg|_{n_i} \]
When implemented, this method is used within the VCS equilibrium solver to calculate the Jacobian elements, which accelerates convergence of the algorithm.
ld | Number of rows in the matrix |
dlnActCoeffdlnN | Output vector of derivatives of the log Activity Coefficients. length = m_kk * m_kk |
Reimplemented from GibbsExcessVPSSTP.
Definition at line 468 of file MargulesVPSSTP.cpp.
|
private |
Initialize lengths of local variables after all species have been identified.
Definition at line 234 of file MargulesVPSSTP.cpp.
|
private |
Update the activity coefficients.
This function will be called to update the internally stored natural logarithm of the activity coefficients
Definition at line 265 of file MargulesVPSSTP.cpp.
|
private |
Update the derivative of the log of the activity coefficients wrt T.
This function will be called to update the internally stored derivative of the natural logarithm of the activity coefficients wrt temperature.
Definition at line 287 of file MargulesVPSSTP.cpp.
|
private |
Update the derivative of the log of the activity coefficients wrt log(mole fraction)
This function will be called to update the internally stored derivative of the natural logarithm of the activity coefficients wrt logarithm of the mole fractions.
Definition at line 432 of file MargulesVPSSTP.cpp.
|
private |
Update the derivative of the log of the activity coefficients wrt log(moles) - diagonal only.
This function will be called to update the internally stored diagonal entries for the derivative of the natural logarithm of the activity coefficients wrt logarithm of the moles.
Definition at line 361 of file MargulesVPSSTP.cpp.
|
private |
Update the derivative of the log of the activity coefficients wrt log(moles_m)
This function will be called to update the internally stored derivative of the natural logarithm of the activity coefficients wrt logarithm of the mole number of species
Definition at line 393 of file MargulesVPSSTP.cpp.
|
protected |
number of binary interaction expressions
Definition at line 418 of file MargulesVPSSTP.h.
|
mutableprotected |
Enthalpy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 422 of file MargulesVPSSTP.h.
|
mutableprotected |
Enthalpy term for the ternary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 426 of file MargulesVPSSTP.h.
|
mutableprotected |
Entropy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 430 of file MargulesVPSSTP.h.
|
mutableprotected |
Entropy term for the ternary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 434 of file MargulesVPSSTP.h.
|
mutableprotected |
Enthalpy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 438 of file MargulesVPSSTP.h.
|
mutableprotected |
Enthalpy term for the ternary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 442 of file MargulesVPSSTP.h.
|
mutableprotected |
Entropy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 446 of file MargulesVPSSTP.h.
|
mutableprotected |
Entropy term for the ternary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 450 of file MargulesVPSSTP.h.
|
protected |
vector of species indices representing species A in the interaction
Each Margules excess Gibbs free energy term involves two species, A and B. This vector identifies species A.
Definition at line 457 of file MargulesVPSSTP.h.
|
protected |
vector of species indices representing species B in the interaction
Each Margules excess Gibbs free energy term involves two species, A and B. This vector identifies species B.
Definition at line 464 of file MargulesVPSSTP.h.
|
protected |
form of the Margules interaction expression
Currently there is only one form.
Definition at line 470 of file MargulesVPSSTP.h.
|
protected |
form of the temperature dependence of the Margules interaction expression
Currently there is only one form -> constant wrt temperature.
Definition at line 476 of file MargulesVPSSTP.h.