Cantera  3.0.0
Loading...
Searching...
No Matches
MixtureFugacityTP Class Reference

This is a filter class for ThermoPhase that implements some preparatory steps for efficiently handling mixture of gases that whose standard states are defined as ideal gases, but which describe also non-ideal solutions. More...

#include <MixtureFugacityTP.h>

Inheritance diagram for MixtureFugacityTP:
[legend]

Detailed Description

This is a filter class for ThermoPhase that implements some preparatory steps for efficiently handling mixture of gases that whose standard states are defined as ideal gases, but which describe also non-ideal solutions.

In addition a multicomponent liquid phase below the critical temperature of the mixture is also allowed. The main subclass is currently a mixture Redlich- Kwong class.

Several concepts are introduced. The first concept is there are temporary variables for holding the species standard state values of Cp, H, S, G, and V at the last temperature and pressure called. These functions are not recalculated if a new call is made using the previous temperature and pressure.

The other concept is that the current state of the mixture is tracked. The state variable is either GAS, LIQUID, or SUPERCRIT fluid. Additionally, the variable LiquidContent is used and may vary between 0 and 1.

Typically, only one liquid phase is allowed to be formed within these classes. Additionally, there is an inherent contradiction between three phase models and the ThermoPhase class. The ThermoPhase class is really only meant to represent a single instantiation of a phase. The three phase models may be in equilibrium with multiple phases of the fluid in equilibrium with each other. This has yet to be resolved.

This class is usually used for non-ideal gases.

Definition at line 65 of file MixtureFugacityTP.h.

Public Member Functions

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.
 
Constructors and Duplicators for MixtureFugacityTP
 MixtureFugacityTP ()=default
 Constructor.
 
Utilities
string type () const override
 String indicating the thermodynamic model implemented.
 
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 setForcedSolutionBranch (int solnBranch)
 Set the solution branch to force the ThermoPhase to exist on one branch or another.
 
int forcedSolutionBranch () const
 Report the solution branch which the solution is restricted to.
 
int reportSolnBranchActual () const
 Report the solution branch which the solution is actually on.
 
Molar Thermodynamic properties
double enthalpy_mole () const override
 Molar enthalpy. Units: J/kmol.
 
double entropy_mole () const override
 Molar entropy. Units: J/kmol/K.
 
Partial Molar Properties of the Solution
void getChemPotentials_RT (double *mu) const override
 Get the array of non-dimensional species chemical potentials These are partial molar Gibbs free energies.
 
Properties of the Standard State of the Species in the Solution

Within MixtureFugacityTP, these properties are calculated via a common routine, _updateStandardStateThermo(), which must be overloaded in inherited objects.

The values are cached within this object, and are not recalculated unless the temperature or pressure changes.

void getStandardChemPotentials (double *mu) const override
 Get the array of chemical potentials at unit activity.
 
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 Enthalpy 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 at their standard states of solution at the current T and P of the solution.
 
void getPureGibbs (double *gpure) const override
 Get the pure Gibbs free energies of each species.
 
void getIntEnergy_RT (double *urt) const override
 Returns the vector of nondimensional internal Energies of the standard state at the current temperature and pressure of the solution for each species.
 
void getCp_R (double *cpr) const override
 Get the nondimensional Heat Capacities at constant pressure for the standard state of the species at the current T and P.
 
void getStandardVolumes (double *vol) const override
 Get the molar volumes of each species in their standard states at the current T and P of the solution.
 
Initialization Methods - For Internal use

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().

bool addSpecies (shared_ptr< Species > spec) override
 Add a Species to this Phase.
 
- 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.
 
virtual double minTemp (size_t k=npos) const
 Minimum temperature for which the thermodynamic data for the species or phase are valid.
 
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.
 
virtual double maxTemp (size_t k=npos) const
 Maximum temperature for which the thermodynamic data for the species are valid.
 
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 cp_mole () const
 Molar heat capacity at constant pressure. Units: J/kmol/K.
 
virtual double cv_mole () const
 Molar heat capacity at constant volume. Units: J/kmol/K.
 
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.
 
virtual Units standardConcentrationUnits () const
 Returns the units of the "standard concentration" for this phase.
 
virtual double standardConcentration (size_t k=0) const
 Return the standard concentration for the kth species.
 
virtual double logStandardConc (size_t k=0) const
 Natural logarithm of the standard concentration of the kth species.
 
virtual void getActivities (double *a) const
 Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration.
 
virtual void getActivityCoefficients (double *ac) const
 Get the array of non-dimensional molar-based activity coefficients at the current solution temperature, pressure, and solution concentration.
 
virtual void getLnActivityCoefficients (double *lnac) const
 Get the array of non-dimensional molar-based ln activity coefficients at the current solution temperature, pressure, and solution concentration.
 
virtual void getChemPotentials (double *mu) const
 Get the species chemical potentials. Units: J/kmol.
 
void getElectrochemPotentials (double *mu) const
 Get the species electrochemical potentials.
 
virtual void getPartialMolarEnthalpies (double *hbar) const
 Returns an array of partial molar enthalpies for the species in the mixture.
 
virtual void getPartialMolarEntropies (double *sbar) const
 Returns an array of partial molar entropies of the species in the solution.
 
virtual void getPartialMolarIntEnergies (double *ubar) const
 Return an array of partial molar internal energies for the species in the mixture.
 
virtual void getPartialMolarCp (double *cpbar) const
 Return an array of partial molar heat capacities for the species in the mixture.
 
virtual void getPartialMolarVolumes (double *vbar) const
 Return an array of partial molar volumes 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_TP (double t, double p)
 Set the temperature (K) and pressure (Pa)
 
virtual void setState_PX (double p, double *x)
 Set the pressure (Pa) and mole fractions.
 
virtual void setState_PY (double p, double *y)
 Set the internally stored pressure (Pa) and mass fractions.
 
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)
 
void setState_RP (double rho, double p)
 Set the density (kg/m**3) and pressure (Pa) at constant composition.
 
virtual void setState_DP (double rho, double p)
 Set the density (kg/m**3) and pressure (Pa) at constant composition.
 
virtual void setState_RPX (double rho, double p, const double *x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions.
 
virtual void setState_RPX (double rho, double p, const Composition &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions.
 
virtual void setState_RPX (double rho, double p, const string &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions.
 
virtual void setState_RPY (double rho, double p, const double *y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions.
 
virtual void setState_RPY (double rho, double p, const Composition &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions.
 
virtual void setState_RPY (double rho, double p, const string &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions.
 
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 satTemperature (double p) const
 Return the saturation temperature given the pressure.
 
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).
 
bool addSpecies (shared_ptr< Species > spec) override
 Add a Species to this Phase.
 
void modifySpecies (size_t k, shared_ptr< Species > spec) override
 Modify the thermodynamic data associated with a species.
 
virtual MultiSpeciesThermospeciesThermo (int k=-1)
 Return a changeable reference to the calculation manager for species reference-state thermodynamic properties.
 
virtual const MultiSpeciesThermospeciesThermo (int k=-1) const
 
void initThermoFile (const string &inputFile, const string &id)
 Initialize a ThermoPhase object using an input file.
 
virtual void initThermo ()
 Initialize the ThermoPhase object after all species have been set up.
 
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.
 
virtual void getSpeciesParameters (const string &name, AnyMap &speciesNode) const
 Get phase-specific parameters of a Species object such that an identical one could be reconstructed and added to this phase.
 
const AnyMapinput () const
 Access input data associated with the phase description.
 
AnyMapinput ()
 
void invalidateCache () override
 Invalidate any cached values which are normally updated only when a change in state is detected.
 
virtual void getdlnActCoeffds (const double dTds, const double *const dXds, double *dlnActCoeffds) const
 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.
 
virtual void getdlnActCoeffdlnX_diag (double *dlnActCoeffdlnX_diag) const
 Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only.
 
virtual void getdlnActCoeffdlnN_diag (double *dlnActCoeffdlnN_diag) const
 Get the array of log species mole number derivatives of the log activity coefficients.
 
virtual void getdlnActCoeffdlnN (const size_t ld, double *const dlnActCoeffdlnN)
 Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers.
 
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
 
virtual void reportCSV (std::ofstream &csvFile) const
 returns a summary of the state of the phase to a comma separated file.
 
- Public Member Functions inherited from Phase
 Phase ()=default
 Default constructor.
 
 Phase (const Phase &)=delete
 
Phaseoperator= (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 (vector< double > &weights) const
 Copy the vector of molecular weights into vector weights.
 
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.
 
const double * moleFractdivMMW () const
 Returns a const pointer to the start of the moleFraction/MW array.
 
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.
 
void getAtoms (size_t k, double *atomArray) const
 Get a vector containing the atomic composition of 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.
 
string speciesSPName (int k) const
 Returns the expanded species name of a species, including the phase name This is guaranteed to be unique within a Cantera problem.
 
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_TRX (double t, double dens, const double *x)
 Set the internally stored temperature (K), density, and mole fractions.
 
void setState_TRX (double t, double dens, const Composition &x)
 Set the internally stored temperature (K), density, and mole fractions.
 
void setState_TRY (double t, double dens, const double *y)
 Set the internally stored temperature (K), density, and mass fractions.
 
void setState_TRY (double t, double dens, const Composition &y)
 Set the internally stored temperature (K), density, and mass fractions.
 
void setState_TNX (double t, double n, const double *x)
 Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions.
 
void setState_TR (double t, double rho)
 Set the internally stored temperature (K) and density (kg/m^3)
 
void setState_TD (double t, double rho)
 Set the internally stored temperature (K) and density (kg/m^3)
 
void setState_TX (double t, double *x)
 Set the internally stored temperature (K) and mole fractions.
 
void setState_TY (double t, double *y)
 Set the internally stored temperature (K) and mass fractions.
 
void setState_RX (double rho, double *x)
 Set the density (kg/m^3) and mole fractions.
 
void setState_RY (double rho, double *y)
 Set the density (kg/m^3) and mass fractions.
 
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 pressure () const
 Return the thermodynamic pressure (Pa).
 
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 setMolarDensity (const double molarDensity)
 Set the internally stored molar density (kmol/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< Speciesspecies (const string &name) const
 Return the Species object for the named species.
 
shared_ptr< Speciesspecies (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 Member Functions

void compositionChanged () override
 Apply changes to the state which are needed after the composition changes.
 
virtual void _updateReferenceStateThermo () const
 Updates the reference state thermodynamic functions at the current T of the solution.
 
Critical State Properties.
double critTemperature () const override
 Critical temperature (K).
 
double critPressure () const override
 Critical pressure (Pa).
 
double critVolume () const override
 Critical volume (m3/kmol).
 
double critCompressibility () const override
 Critical compressibility (unitless).
 
double critDensity () const override
 Critical density (kg/m3).
 
virtual void calcCriticalConditions (double &pc, double &tc, double &vc) const
 
int solveCubic (double T, double pres, double a, double b, double aAlpha, double Vroot[3], double an, double bn, double cn, double dn, double tc, double vc) const
 Solve the cubic equation of state.
 
- Protected Member Functions inherited from ThermoPhase
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.
 
virtual void getCsvReportData (vector< string > &names, vector< vector< double > > &data) const
 Fills names and data with the column names and species thermo properties to be included in the output of the reportCSV method.
 
- 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.
 

Protected Attributes

vector< double > m_tmpV
 Temporary storage - length = m_kk.
 
vector< double > moleFractions_
 Storage for the current values of the mole fractions of the species.
 
int iState_ = FLUID_GAS
 Current state of the fluid.
 
int forcedState_ = FLUID_UNDEFINED
 Force the system to be on a particular side of the spinodal curve.
 
vector< double > m_h0_RT
 Temporary storage for dimensionless reference state enthalpies.
 
vector< double > m_cp0_R
 Temporary storage for dimensionless reference state heat capacities.
 
vector< double > m_g0_RT
 Temporary storage for dimensionless reference state Gibbs energies.
 
vector< double > m_s0_R
 Temporary storage for dimensionless reference state entropies.
 
- 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.
 

Thermodynamic Values for the Species Reference States

There are also temporary variables for holding the species reference- state values of Cp, H, S, and V at the last temperature and reference pressure called.

These functions are not recalculated if a new call is made using the previous temperature. All calculations are done within the routine _updateRefStateThermo().

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.
 
const vector< double > & gibbs_RT_ref () const
 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.
 

Special Functions for fugacity classes

virtual double liquidVolEst (double TKelvin, double &pres) const
 Estimate for the molar volume of the liquid.
 
virtual double densityCalc (double TKelvin, double pressure, int phaseRequested, double rhoguess)
 Calculates the density given the temperature and the pressure and a guess at the density.
 
int phaseState (bool checkState=false) const
 Returns the Phase State flag for the current state of the object.
 
virtual double densSpinodalLiquid () const
 Return the value of the density at the liquid spinodal point (on the liquid side) for the current temperature.
 
virtual double densSpinodalGas () const
 Return the value of the density at the gas spinodal point (on the gas side) for the current temperature.
 
double calculatePsat (double TKelvin, double &molarVolGas, double &molarVolLiquid)
 Calculate the saturation pressure at the current mixture content for the given temperature.
 
double satPressure (double TKelvin) override
 Calculate the saturation pressure at the current mixture content for the given temperature.
 
void getActivityConcentrations (double *c) const override
 This method returns an array of generalized concentrations.
 
double z () const
 Calculate the value of z.
 
virtual double sresid () const
 Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
 
virtual double hresid () const
 Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.
 
virtual double psatEst (double TKelvin) const
 Estimate for the saturation pressure.
 
int corr0 (double TKelvin, double pres, double &densLiq, double &densGas, double &liqGRT, double &gasGRT)
 Utility routine in the calculation of the saturation pressure.
 
virtual double dpdVCalc (double TKelvin, double molarVol, double &presCalc) const
 Calculate the pressure and the pressure derivative given the temperature and the molar volume.
 
virtual void updateMixingExpressions ()
 

Constructor & Destructor Documentation

◆ MixtureFugacityTP()

MixtureFugacityTP ( )
default

Constructor.

Member Function Documentation

◆ type()

string type ( ) const
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.

Since
Starting in Cantera 3.0, the name returned by this method corresponds to the canonical name used in the YAML input format.

Reimplemented from Phase.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 78 of file MixtureFugacityTP.h.

◆ standardStateConvention()

int standardStateConvention ( ) const
overridevirtual

This method returns the convention used in specification of the standard state, of which there are currently two, temperature based, and variable pressure based.

Currently, there are two standard state conventions:

  • Temperature-based activities cSS_CONVENTION_TEMPERATURE 0
    • default
  • Variable Pressure and Temperature -based activities cSS_CONVENTION_VPSS 1
  • Thermodynamics is set via slave ThermoPhase objects with nothing being carried out at this ThermoPhase object level cSS_CONVENTION_SLAVE 2

Reimplemented from ThermoPhase.

Definition at line 21 of file MixtureFugacityTP.cpp.

◆ setForcedSolutionBranch()

void setForcedSolutionBranch ( int  solnBranch)

Set the solution branch to force the ThermoPhase to exist on one branch or another.

Parameters
solnBranchBranch that the solution is restricted to. the value -1 means gas. The value -2 means unrestricted. Values of zero or greater refer to species dominated condensed phases.

Definition at line 26 of file MixtureFugacityTP.cpp.

◆ forcedSolutionBranch()

int forcedSolutionBranch ( ) const

Report the solution branch which the solution is restricted to.

Returns
Branch that the solution is restricted to. the value -1 means gas. The value -2 means unrestricted. Values of zero or greater refer to species dominated condensed phases.

Definition at line 31 of file MixtureFugacityTP.cpp.

◆ reportSolnBranchActual()

int reportSolnBranchActual ( ) const

Report the solution branch which the solution is actually on.

Returns
Branch that the solution is restricted to. the value -1 means gas. The value -2 means superfluid.. Values of zero or greater refer to species dominated condensed phases.

Definition at line 36 of file MixtureFugacityTP.cpp.

◆ enthalpy_mole()

double enthalpy_mole ( ) const
overridevirtual

Molar enthalpy. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 42 of file MixtureFugacityTP.cpp.

◆ entropy_mole()

double entropy_mole ( ) const
overridevirtual

Molar entropy. Units: J/kmol/K.

Reimplemented from ThermoPhase.

Definition at line 50 of file MixtureFugacityTP.cpp.

◆ getChemPotentials_RT()

void getChemPotentials_RT ( double *  mu) const
overridevirtual

Get the array of non-dimensional species chemical potentials These are partial molar Gibbs free energies.

\( \mu_k / \hat R T \). Units: unitless

We close the loop on this function, here, calling getChemPotentials() and then dividing by RT. No need for child classes to handle.

Parameters
muOutput vector of non-dimensional species chemical potentials Length: m_kk.
Deprecated:
To be removed after Cantera 3.0. Use getChemPotentials() instead.

Reimplemented from ThermoPhase.

Reimplemented in RedlichKwongMFTP.

Definition at line 60 of file MixtureFugacityTP.cpp.

◆ getStandardChemPotentials()

void getStandardChemPotentials ( double *  mu) const
overridevirtual

Get the array of chemical potentials at unit activity.

These are the standard state chemical potentials \( \mu^0_k(T,P) \). The values are evaluated at the current temperature and pressure.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution.

Parameters
muOutput vector of standard state chemical potentials. length = m_kk. units are J / kmol.

Reimplemented from ThermoPhase.

Definition at line 72 of file MixtureFugacityTP.cpp.

◆ getEnthalpy_RT()

void getEnthalpy_RT ( double *  hrt) const
overridevirtual

Get the nondimensional Enthalpy functions for the species at their standard states at the current T and P of the solution.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution.

Parameters
hrtOutput vector of standard state enthalpies. length = m_kk. units are unitless.

Reimplemented from ThermoPhase.

Definition at line 81 of file MixtureFugacityTP.cpp.

◆ getEntropy_R()

void getEntropy_R ( double *  sr) const
overridevirtual

Get the array of nondimensional Enthalpy functions for the standard state species at the current T and P of the solution.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution.

Parameters
srOutput vector of nondimensional standard state entropies. length = m_kk.

Reimplemented from ThermoPhase.

Definition at line 86 of file MixtureFugacityTP.cpp.

◆ getGibbs_RT()

void getGibbs_RT ( double *  grt) const
overridevirtual

Get the nondimensional Gibbs functions for the species at their standard states of solution at the current T and P of the solution.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution.

Parameters
grtOutput vector of nondimensional standard state Gibbs free energies. length = m_kk.

Reimplemented from ThermoPhase.

Definition at line 95 of file MixtureFugacityTP.cpp.

◆ getPureGibbs()

void getPureGibbs ( double *  gpure) const
overridevirtual

Get the pure Gibbs free energies of each species.

Species are assumed to be in their standard states.

This is the same as getStandardChemPotentials().

Parameters
[out]gpureArray of standard state Gibbs free energies. length = m_kk. units are J/kmol.

Reimplemented from ThermoPhase.

Definition at line 104 of file MixtureFugacityTP.cpp.

◆ getIntEnergy_RT()

void getIntEnergy_RT ( double *  urt) const
overridevirtual

Returns the vector of nondimensional internal Energies of the standard state at the current temperature and pressure of the solution for each species.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution.

\[ u^{ss}_k(T,P) = h^{ss}_k(T) - P * V^{ss}_k \]

Parameters
urtOutput vector of nondimensional standard state internal energies. length = m_kk.

Reimplemented from ThermoPhase.

Definition at line 113 of file MixtureFugacityTP.cpp.

◆ getCp_R()

void getCp_R ( double *  cpr) const
overridevirtual

Get the nondimensional Heat Capacities at constant pressure for the standard state of the species at the current T and P.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution.

Parameters
cprOutput vector containing the the nondimensional Heat Capacities at constant pressure for the standard state of the species. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 121 of file MixtureFugacityTP.cpp.

◆ getStandardVolumes()

void getStandardVolumes ( double *  vol) const
overridevirtual

Get the molar volumes of each species in their standard states at the current T and P of the solution.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution.

units = m^3 / kmol

Parameters
volOutput vector of species volumes. length = m_kk. units = m^3 / kmol

Reimplemented from ThermoPhase.

Definition at line 126 of file MixtureFugacityTP.cpp.

◆ setTemperature()

void setTemperature ( const double  temp)
overridevirtual

Set the temperature of the phase.

Currently this passes down to setState_TP(). It does not make sense to calculate the standard state without first setting T and P.

Parameters
tempTemperature (kelvin)

Reimplemented from Phase.

Definition at line 191 of file MixtureFugacityTP.cpp.

◆ setPressure()

void setPressure ( double  p)
overridevirtual

Set the internally stored pressure (Pa) at constant temperature and composition.

Currently this passes down to setState_TP(). It does not make sense to calculate the standard state without first setting T and P.

Parameters
pinput Pressure (Pa)

Reimplemented from Phase.

Definition at line 200 of file MixtureFugacityTP.cpp.

◆ compositionChanged()

void compositionChanged ( )
overrideprotectedvirtual

Apply changes to the state which are needed after the composition changes.

This function is called after any call to setMassFractions(), setMoleFractions(), or similar. For phases which need to execute a callback after any change to the composition, it should be done by overriding this function rather than overriding all of the composition- setting functions. Derived class implementations of compositionChanged() should call the parent class method as well.

Reimplemented from Phase.

Definition at line 264 of file MixtureFugacityTP.cpp.

◆ _updateReferenceStateThermo()

void _updateReferenceStateThermo ( ) const
protectedvirtual

Updates the reference state thermodynamic functions at the current T of the solution.

This function must be called for every call to functions in this class. It checks to see whether the temperature has changed and thus the ss thermodynamics functions for all of the species must be recalculated.

This function is responsible for updating the following internal members:

  • m_h0_RT;
  • m_cp0_R;
  • m_g0_RT;
  • m_s0_R;

Definition at line 734 of file MixtureFugacityTP.cpp.

◆ getEnthalpy_RT_ref()

void getEnthalpy_RT_ref ( double *  hrt) const
overridevirtual

Returns the vector of nondimensional enthalpies of the reference state at the current temperature of the solution and the reference pressure for the species.

Parameters
hrtOutput vector containing the nondimensional reference state enthalpies. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 135 of file MixtureFugacityTP.cpp.

◆ getGibbs_RT_ref()

void getGibbs_RT_ref ( double *  grt) const
overridevirtual

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.

Parameters
grtOutput vector containing the nondimensional reference state Gibbs Free energies. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 140 of file MixtureFugacityTP.cpp.

◆ gibbs_RT_ref()

const vector< double > & gibbs_RT_ref ( ) const
protected

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.

Returns
Output vector contains the nondimensional Gibbs free energies of the reference state of the species length = m_kk, units = dimensionless.

Definition at line 151 of file MixtureFugacityTP.cpp.

◆ getGibbs_ref()

void getGibbs_ref ( double *  g) const
overridevirtual

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.

Parameters
gOutput vector containing the reference state Gibbs Free energies. Length: m_kk. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 145 of file MixtureFugacityTP.cpp.

◆ getEntropy_R_ref()

void getEntropy_R_ref ( double *  er) const
overridevirtual

Returns the vector of nondimensional entropies of the reference state at the current temperature of the solution and the reference pressure for each species.

Parameters
erOutput vector containing the nondimensional reference state entropies. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 156 of file MixtureFugacityTP.cpp.

◆ getCp_R_ref()

void getCp_R_ref ( double *  cprt) const
overridevirtual

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.

Parameters
cprtOutput vector of nondimensional reference state heat capacities at constant pressure for the species. Length: m_kk

Reimplemented from ThermoPhase.

Definition at line 161 of file MixtureFugacityTP.cpp.

◆ getStandardVolumes_ref()

void getStandardVolumes_ref ( double *  vol) const
overridevirtual

Get the molar volumes of the species reference states at the current T and P_ref of the solution.

units = m^3 / kmol

Parameters
volOutput vector containing the standard state volumes. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 166 of file MixtureFugacityTP.cpp.

◆ addSpecies()

bool addSpecies ( shared_ptr< Species spec)
overridevirtual

Add a Species to this Phase.

Returns true if the species was successfully added, or false if the species was ignored.

Derived classes which need to size arrays according to the number of species should overload this method. The derived class implementation should call the base class method, and, if this returns true (indicating that the species has been added), adjust their array sizes accordingly.

See also
ignoreUndefinedElements addUndefinedElements throwUndefinedElements

Reimplemented from Phase.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 173 of file MixtureFugacityTP.cpp.

◆ z()

double z ( ) const
protected

Calculate the value of z.

\[ z = \frac{P v}{R T} \]

returns the value of z

Definition at line 280 of file MixtureFugacityTP.cpp.

◆ sresid()

double sresid ( ) const
protectedvirtual

Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.

Here we use the current state conditions

Returns
the change in entropy in units of J kmol-1 K-1.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 285 of file MixtureFugacityTP.cpp.

◆ hresid()

double hresid ( ) const
protectedvirtual

Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.

Here we use the current state conditions

Returns
the change in entropy in units of J kmol-1.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 290 of file MixtureFugacityTP.cpp.

◆ psatEst()

double psatEst ( double  TKelvin) const
protectedvirtual

Estimate for the saturation pressure.

Note: this is only used as a starting guess for later routines that actually calculate an accurate value for the saturation pressure.

Parameters
TKelvintemperature in kelvin
Returns
the estimated saturation pressure at the given temperature

Definition at line 295 of file MixtureFugacityTP.cpp.

◆ liquidVolEst()

double liquidVolEst ( double  TKelvin,
double &  pres 
) const
virtual

Estimate for the molar volume of the liquid.

Note: this is only used as a starting guess for later routines that actually calculate an accurate value for the liquid molar volume. This routine doesn't change the state of the system.

Parameters
TKelvintemperature in kelvin
presPressure in Pa. This is used as an initial guess. If the routine needs to change the pressure to find a stable liquid state, the new pressure is returned in this variable.
Returns
the estimate of the liquid volume. If the liquid can't be found, this routine returns -1.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 306 of file MixtureFugacityTP.cpp.

◆ densityCalc()

double densityCalc ( double  TKelvin,
double  pressure,
int  phaseRequested,
double  rhoguess 
)
virtual

Calculates the density given the temperature and the pressure and a guess at the density.

Note, below T_c, this is a multivalued function. We do not cross the vapor dome in this. This is protected because it is called during setState_TP() routines. Infinite loops would result if it were not protected.

Parameters
TKelvinTemperature in Kelvin
pressurePressure in Pascals (Newton/m**2)
phaseRequestedint representing the phase whose density we are requesting. If we put a gas or liquid phase here, we will attempt to find a volume in that part of the volume space, only, in this routine. A value of FLUID_UNDEFINED means that we will accept anything.
rhoguessGuessed density of the fluid. A value of -1.0 indicates that there is no guessed density
Returns
We return the density of the fluid at the requested phase. If we have not found any acceptable density we return a -1. If we have found an acceptable density at a different phase, we return a -2.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 311 of file MixtureFugacityTP.cpp.

◆ corr0()

int corr0 ( double  TKelvin,
double  pres,
double &  densLiq,
double &  densGas,
double &  liqGRT,
double &  gasGRT 
)
protected

Utility routine in the calculation of the saturation pressure.

Parameters
TKelvintemperature (kelvin)
prespressure (Pascal)
[out]densLiqdensity of liquid
[out]densGasdensity of gas
[out]liqGRTdeltaG/RT of liquid
[out]gasGRTdeltaG/RT of gas

Definition at line 452 of file MixtureFugacityTP.cpp.

◆ phaseState()

int phaseState ( bool  checkState = false) const

Returns the Phase State flag for the current state of the object.

Parameters
checkStateIf true, this function does a complete check to see where in parameters space we are

There are three values:

  • WATER_GAS below the critical temperature but below the critical density
  • WATER_LIQUID below the critical temperature but above the critical density
  • WATER_SUPERCRIT above the critical temperature

Definition at line 481 of file MixtureFugacityTP.cpp.

◆ densSpinodalLiquid()

double densSpinodalLiquid ( ) const
virtual

Return the value of the density at the liquid spinodal point (on the liquid side) for the current temperature.

Returns
the density with units of kg m-3

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 522 of file MixtureFugacityTP.cpp.

◆ densSpinodalGas()

double densSpinodalGas ( ) const
virtual

Return the value of the density at the gas spinodal point (on the gas side) for the current temperature.

Returns
the density with units of kg m-3

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 527 of file MixtureFugacityTP.cpp.

◆ calculatePsat()

double calculatePsat ( double  TKelvin,
double &  molarVolGas,
double &  molarVolLiquid 
)

Calculate the saturation pressure at the current mixture content for the given temperature.

Parameters
TKelvin(input) Temperature (Kelvin)
molarVolGas(return) Molar volume of the gas
molarVolLiquid(return) Molar volume of the liquid
Returns
the saturation pressure at the given temperature

Definition at line 539 of file MixtureFugacityTP.cpp.

◆ satPressure()

double satPressure ( double  TKelvin)
overridevirtual

Calculate the saturation pressure at the current mixture content for the given temperature.

Parameters
TKelvinTemperature (Kelvin)
Returns
The saturation pressure at the given temperature

Reimplemented from ThermoPhase.

Definition at line 532 of file MixtureFugacityTP.cpp.

◆ getActivityConcentrations()

void getActivityConcentrations ( double *  c) const
overridevirtual

This method returns an array of generalized concentrations.

\( C^a_k \) are defined such that \( a_k = C^a_k / C^0_k, \) where \( C^0_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. Note that they may or may not have units of concentration — they might be partial pressures, mole fractions, or surface coverages, for example.

Parameters
cOutput array of generalized concentrations. The units depend upon the implementation of the reaction rate expressions within the phase.

Reimplemented from ThermoPhase.

Definition at line 271 of file MixtureFugacityTP.cpp.

◆ dpdVCalc()

double dpdVCalc ( double  TKelvin,
double  molarVol,
double &  presCalc 
) const
protectedvirtual

Calculate the pressure and the pressure derivative given the temperature and the molar volume.

Temperature and mole number are held constant

Parameters
TKelvintemperature in kelvin
molarVolmolar volume ( m3/kmol)
presCalcReturns the pressure.
Returns
the derivative of the pressure wrt the molar volume

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 729 of file MixtureFugacityTP.cpp.

◆ updateMixingExpressions()

void updateMixingExpressions ( )
protectedvirtual

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 448 of file MixtureFugacityTP.cpp.

◆ critTemperature()

double critTemperature ( ) const
overrideprotectedvirtual

Critical temperature (K).

Reimplemented from ThermoPhase.

Definition at line 756 of file MixtureFugacityTP.cpp.

◆ critPressure()

double critPressure ( ) const
overrideprotectedvirtual

Critical pressure (Pa).

Reimplemented from ThermoPhase.

Definition at line 763 of file MixtureFugacityTP.cpp.

◆ critVolume()

double critVolume ( ) const
overrideprotectedvirtual

Critical volume (m3/kmol).

Reimplemented from ThermoPhase.

Definition at line 770 of file MixtureFugacityTP.cpp.

◆ critCompressibility()

double critCompressibility ( ) const
overrideprotectedvirtual

Critical compressibility (unitless).

Reimplemented from ThermoPhase.

Definition at line 777 of file MixtureFugacityTP.cpp.

◆ critDensity()

double critDensity ( ) const
overrideprotectedvirtual

Critical density (kg/m3).

Reimplemented from ThermoPhase.

Definition at line 784 of file MixtureFugacityTP.cpp.

◆ calcCriticalConditions()

void calcCriticalConditions ( double &  pc,
double &  tc,
double &  vc 
) const
protectedvirtual

Definition at line 792 of file MixtureFugacityTP.cpp.

◆ solveCubic()

int solveCubic ( double  T,
double  pres,
double  a,
double  b,
double  aAlpha,
double  Vroot[3],
double  an,
double  bn,
double  cn,
double  dn,
double  tc,
double  vc 
) const
protected

Solve the cubic equation of state.

Returns the number of solutions found. For the gas phase solution, it returns a positive number (1 or 2). If it only finds the liquid branch solution, it will return -1 or -2 instead of 1 or 2. If it returns 0, then there is an error. The cubic equation is solved using Nickalls' method [28].

Parameters
Ttemperature (kelvin)
prespressure (Pa)
a"a" parameter in the non-ideal EoS [Pa-m^6/kmol^2]
b"b" parameter in the non-ideal EoS [m^3/kmol]
aAlphaa*alpha (temperature dependent function for P-R EoS, 1 for R-K EoS)
VrootRoots of the cubic equation for molar volume (m3/kmol)
anconstant used in cubic equation
bnconstant used in cubic equation
cnconstant used in cubic equation
dnconstant used in cubic equation
tcCritical temperature (kelvin)
vcCritical volume
Returns
the number of solutions found

Definition at line 797 of file MixtureFugacityTP.cpp.

Member Data Documentation

◆ m_tmpV

vector<double> m_tmpV
mutableprotected

Temporary storage - length = m_kk.

Definition at line 288 of file MixtureFugacityTP.h.

◆ moleFractions_

vector<double> moleFractions_
protected

Storage for the current values of the mole fractions of the species.

This vector is kept up-to-date when some the setState functions are called.

Definition at line 535 of file MixtureFugacityTP.h.

◆ iState_

int iState_ = FLUID_GAS
protected

Current state of the fluid.

There are three possible states of the fluid:

  • FLUID_GAS
  • FLUID_LIQUID
  • FLUID_SUPERCRIT

Definition at line 544 of file MixtureFugacityTP.h.

◆ forcedState_

int forcedState_ = FLUID_UNDEFINED
protected

Force the system to be on a particular side of the spinodal curve.

Definition at line 547 of file MixtureFugacityTP.h.

◆ m_h0_RT

vector<double> m_h0_RT
mutableprotected

Temporary storage for dimensionless reference state enthalpies.

Definition at line 550 of file MixtureFugacityTP.h.

◆ m_cp0_R

vector<double> m_cp0_R
mutableprotected

Temporary storage for dimensionless reference state heat capacities.

Definition at line 553 of file MixtureFugacityTP.h.

◆ m_g0_RT

vector<double> m_g0_RT
mutableprotected

Temporary storage for dimensionless reference state Gibbs energies.

Definition at line 556 of file MixtureFugacityTP.h.

◆ m_s0_R

vector<double> m_s0_R
mutableprotected

Temporary storage for dimensionless reference state entropies.

Definition at line 559 of file MixtureFugacityTP.h.


The documentation for this class was generated from the following files: