Cantera
3.1.0a1
|
Implementation of a multi-species Redlich-Kwong equation of state. More...
#include <RedlichKwongMFTP.h>
Implementation of a multi-species Redlich-Kwong equation of state.
Definition at line 19 of file RedlichKwongMFTP.h.
Public Member Functions | |
RedlichKwongMFTP (const string &infile="", const string &id="") | |
Construct a RedlichKwongMFTP object from an input file. More... | |
string | type () const override |
String indicating the thermodynamic model implemented. More... | |
double | standardConcentration (size_t k=0) const override |
Returns the standard concentration \( C^0_k \), which is used to normalize the generalized concentration. More... | |
void | getActivityCoefficients (double *ac) const override |
Get the array of non-dimensional activity coefficients at the current solution temperature, pressure, and solution concentration. More... | |
double | liquidVolEst (double TKelvin, double &pres) const override |
Estimate for the molar volume of the liquid. More... | |
double | densityCalc (double T, double pressure, int phase, double rhoguess) override |
Calculates the density given the temperature and the pressure and a guess at the density. More... | |
double | densSpinodalLiquid () const override |
Return the value of the density at the liquid spinodal point (on the liquid side) for the current temperature. More... | |
double | densSpinodalGas () const override |
Return the value of the density at the gas spinodal point (on the gas side) for the current temperature. More... | |
double | dpdVCalc (double TKelvin, double molarVol, double &presCalc) const override |
Calculate the pressure and the pressure derivative given the temperature and the molar volume. More... | |
double | isothermalCompressibility () const override |
Returns the isothermal compressibility. Units: 1/Pa. More... | |
double | thermalExpansionCoeff () const override |
Return the volumetric thermal expansion coefficient. Units: 1/K. More... | |
double | soundSpeed () const override |
Return the speed of sound. Units: m/s. More... | |
void | pressureDerivatives () const |
Calculate dpdV and dpdT at the current conditions. More... | |
void | updateMixingExpressions () override |
Update the a and b parameters. More... | |
void | calculateAB (double temp, double &aCalc, double &bCalc) const |
Calculate the a and the b parameters given the temperature. More... | |
double | da_dt () const |
void | calcCriticalConditions (double &pc, double &tc, double &vc) const override |
int | solveCubic (double T, double pres, double a, double b, double Vroot[3]) const |
Prepare variables and call the function to solve the cubic equation of state. More... | |
Molar Thermodynamic properties | |
double | cp_mole () const override |
Molar heat capacity at constant pressure. Units: J/kmol/K. More... | |
double | cv_mole () const override |
Molar heat capacity at constant volume. Units: J/kmol/K. More... | |
Mechanical Properties | |
double | pressure () const override |
Return the thermodynamic pressure (Pa). More... | |
Partial Molar Properties of the Solution | |
void | getChemPotentials (double *mu) const override |
Get the species chemical potentials. Units: J/kmol. More... | |
void | getPartialMolarEnthalpies (double *hbar) const override |
Returns an array of partial molar enthalpies for the species in the mixture. More... | |
void | getPartialMolarEntropies (double *sbar) const override |
Returns an array of partial molar entropies of the species in the solution. More... | |
void | getPartialMolarIntEnergies (double *ubar) const override |
Return an array of partial molar internal energies for the species in the mixture. More... | |
void | getPartialMolarCp (double *cpbar) const override |
Return an array of partial molar heat capacities for the species in the mixture. More... | |
void | getPartialMolarVolumes (double *vbar) const override |
Return an array of partial molar volumes for the species in the mixture. More... | |
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. More... | |
void | initThermo () override |
Initialize the ThermoPhase object after all species have been set up. More... | |
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. More... | |
void | setSpeciesCoeffs (const string &species, double a0, double a1, double b) |
Set the pure fluid interaction parameters for a species. More... | |
void | setBinaryCoeffs (const string &species_i, const string &species_j, double a0, double a1) |
Set values for the interaction parameter between two species. More... | |
Public Member Functions inherited from MixtureFugacityTP | |
void | setTemperature (const double temp) override |
Set the temperature of the phase. More... | |
void | setPressure (double p) override |
Set the internally stored pressure (Pa) at constant temperature and composition. More... | |
MixtureFugacityTP ()=default | |
Constructor. More... | |
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. More... | |
void | setForcedSolutionBranch (int solnBranch) |
Set the solution branch to force the ThermoPhase to exist on one branch or another. More... | |
int | forcedSolutionBranch () const |
Report the solution branch which the solution is restricted to. More... | |
int | reportSolnBranchActual () const |
Report the solution branch which the solution is actually on. More... | |
double | enthalpy_mole () const override |
Molar enthalpy. Units: J/kmol. More... | |
double | entropy_mole () const override |
Molar entropy. Units: J/kmol/K. More... | |
void | getStandardChemPotentials (double *mu) const override |
Get the array of chemical potentials at unit activity. More... | |
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. More... | |
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. More... | |
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. More... | |
void | getPureGibbs (double *gpure) const override |
Get the pure Gibbs free energies of each species. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
int | phaseState (bool checkState=false) const |
Returns the Phase State flag for the current state of the object. More... | |
double | calculatePsat (double TKelvin, double &molarVolGas, double &molarVolLiquid) |
Calculate the saturation pressure at the current mixture content for the given temperature. More... | |
double | satPressure (double TKelvin) override |
Calculate the saturation pressure at the current mixture content for the given temperature. More... | |
void | getActivityConcentrations (double *c) const override |
This method returns an array of generalized concentrations. More... | |
Public Member Functions inherited from ThermoPhase | |
ThermoPhase ()=default | |
Constructor. More... | |
double | RT () const |
Return the Gas Constant multiplied by the current temperature. More... | |
double | equivalenceRatio () const |
Compute the equivalence ratio for the current mixture from available oxygen and required oxygen. More... | |
string | type () const override |
String indicating the thermodynamic model implemented. More... | |
virtual bool | isIdeal () const |
Boolean indicating whether phase is ideal. More... | |
virtual string | phaseOfMatter () const |
String indicating the mechanical phase of the matter in this Phase. More... | |
virtual double | refPressure () const |
Returns the reference pressure in Pa. More... | |
virtual double | minTemp (size_t k=npos) const |
Minimum temperature for which the thermodynamic data for the species or phase are valid. More... | |
double | Hf298SS (const size_t k) const |
Report the 298 K Heat of Formation of the standard state of one species (J kmol-1) More... | |
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) More... | |
virtual void | resetHf298 (const size_t k=npos) |
Restore the original heat of formation of one or more species. More... | |
virtual double | maxTemp (size_t k=npos) const |
Maximum temperature for which the thermodynamic data for the species are valid. More... | |
bool | chargeNeutralityNecessary () const |
Returns the chargeNeutralityNecessity boolean. More... | |
virtual double | intEnergy_mole () const |
Molar internal energy. Units: J/kmol. More... | |
virtual double | gibbs_mole () const |
Molar Gibbs function. Units: J/kmol. More... | |
void | setElectricPotential (double v) |
Set the electric potential of this phase (V). More... | |
double | electricPotential () const |
Returns the electric potential of this phase (V). More... | |
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. More... | |
virtual Units | standardConcentrationUnits () const |
Returns the units of the "standard concentration" for this phase. More... | |
virtual double | logStandardConc (size_t k=0) const |
Natural logarithm of the standard concentration of the kth species. More... | |
virtual void | getActivities (double *a) const |
Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration. More... | |
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. More... | |
void | getElectrochemPotentials (double *mu) const |
Get the species electrochemical potentials. More... | |
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. More... | |
double | enthalpy_mass () const |
Specific enthalpy. Units: J/kg. More... | |
double | intEnergy_mass () const |
Specific internal energy. Units: J/kg. More... | |
double | entropy_mass () const |
Specific entropy. Units: J/kg/K. More... | |
double | gibbs_mass () const |
Specific Gibbs function. Units: J/kg. More... | |
double | cp_mass () const |
Specific heat at constant pressure. Units: J/kg/K. More... | |
double | cv_mass () const |
Specific heat at constant volume. Units: J/kg/K. More... | |
virtual void | setState_TPX (double t, double p, const double *x) |
Set the temperature (K), pressure (Pa), and mole fractions. More... | |
virtual void | setState_TPX (double t, double p, const Composition &x) |
Set the temperature (K), pressure (Pa), and mole fractions. More... | |
virtual void | setState_TPX (double t, double p, const string &x) |
Set the temperature (K), pressure (Pa), and mole fractions. More... | |
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. More... | |
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. More... | |
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. More... | |
virtual void | setState_TP (double t, double p) |
Set the temperature (K) and pressure (Pa) More... | |
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. More... | |
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). More... | |
virtual void | setState_SP (double s, double p, double tol=1e-9) |
Set the specific entropy (J/kg/K) and pressure (Pa). More... | |
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). More... | |
virtual void | setState_ST (double s, double t, double tol=1e-9) |
Set the specific entropy (J/kg/K) and temperature (K). More... | |
virtual void | setState_TV (double t, double v, double tol=1e-9) |
Set the temperature (K) and specific volume (m^3/kg). More... | |
virtual void | setState_PV (double p, double v, double tol=1e-9) |
Set the pressure (Pa) and specific volume (m^3/kg). More... | |
virtual void | setState_UP (double u, double p, double tol=1e-9) |
Set the specific internal energy (J/kg) and pressure (Pa). More... | |
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) More... | |
virtual void | setState_TH (double t, double h, double tol=1e-9) |
Set the temperature (K) and the specific enthalpy (J/kg) More... | |
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) More... | |
virtual void | setState_DP (double rho, double p) |
Set the density (kg/m**3) and pressure (Pa) at constant composition. More... | |
virtual void | setState (const AnyMap &state) |
Set the state using an AnyMap containing any combination of properties supported by the thermodynamic model. More... | |
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) More... | |
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) More... | |
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) More... | |
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. More... | |
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. More... | |
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. More... | |
void | setEquivalenceRatio (double phi, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. More... | |
void | setEquivalenceRatio (double phi, const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. More... | |
void | setEquivalenceRatio (double phi, const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
virtual void | setToEquilState (const double *mu_RT) |
This method is used by the ChemEquil equilibrium solver. More... | |
virtual bool | compatibleWithMultiPhase () const |
Indicates whether this phase type can be used with class MultiPhase for equilibrium calculations. More... | |
virtual double | satTemperature (double p) const |
Return the saturation temperature given the pressure. More... | |
virtual double | vaporFraction () const |
Return the fraction of vapor at the current conditions. More... | |
virtual void | setState_Tsat (double t, double x) |
Set the state to a saturated system at a particular temperature. More... | |
virtual void | setState_Psat (double p, double x) |
Set the state to a saturated system at a particular pressure. More... | |
void | setState_TPQ (double T, double P, double Q) |
Set the temperature, pressure, and vapor fraction (quality). More... | |
bool | addSpecies (shared_ptr< Species > spec) override |
Add a Species to this Phase. More... | |
void | modifySpecies (size_t k, shared_ptr< Species > spec) override |
Modify the thermodynamic data associated with a species. More... | |
virtual MultiSpeciesThermo & | speciesThermo (int k=-1) |
Return a changeable reference to the calculation manager for species reference-state thermodynamic properties. More... | |
virtual const MultiSpeciesThermo & | speciesThermo (int k=-1) const |
void | initThermoFile (const string &inputFile, const string &id) |
Initialize a ThermoPhase object using an input file. More... | |
virtual void | setParameters (const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap()) |
Set equation of state parameters from an AnyMap phase description. More... | |
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. More... | |
const AnyMap & | input () const |
Access input data associated with the phase description. More... | |
AnyMap & | input () |
void | invalidateCache () override |
Invalidate any cached values which are normally updated only when a change in state is detected. More... | |
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. More... | |
virtual void | getdlnActCoeffdlnX_diag (double *dlnActCoeffdlnX_diag) const |
Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only. More... | |
virtual void | getdlnActCoeffdlnN_diag (double *dlnActCoeffdlnN_diag) const |
Get the array of log species mole number derivatives of the log activity coefficients. More... | |
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. More... | |
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 More... | |
Public Member Functions inherited from Phase | |
Phase ()=default | |
Default constructor. More... | |
Phase (const Phase &)=delete | |
Phase & | operator= (const Phase &)=delete |
virtual bool | isPure () const |
Return whether phase represents a pure (single species) substance. More... | |
virtual bool | hasPhaseTransition () const |
Return whether phase represents a substance with phase transitions. More... | |
virtual bool | isCompressible () const |
Return whether phase represents a compressible substance. More... | |
virtual map< string, size_t > | nativeState () const |
Return a map of properties defining the native state of a substance. More... | |
string | nativeMode () const |
Return string acronym representing the native state of a Phase. More... | |
virtual vector< string > | fullStates () const |
Return a vector containing full states defining a phase. More... | |
virtual vector< string > | partialStates () const |
Return a vector of settable partial property sets within a phase. More... | |
virtual size_t | stateSize () const |
Return size of vector defining internal state of the phase. More... | |
void | saveState (vector< double > &state) const |
Save the current internal state of the phase. More... | |
virtual void | saveState (size_t lenstate, double *state) const |
Write to array 'state' the current internal state. More... | |
void | restoreState (const vector< double > &state) |
Restore a state saved on a previous call to saveState. More... | |
virtual void | restoreState (size_t lenstate, const double *state) |
Restore the state of the phase from a previously saved state vector. More... | |
double | molecularWeight (size_t k) const |
Molecular weight of species k . More... | |
void | getMolecularWeights (double *weights) const |
Copy the vector of molecular weights into array weights. More... | |
const vector< double > & | molecularWeights () const |
Return a const reference to the internal vector of molecular weights. More... | |
const vector< double > & | inverseMolecularWeights () const |
Return a const reference to the internal vector of molecular weights. More... | |
void | getCharges (double *charges) const |
Copy the vector of species charges into array charges. More... | |
virtual void | setMolesNoTruncate (const double *const N) |
Set the state of the object with moles in [kmol]. More... | |
double | elementalMassFraction (const size_t m) const |
Elemental mass fraction of element m. More... | |
double | elementalMoleFraction (const size_t m) const |
Elemental mole fraction of element m. More... | |
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. More... | |
double | chargeDensity () const |
Charge density [C/m^3]. More... | |
size_t | nDim () const |
Returns the number of spatial dimensions (1, 2, or 3) More... | |
void | setNDim (size_t ndim) |
Set the number of spatial dimensions (1, 2, or 3). More... | |
virtual bool | ready () const |
Returns a bool indicating whether the object is ready for use. More... | |
int | stateMFNumber () const |
Return the State Mole Fraction Number. More... | |
bool | caseSensitiveSpecies () const |
Returns true if case sensitive species names are enforced. More... | |
void | setCaseSensitiveSpecies (bool cflag=true) |
Set flag that determines whether case sensitive species are enforced in look-up operations, for example speciesIndex. More... | |
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. More... | |
void | massFractionsToMoleFractions (const double *Y, double *X) const |
Converts a mixture composition from mole fractions to mass fractions. More... | |
void | moleFractionsToMassFractions (const double *X, double *Y) const |
Converts a mixture composition from mass fractions to mole fractions. More... | |
string | name () const |
Return the name of the phase. More... | |
void | setName (const string &nm) |
Sets the string name for the phase. More... | |
string | elementName (size_t m) const |
Name of the element with index m. More... | |
size_t | elementIndex (const string &name) const |
Return the index of element named 'name'. More... | |
const vector< string > & | elementNames () const |
Return a read-only reference to the vector of element names. More... | |
double | atomicWeight (size_t m) const |
Atomic weight of element m. More... | |
double | entropyElement298 (size_t m) const |
Entropy of the element in its standard state at 298 K and 1 bar. More... | |
int | atomicNumber (size_t m) const |
Atomic number of element m. More... | |
int | elementType (size_t m) const |
Return the element constraint type Possible types include: More... | |
int | changeElementType (int m, int elem_type) |
Change the element type of the mth constraint Reassigns an element type. More... | |
const vector< double > & | atomicWeights () const |
Return a read-only reference to the vector of atomic weights. More... | |
size_t | nElements () const |
Number of elements. More... | |
void | checkElementIndex (size_t m) const |
Check that the specified element index is in range. More... | |
void | checkElementArraySize (size_t mm) const |
Check that an array size is at least nElements(). More... | |
double | nAtoms (size_t k, size_t m) const |
Number of atoms of element m in species k . More... | |
size_t | speciesIndex (const string &name) const |
Returns the index of a species named 'name' within the Phase object. More... | |
string | speciesName (size_t k) const |
Name of the species with index k. More... | |
const vector< string > & | speciesNames () const |
Return a const reference to the vector of species names. More... | |
size_t | nSpecies () const |
Returns the number of species in the phase. More... | |
void | checkSpeciesIndex (size_t k) const |
Check that the specified species index is in range. More... | |
void | checkSpeciesArraySize (size_t kk) const |
Check that an array size is at least nSpecies(). More... | |
void | setMoleFractionsByName (const Composition &xMap) |
Set the species mole fractions by name. More... | |
void | setMoleFractionsByName (const string &x) |
Set the mole fractions of a group of species by name. More... | |
void | setMassFractionsByName (const Composition &yMap) |
Set the species mass fractions by name. More... | |
void | setMassFractionsByName (const string &x) |
Set the species mass fractions by name. More... | |
void | setState_TD (double t, double rho) |
Set the internally stored temperature (K) and density (kg/m^3) More... | |
Composition | getMoleFractionsByName (double threshold=0.0) const |
Get the mole fractions by name. More... | |
double | moleFraction (size_t k) const |
Return the mole fraction of a single species. More... | |
double | moleFraction (const string &name) const |
Return the mole fraction of a single species. More... | |
Composition | getMassFractionsByName (double threshold=0.0) const |
Get the mass fractions by name. More... | |
double | massFraction (size_t k) const |
Return the mass fraction of a single species. More... | |
double | massFraction (const string &name) const |
Return the mass fraction of a single species. More... | |
void | getMoleFractions (double *const x) const |
Get the species mole fraction vector. More... | |
virtual void | setMoleFractions (const double *const x) |
Set the mole fractions to the specified values. More... | |
virtual void | setMoleFractions_NoNorm (const double *const x) |
Set the mole fractions to the specified values without normalizing. More... | |
void | getMassFractions (double *const y) const |
Get the species mass fractions. More... | |
const double * | massFractions () const |
Return a const pointer to the mass fraction array. More... | |
virtual void | setMassFractions (const double *const y) |
Set the mass fractions to the specified values and normalize them. More... | |
virtual void | setMassFractions_NoNorm (const double *const y) |
Set the mass fractions to the specified values without normalizing. More... | |
virtual void | getConcentrations (double *const c) const |
Get the species concentrations (kmol/m^3). More... | |
virtual double | concentration (const size_t k) const |
Concentration of species k. More... | |
virtual void | setConcentrations (const double *const conc) |
Set the concentrations to the specified values within the phase. More... | |
virtual void | setConcentrationsNoNorm (const double *const conc) |
Set the concentrations without ignoring negative concentrations. More... | |
double | temperature () const |
Temperature (K). More... | |
virtual double | electronTemperature () const |
Electron Temperature (K) More... | |
virtual double | density () const |
Density (kg/m^3). More... | |
virtual double | molarDensity () const |
Molar density (kmol/m^3). More... | |
virtual double | molarVolume () const |
Molar volume (m^3/kmol). More... | |
virtual void | setDensity (const double density_) |
Set the internally stored density (kg/m^3) of the phase. More... | |
virtual void | setElectronTemperature (double etemp) |
Set the internally stored electron temperature of the phase (K). More... | |
double | mean_X (const double *const Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. More... | |
double | mean_X (const vector< double > &Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. More... | |
double | meanMolecularWeight () const |
The mean molecular weight. Units: (kg/kmol) More... | |
double | sum_xlogx () const |
Evaluate \( \sum_k X_k \ln X_k \). More... | |
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. More... | |
void | addSpeciesAlias (const string &name, const string &alias) |
Add a species alias (that is, a user-defined alternative species name). More... | |
virtual vector< string > | findIsomers (const Composition &compMap) const |
Return a vector with isomers names matching a given composition map. More... | |
virtual vector< string > | findIsomers (const string &comp) const |
Return a vector with isomers names matching a given composition string. More... | |
shared_ptr< Species > | species (const string &name) const |
Return the Species object for the named species. More... | |
shared_ptr< Species > | species (size_t k) const |
Return the Species object for species whose index is k. More... | |
void | ignoreUndefinedElements () |
Set behavior when adding a species containing undefined elements to just skip the species. More... | |
void | addUndefinedElements () |
Set behavior when adding a species containing undefined elements to add those elements to the phase. More... | |
void | throwUndefinedElements () |
Set the behavior when adding a species containing undefined elements to throw an exception. More... | |
Protected Types | |
enum class | CoeffSource { EoS , CritProps , Database } |
Protected Member Functions | |
double | sresid () const override |
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture. More... | |
double | hresid () const override |
Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture. More... | |
Protected Member Functions inherited from MixtureFugacityTP | |
void | compositionChanged () override |
Apply changes to the state which are needed after the composition changes. More... | |
virtual void | _updateReferenceStateThermo () const |
Updates the reference state thermodynamic functions at the current T of the solution. More... | |
double | critTemperature () const override |
Critical temperature (K). More... | |
double | critPressure () const override |
Critical pressure (Pa). More... | |
double | critVolume () const override |
Critical volume (m3/kmol). More... | |
double | critCompressibility () const override |
Critical compressibility (unitless). More... | |
double | critDensity () const override |
Critical density (kg/m3). More... | |
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. More... | |
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. More... | |
double | z () const |
Calculate the value of z. More... | |
virtual double | psatEst (double TKelvin) const |
Estimate for the saturation pressure. More... | |
int | corr0 (double TKelvin, double pres, double &densLiq, double &densGas, double &liqGRT, double &gasGRT) |
Utility routine in the calculation of the saturation pressure. More... | |
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. More... | |
Protected Member Functions inherited from Phase | |
void | assertCompressible (const string &setter) const |
Ensure that phase is compressible. More... | |
void | assignDensity (const double density_) |
Set the internally stored constant density (kg/m^3) of the phase. More... | |
void | setMolecularWeight (const int k, const double mw) |
Set the molecular weight of a single species to a given value. More... | |
Protected Attributes | |
int | m_formTempParam = 0 |
Form of the temperature parameterization. More... | |
double | m_b_current = 0.0 |
Value of b in the equation of state. More... | |
double | m_a_current = 0.0 |
Value of a in the equation of state. More... | |
vector< double > | a_vec_Curr_ |
vector< double > | b_vec_Curr_ |
Array2D | a_coeff_vec |
map< string, map< string, pair< double, double > > > | m_binaryParameters |
Explicitly-specified binary interaction parameters. More... | |
vector< CoeffSource > | m_coeffSource |
For each species, specifies the source of the a and b coefficients. More... | |
int | NSolns_ = 0 |
double | Vroot_ [3] = {0.0, 0.0, 0.0} |
vector< double > | m_pp |
Temporary storage - length = m_kk. More... | |
vector< double > | m_partialMolarVolumes |
double | dpdV_ = 0.0 |
The derivative of the pressure wrt the volume. More... | |
double | dpdT_ = 0.0 |
The derivative of the pressure wrt the temperature. More... | |
vector< double > | dpdni_ |
Vector of derivatives of pressure wrt mole number. More... | |
Protected Attributes inherited from MixtureFugacityTP | |
vector< double > | m_tmpV |
Temporary storage - length = m_kk. More... | |
vector< double > | moleFractions_ |
Storage for the current values of the mole fractions of the species. More... | |
int | iState_ = FLUID_GAS |
Current state of the fluid. More... | |
int | forcedState_ = FLUID_UNDEFINED |
Force the system to be on a particular side of the spinodal curve. More... | |
vector< double > | m_h0_RT |
Temporary storage for dimensionless reference state enthalpies. More... | |
vector< double > | m_cp0_R |
Temporary storage for dimensionless reference state heat capacities. More... | |
vector< double > | m_g0_RT |
Temporary storage for dimensionless reference state Gibbs energies. More... | |
vector< double > | m_s0_R |
Temporary storage for dimensionless reference state entropies. More... | |
Protected Attributes inherited from ThermoPhase | |
MultiSpeciesThermo | m_spthermo |
Pointer to the calculation manager for species reference-state thermodynamic properties. More... | |
AnyMap | m_input |
Data supplied via setParameters. More... | |
double | m_phi = 0.0 |
Stored value of the electric potential for this phase. Units are Volts. More... | |
bool | m_chargeNeutralityNecessary = false |
Boolean indicating whether a charge neutrality condition is a necessity. More... | |
int | m_ssConvention = cSS_CONVENTION_TEMPERATURE |
Contains the standard state convention. More... | |
double | m_tlast = 0.0 |
last value of the temperature processed by reference state More... | |
Protected Attributes inherited from Phase | |
ValueCache | m_cache |
Cached for saved calculations within each ThermoPhase. More... | |
size_t | m_kk = 0 |
Number of species in the phase. More... | |
size_t | m_ndim = 3 |
Dimensionality of the phase. More... | |
vector< double > | m_speciesComp |
Atomic composition of the species. More... | |
vector< double > | m_speciesCharge |
Vector of species charges. length m_kk. More... | |
map< string, shared_ptr< Species > > | m_species |
UndefElement::behavior | m_undefinedElementBehavior = UndefElement::add |
Flag determining behavior when adding species with an undefined element. More... | |
bool | m_caseSensitiveSpecies = false |
Flag determining whether case sensitive species names are enforced. More... | |
Static Private Attributes | |
static const double | omega_a = 4.27480233540E-01 |
Omega constant for a -> value of a in terms of critical properties. More... | |
static const double | omega_b = 8.66403499650E-02 |
Omega constant for b. More... | |
static const double | omega_vc = 3.33333333333333E-01 |
Omega constant for the critical molar volume. More... | |
|
explicit |
Construct a RedlichKwongMFTP object from an input file.
infile | 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 24 of file RedlichKwongMFTP.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 MixtureFugacityTP.
Definition at line 31 of file RedlichKwongMFTP.h.
|
overridevirtual |
Molar heat capacity at constant pressure. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 100 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Molar heat capacity at constant volume. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 116 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Return the thermodynamic pressure (Pa).
Since the mass density, temperature, and mass fractions are stored, this method uses these values to implement the mechanical equation of state \( P(T, \rho, Y_1, \dots, Y_K) \).
\[ P = \frac{RT}{v-b_{mix}} - \frac{a_{mix}}{T^{0.5} v \left( v + b_{mix} \right) } \]
Reimplemented from Phase.
Definition at line 130 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Returns the standard concentration \( C^0_k \), which is used to normalize the generalized concentration.
This is defined as the concentration by which the generalized concentration is normalized to produce the activity. In many cases, this quantity will be the same for all species in a phase. Since the activity for an ideal gas mixture is simply the mole fraction, for an ideal gas \( C^0_k = P/\hat R T \).
k | Optional parameter indicating the species. The default is to assume this refers to species 0. |
Reimplemented from ThermoPhase.
Definition at line 141 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Get the array of non-dimensional activity coefficients at the current solution temperature, pressure, and solution concentration.
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. The activities are based on this standard state.
ac | Output vector of activity coefficients. Length: m_kk. |
Reimplemented from ThermoPhase.
Definition at line 147 of file RedlichKwongMFTP.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 179 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Returns an array of partial molar enthalpies for the species in the mixture.
Units (J/kmol)
hbar | Output vector of species partial molar enthalpies. Length: m_kk. units are J/kmol. |
Reimplemented from ThermoPhase.
Definition at line 213 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Returns an array of partial molar entropies of the species in the solution.
Units: J/kmol/K.
sbar | Output vector of species partial molar entropies. Length = m_kk. units are J/kmol/K. |
Reimplemented from ThermoPhase.
Definition at line 258 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Return an array of partial molar internal energies for the species in the mixture.
Units: J/kmol.
ubar | Output vector of species partial molar internal energies. Length = m_kk. units are J/kmol. |
Reimplemented from ThermoPhase.
Definition at line 309 of file RedlichKwongMFTP.cpp.
|
inlineoverridevirtual |
Return an array of partial molar heat capacities for the species in the mixture.
Units: J/kmol/K
cpbar | Output vector of species partial molar heat capacities at constant pressure. Length = m_kk. units are J/kmol/K. |
Reimplemented from ThermoPhase.
Definition at line 93 of file RedlichKwongMFTP.h.
|
overridevirtual |
Return an array of partial molar volumes for the species in the mixture.
Units: m^3/kmol.
vbar | Output vector of species partial molar volumes. Length = m_kk. units are m^3/kmol. |
Reimplemented from ThermoPhase.
Definition at line 320 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
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.
Reimplemented from MixtureFugacityTP.
Definition at line 353 of file RedlichKwongMFTP.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 372 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Get phase-specific parameters of a Species object such that an identical one could be reconstructed and added to this phase.
name | Name of the species |
speciesNode | Mapping to be populated with parameters |
Reimplemented from ThermoPhase.
Definition at line 471 of file RedlichKwongMFTP.cpp.
void setSpeciesCoeffs | ( | const string & | species, |
double | a0, | ||
double | a1, | ||
double | b | ||
) |
Set the pure fluid interaction parameters for a species.
The "a" parameter for species i in the Redlich-Kwong model is assumed to be a linear function of temperature:
\[ a = a_0 + a_1 T \]
species | Name of the species |
a0 | constant term in the expression for the "a" parameter of the specified species [Pa-m^6/kmol^2] |
a1 | temperature-proportional term in the expression for the "a" parameter of the specified species [Pa-m^6/kmol^2/K] |
b | "b" parameter in the Redlich-Kwong model [m^3/kmol] |
Definition at line 29 of file RedlichKwongMFTP.cpp.
void setBinaryCoeffs | ( | const string & | species_i, |
const string & | species_j, | ||
double | a0, | ||
double | a1 | ||
) |
Set values for the interaction parameter between two species.
The "a" parameter for interactions between species i and j is assumed by default to be computed as:
\[ a_{ij} = \sqrt(a_{i,0} a_{j,0}) + \sqrt(a_{i,1} a_{j,1}) T \]
This function overrides the defaults with the specified parameters:
\[ a_{ij} = a_{ij,0} + a_{ij,1} T \]
species_i | Name of one species |
species_j | Name of the other species |
a0 | constant term in the "a" expression [Pa-m^6/kmol^2] |
a1 | temperature-proportional term in the "a" expression [Pa-m^6/kmol^2/K] |
Definition at line 71 of file RedlichKwongMFTP.cpp.
|
overrideprotectedvirtual |
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
Here we use the current state conditions
Reimplemented from MixtureFugacityTP.
Definition at line 520 of file RedlichKwongMFTP.cpp.
|
overrideprotectedvirtual |
Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.
Here we use the current state conditions
Reimplemented from MixtureFugacityTP.
Definition at line 536 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
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.
TKelvin | temperature in kelvin |
pres | Pressure 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. |
Reimplemented from MixtureFugacityTP.
Definition at line 551 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
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.
TKelvin | Temperature in Kelvin |
pressure | Pressure in Pascals (Newton/m**2) |
phaseRequested | int 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. |
rhoguess | Guessed density of the fluid. A value of -1.0 indicates that there is no guessed density |
Reimplemented from MixtureFugacityTP.
Definition at line 583 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Return the value of the density at the liquid spinodal point (on the liquid side) for the current temperature.
Reimplemented from MixtureFugacityTP.
Definition at line 638 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Return the value of the density at the gas spinodal point (on the gas side) for the current temperature.
Reimplemented from MixtureFugacityTP.
Definition at line 660 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Calculate the pressure and the pressure derivative given the temperature and the molar volume.
Temperature and mole number are held constant
TKelvin | temperature in kelvin |
molarVol | molar volume ( m3/kmol) |
presCalc | Returns the pressure. |
Reimplemented from MixtureFugacityTP.
Definition at line 682 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Returns the isothermal compressibility. Units: 1/Pa.
The isothermal compressibility is defined as
\[ \kappa_T = -\frac{1}{v}\left(\frac{\partial v}{\partial P}\right)_T \]
or
\[ \kappa_T = \frac{1}{\rho}\left(\frac{\partial \rho}{\partial P}\right)_T \]
Reimplemented from ThermoPhase.
Definition at line 695 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Return the volumetric thermal expansion coefficient. Units: 1/K.
The thermal expansion coefficient is defined as
\[ \beta = \frac{1}{v}\left(\frac{\partial v}{\partial T}\right)_P \]
Reimplemented from ThermoPhase.
Definition at line 701 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Return the speed of sound. Units: m/s.
The speed of sound is defined as
\[ c = \sqrt{\left(\frac{\partial P}{\partial\rho}\right)_s} \]
Reimplemented from ThermoPhase.
Definition at line 707 of file RedlichKwongMFTP.cpp.
void pressureDerivatives | ( | ) | const |
Calculate dpdV and dpdT at the current conditions.
These are stored internally.
Definition at line 713 of file RedlichKwongMFTP.cpp.
|
overridevirtual |
Update the a and b parameters.
The a and the b parameters depend on the mole fraction and the temperature. This function updates the internal numbers based on the state of the object.
Reimplemented from MixtureFugacityTP.
Definition at line 728 of file RedlichKwongMFTP.cpp.
void calculateAB | ( | double | temp, |
double & | aCalc, | ||
double & | bCalc | ||
) | const |
Calculate the a and the b parameters given the temperature.
This function doesn't change the internal state of the object, so it is a const function. It does use the stored mole fractions in the object.
temp | Temperature (TKelvin) |
aCalc | (output) Returns the a value |
bCalc | (output) Returns the b value. |
Definition at line 765 of file RedlichKwongMFTP.cpp.
int solveCubic | ( | double | T, |
double | pres, | ||
double | a, | ||
double | b, | ||
double | Vroot[3] | ||
) | const |
Prepare variables and call the function to solve the cubic equation of state.
Definition at line 857 of file RedlichKwongMFTP.cpp.
|
protected |
Form of the temperature parameterization.
Definition at line 203 of file RedlichKwongMFTP.h.
|
protected |
Value of b in the equation of state.
m_b is a function of the temperature and the mole fraction.
Definition at line 209 of file RedlichKwongMFTP.h.
|
protected |
Value of a in the equation of state.
a_b is a function of the temperature and the mole fraction.
Definition at line 215 of file RedlichKwongMFTP.h.
|
protected |
Explicitly-specified binary interaction parameters.
Definition at line 223 of file RedlichKwongMFTP.h.
|
protected |
For each species, specifies the source of the a and b coefficients.
Definition at line 227 of file RedlichKwongMFTP.h.
|
mutableprotected |
Temporary storage - length = m_kk.
Definition at line 234 of file RedlichKwongMFTP.h.
|
mutableprotected |
The derivative of the pressure wrt the volume.
Calculated at the current conditions. temperature and mole number kept constant
Definition at line 244 of file RedlichKwongMFTP.h.
|
mutableprotected |
The derivative of the pressure wrt the temperature.
Calculated at the current conditions. Total volume and mole number kept constant
Definition at line 251 of file RedlichKwongMFTP.h.
|
mutableprotected |
Vector of derivatives of pressure wrt mole number.
Calculated at the current conditions. Total volume, temperature and other mole number kept constant
Definition at line 258 of file RedlichKwongMFTP.h.
|
staticprivate |
Omega constant for a -> value of a in terms of critical properties.
this was calculated from a small nonlinear solve
Definition at line 265 of file RedlichKwongMFTP.h.
|
staticprivate |
Omega constant for b.
Definition at line 268 of file RedlichKwongMFTP.h.
|
staticprivate |
Omega constant for the critical molar volume.
Definition at line 271 of file RedlichKwongMFTP.h.