Cantera  3.1.0a2
Loading...
Searching...
No Matches

Class for pressure dependent standard states corresponding to ionic solutes in electrolyte water. More...

#include <PDSS_HKFT.h>

Inheritance diagram for PDSS_HKFT:
[legend]

Detailed Description

Class for pressure dependent standard states corresponding to ionic solutes in electrolyte water.

Definition at line 27 of file PDSS_HKFT.h.

Public Member Functions

 PDSS_HKFT ()
 Default Constructor.
 
Molar Thermodynamic Properties of the Solution
double enthalpy_mole () const override
 Return the molar enthalpy in units of J kmol-1.
 
double intEnergy_mole () const override
 Return the molar internal Energy in units of J kmol-1.
 
double entropy_mole () const override
 Return the molar entropy in units of J kmol-1 K-1.
 
double gibbs_mole () const override
 Return the molar Gibbs free energy in units of J kmol-1.
 
double cp_mole () const override
 Return the molar const pressure heat capacity in units of J kmol-1 K-1.
 
double molarVolume () const override
 Return the molar volume at standard state.
 
double density () const override
 Return the standard state density at standard state.
 
Properties of the Reference State of the Species in the Solution
double refPressure () const
 
double gibbs_RT_ref () const override
 Return the molar Gibbs free energy divided by RT at reference pressure.
 
double enthalpy_RT_ref () const override
 Return the molar enthalpy divided by RT at reference pressure.
 
double entropy_R_ref () const override
 Return the molar entropy divided by R at reference pressure.
 
double cp_R_ref () const override
 Return the molar heat capacity divided by R at reference pressure.
 
double molarVolume_ref () const override
 Return the molar volume at reference pressure.
 
Mechanical Equation of State Properties
void setState_TP (double temp, double pres) override
 Set the internal temperature and pressure.
 
Initialization of the Object
void setParent (VPStandardStateTP *phase, size_t k) override
 Set the parent VPStandardStateTP object of this PDSS object.
 
void initThermo () override
 Initialization routine.
 
void setDeltaH0 (double dh0)
 Set enthalpy of formation at Pr, Tr [J/kmol].
 
void setDeltaG0 (double dg0)
 Set Gibbs free energy of formation at Pr, Tr [J/kmol].
 
void setS0 (double s0)
 Set entropy of formation at Pr, Tr [J/kmol/K].
 
void set_a (double *a)
 Set "a" coefficients (array of 4 elements).
 
void set_c (double *c)
 Set "c" coefficients (array of 2 elements).
 
void setOmega (double omega)
 Set omega [J/kmol].
 
void getParameters (AnyMap &eosNode) const override
 Store the parameters needed to reconstruct a copy of this PDSS object.
 
- Public Member Functions inherited from PDSS_Molar
double enthalpy_RT () const override
 Return the standard state molar enthalpy divided by RT.
 
double entropy_R () const override
 Return the standard state entropy divided by RT.
 
double gibbs_RT () const override
 Return the molar Gibbs free energy divided by RT.
 
double cp_R () const override
 Return the molar const pressure heat capacity divided by RT.
 
- Public Member Functions inherited from PDSS
virtual void setTemperature (double temp)
 Set the internal temperature.
 
virtual double temperature () const
 Return the current stored temperature.
 
virtual void setState_TP (double temp, double pres)
 Set the internal temperature and pressure.
 
virtual double critTemperature () const
 critical temperature
 
virtual double critPressure () const
 critical pressure
 
virtual double critDensity () const
 critical density
 
virtual double satPressure (double T)
 saturation pressure
 
double molecularWeight () const
 Return the molecular weight of the species in units of kg kmol-1.
 
void setMolecularWeight (double mw)
 Set the molecular weight of the species.
 
 PDSS ()=default
 Default Constructor.
 
 PDSS (const PDSS &b)=delete
 
PDSSoperator= (const PDSS &b)=delete
 
virtual double cv_mole () const
 Return the molar const volume heat capacity in units of J kmol-1 K-1.
 
double refPressure () const
 Return the reference pressure for this phase.
 
double minTemp () const
 return the minimum temperature
 
double maxTemp () const
 return the minimum temperature
 
virtual double pressure () const
 Returns the pressure (Pa)
 
virtual void setPressure (double pres)
 Sets the pressure in the object.
 
virtual double thermalExpansionCoeff () const
 Return the volumetric thermal expansion coefficient. Units: 1/K.
 
void setReferenceThermo (shared_ptr< SpeciesThermoInterpType > stit)
 Set the SpeciesThermoInterpType object used to calculate reference state properties.
 
void setParameters (const AnyMap &node)
 Set model parameters from an AnyMap phase description, for example from the equation-of-state field of a species definition.
 

Private Member Functions

double deltaG () const
 Main routine that actually calculates the Gibbs free energy difference between the reference state at Tr, Pr and T,P.
 
double deltaS () const
 Main routine that actually calculates the entropy difference between the reference state at Tr, Pr and T,P.
 
double ag (const double temp, const int ifunc=0) const
 Internal formula for the calculation of a_g()
 
double bg (const double temp, const int ifunc=0) const
 Internal formula for the calculation of b_g()
 
double g (const double temp, const double pres, const int ifunc=0) const
 function g appearing in the formulation
 
double f (const double temp, const double pres, const int ifunc=0) const
 Difference function f appearing in the formulation.
 
double gstar (const double temp, const double pres, const int ifunc=0) const
 Evaluate the Gstar value appearing in the HKFT formulation.
 
double LookupGe (const string &elemName)
 Function to look up Element Free Energies.
 
void convertDGFormation ()
 Translate a Gibbs free energy of formation value to a NIST-based Chemical potential.
 

Private Attributes

VPStandardStateTPm_tp
 Parent VPStandardStateTP (ThermoPhase) object.
 
size_t m_spindex
 Index of this species within the parent phase.
 
PDSS_Waterm_waterSS = nullptr
 Water standard state calculator.
 
UnitSystem m_units
 
double m_densWaterSS = -1.0
 density of standard-state water. internal temporary variable
 
unique_ptr< WaterPropsm_waterProps
 Pointer to the water property calculator.
 
double m_deltaG_formation_tr_pr = NAN
 Input value of deltaG of Formation at Tr and Pr (cal gmol-1)
 
double m_deltaH_formation_tr_pr = NAN
 Input value of deltaH of Formation at Tr and Pr (cal gmol-1)
 
double m_Mu0_tr_pr = 0.0
 Value of the Absolute Gibbs Free Energy NIST scale at T_r and P_r.
 
double m_Entrop_tr_pr = NAN
 Input value of S_j at Tr and Pr (cal gmol-1 K-1)
 
double m_a1 = 0.0
 Input a1 coefficient (cal gmol-1 bar-1)
 
double m_a2 = 0.0
 Input a2 coefficient (cal gmol-1)
 
double m_a3 = 0.0
 Input a3 coefficient (cal K gmol-1 bar-1)
 
double m_a4 = 0.0
 Input a4 coefficient (cal K gmol-1)
 
double m_c1 = 0.0
 Input c1 coefficient (cal gmol-1 K-1)
 
double m_c2 = 0.0
 Input c2 coefficient (cal K gmol-1)
 
double m_omega_pr_tr = 0.0
 Input omega_pr_tr coefficient(cal gmol-1)
 
double m_Y_pr_tr = 0.0
 y = dZdT = 1/(esp*esp) desp/dT at 298.15 and 1 bar
 
double m_Z_pr_tr = 0.0
 Z = -1 / relEpsilon at 298.15 and 1 bar.
 
double m_presR_bar = OneAtm * 1.0E-5
 Reference pressure is 1 atm in units of bar= 1.0132.
 
double m_domega_jdT_prtr = 0.0
 small value that is not quite zero
 
double m_charge_j = 0.0
 Charge of the ion.
 

Static Private Attributes

static int s_InputInconsistencyErrorExit = 1
 Static variable determining error exiting.
 

Additional Inherited Members

- Protected Attributes inherited from PDSS
double m_temp = -1.0
 Current temperature used by the PDSS object.
 
double m_pres = -1.0
 State of the system - pressure.
 
double m_p0 = -1.0
 Reference state pressure of the species.
 
double m_minTemp = -1.0
 Minimum temperature.
 
double m_maxTemp = 10000.0
 Maximum temperature.
 
double m_mw = 0.0
 Molecular Weight of the species.
 
AnyMap m_input
 Input data supplied via setParameters.
 
shared_ptr< SpeciesThermoInterpTypem_spthermo
 Pointer to the species thermodynamic property manager.
 

Constructor & Destructor Documentation

◆ PDSS_HKFT()

PDSS_HKFT ( )

Default Constructor.

Definition at line 23 of file PDSS_HKFT.cpp.

Member Function Documentation

◆ enthalpy_mole()

double enthalpy_mole ( ) const
overridevirtual

Return the molar enthalpy in units of J kmol-1.

Returns
the species standard state enthalpy in J kmol-1 at the current temperature and pressure.

Reimplemented from PDSS.

Definition at line 29 of file PDSS_HKFT.cpp.

◆ intEnergy_mole()

double intEnergy_mole ( ) const
overridevirtual

Return the molar internal Energy in units of J kmol-1.

Returns
The species standard state internal Energy in J kmol-1 at the current temperature and pressure.

Reimplemented from PDSS.

Definition at line 36 of file PDSS_HKFT.cpp.

◆ entropy_mole()

double entropy_mole ( ) const
overridevirtual

Return the molar entropy in units of J kmol-1 K-1.

Returns
The species standard state entropy in J kmol-1 K-1 at the current temperature and pressure.

Reimplemented from PDSS.

Definition at line 41 of file PDSS_HKFT.cpp.

◆ gibbs_mole()

double gibbs_mole ( ) const
overridevirtual

Return the molar Gibbs free energy in units of J kmol-1.

Returns
The species standard state Gibbs free energy in J kmol-1 at the current temperature and pressure.

Reimplemented from PDSS.

Definition at line 46 of file PDSS_HKFT.cpp.

◆ cp_mole()

double cp_mole ( ) const
overridevirtual

Return the molar const pressure heat capacity in units of J kmol-1 K-1.

Returns
The species standard state Cp in J kmol-1 K-1 at the current temperature and pressure.

Reimplemented from PDSS.

Definition at line 51 of file PDSS_HKFT.cpp.

◆ molarVolume()

double molarVolume ( ) const
overridevirtual

Return the molar volume at standard state.

Returns
The standard state molar volume at the current temperature and pressure. Units are m**3 kmol-1.

Reimplemented from PDSS.

Definition at line 107 of file PDSS_HKFT.cpp.

◆ density()

double density ( ) const
overridevirtual

Return the standard state density at standard state.

Returns
The standard state density at the current temperature and pressure. units are kg m-3

Reimplemented from PDSS.

Definition at line 150 of file PDSS_HKFT.cpp.

◆ refPressure()

double refPressure ( ) const
inline

Definition at line 51 of file PDSS_HKFT.h.

◆ gibbs_RT_ref()

double gibbs_RT_ref ( ) const
overridevirtual

Return the molar Gibbs free energy divided by RT at reference pressure.

Returns
The reference state Gibbs free energy at the current temperature, divided by RT.

Reimplemented from PDSS.

Definition at line 155 of file PDSS_HKFT.cpp.

◆ enthalpy_RT_ref()

double enthalpy_RT_ref ( ) const
overridevirtual

Return the molar enthalpy divided by RT at reference pressure.

Returns
The species reference state enthalpy at the current temperature, divided by RT.

Reimplemented from PDSS.

Definition at line 164 of file PDSS_HKFT.cpp.

◆ entropy_R_ref()

double entropy_R_ref ( ) const
overridevirtual

Return the molar entropy divided by R at reference pressure.

Returns
The species reference state entropy at the current temperature, divided by R.

Reimplemented from PDSS.

Definition at line 173 of file PDSS_HKFT.cpp.

◆ cp_R_ref()

double cp_R_ref ( ) const
overridevirtual

Return the molar heat capacity divided by R at reference pressure.

Returns
The species reference state heat capacity divided by R at the current temperature.

Reimplemented from PDSS.

Definition at line 182 of file PDSS_HKFT.cpp.

◆ molarVolume_ref()

double molarVolume_ref ( ) const
overridevirtual

Return the molar volume at reference pressure.

Returns
The reference state molar volume. units are m**3 kmol-1.

Reimplemented from PDSS.

Definition at line 191 of file PDSS_HKFT.cpp.

◆ setState_TP()

void setState_TP ( double  temp,
double  pres 
)
overridevirtual

Set the internal temperature and pressure.

Parameters
tempTemperature (Kelvin)
prespressure (Pascals)

Reimplemented from PDSS.

Definition at line 200 of file PDSS_HKFT.cpp.

◆ setParent()

void setParent ( VPStandardStateTP phase,
size_t  k 
)
inlineoverridevirtual

Set the parent VPStandardStateTP object of this PDSS object.

This information is only used by certain PDSS subclasses

Parameters
phasePointer to the parent phase
kIndex of this species in the phase

Reimplemented from PDSS.

Definition at line 71 of file PDSS_HKFT.h.

◆ initThermo()

void initThermo ( )
overridevirtual

Initialization routine.

This is a cascading call, where each level should call the the parent level.

Reimplemented from PDSS.

Definition at line 206 of file PDSS_HKFT.cpp.

◆ setDeltaH0()

void setDeltaH0 ( double  dh0)

Set enthalpy of formation at Pr, Tr [J/kmol].

Definition at line 315 of file PDSS_HKFT.cpp.

◆ setDeltaG0()

void setDeltaG0 ( double  dg0)

Set Gibbs free energy of formation at Pr, Tr [J/kmol].

Definition at line 319 of file PDSS_HKFT.cpp.

◆ setS0()

void setS0 ( double  s0)

Set entropy of formation at Pr, Tr [J/kmol/K].

Definition at line 323 of file PDSS_HKFT.cpp.

◆ set_a()

void set_a ( double *  a)

Set "a" coefficients (array of 4 elements).

Units of each coefficient are [J/kmol/Pa, J/kmol, J*K/kmol/Pa, J*K/kmol]

Definition at line 327 of file PDSS_HKFT.cpp.

◆ set_c()

void set_c ( double *  c)

Set "c" coefficients (array of 2 elements).

Units of each coefficient are [J/kmol/K, J*K/kmol]

Definition at line 334 of file PDSS_HKFT.cpp.

◆ setOmega()

void setOmega ( double  omega)

Set omega [J/kmol].

Definition at line 339 of file PDSS_HKFT.cpp.

◆ getParameters()

void getParameters ( AnyMap eosNode) const
overridevirtual

Store the parameters needed to reconstruct a copy of this PDSS object.

Reimplemented from PDSS.

Definition at line 343 of file PDSS_HKFT.cpp.

◆ deltaG()

double deltaG ( ) const
private

Main routine that actually calculates the Gibbs free energy difference between the reference state at Tr, Pr and T,P.

This is Eqn. 59 in Johnson et al. [14].

Definition at line 365 of file PDSS_HKFT.cpp.

◆ deltaS()

double deltaS ( ) const
private

Main routine that actually calculates the entropy difference between the reference state at Tr, Pr and T,P.

This is Eqn. 61 in Johnson et al. [14]. Actually, there appears to be an error in the latter. This is a correction.

Definition at line 398 of file PDSS_HKFT.cpp.

◆ ag()

double ag ( const double  temp,
const int  ifunc = 0 
) const
private

Internal formula for the calculation of a_g()

The output of this is in units of Angstroms

Parameters
tempTemperature (K)
ifuncparameters specifying the desired information
  • 0 function value
  • 1 derivative wrt temperature
  • 2 2nd derivative wrt temperature
  • 3 derivative wrt pressure

Definition at line 438 of file PDSS_HKFT.cpp.

◆ bg()

double bg ( const double  temp,
const int  ifunc = 0 
) const
private

Internal formula for the calculation of b_g()

the output of this is unitless

Parameters
tempTemperature (K)
ifuncparameters specifying the desired information
  • 0 function value
  • 1 derivative wrt temperature
  • 2 2nd derivative wrt temperature
  • 3 derivative wrt pressure

Definition at line 452 of file PDSS_HKFT.cpp.

◆ g()

double g ( const double  temp,
const double  pres,
const int  ifunc = 0 
) const
private

function g appearing in the formulation

Function \( g \) (Eqn. 49) appearing in the Johnson et al. [14] formulation.

Parameters
tempTemperature kelvin
presPressure (pascal)
ifuncparameters specifying the desired information
  • 0 function value
  • 1 derivative wrt temperature
  • 2 2nd derivative wrt temperature
  • 3 derivative wrt pressure

Definition at line 500 of file PDSS_HKFT.cpp.

◆ f()

double f ( const double  temp,
const double  pres,
const int  ifunc = 0 
) const
private

Difference function f appearing in the formulation.

Function \( f \) (Eqn. 52) appearing in the Johnson et al. [14] formulation of \( \omega_j \) (Eqn. 46).

Parameters
tempTemperature kelvin
presPressure (pascal)
ifuncparameters specifying the desired information
  • 0 function value
  • 1 derivative wrt temperature
  • 2 2nd derivative wrt temperature
  • 3 derivative wrt pressure

Definition at line 466 of file PDSS_HKFT.cpp.

◆ gstar()

double gstar ( const double  temp,
const double  pres,
const int  ifunc = 0 
) const
private

Evaluate the Gstar value appearing in the HKFT formulation.

Parameters
tempTemperature kelvin
presPressure (pascal)
ifuncparameters specifying the desired information
  • 0 function value
  • 1 derivative wrt temperature
  • 2 2nd derivative wrt temperature
  • 3 derivative wrt pressure

Definition at line 553 of file PDSS_HKFT.cpp.

◆ LookupGe()

double LookupGe ( const string &  elemName)
private

Function to look up Element Free Energies.

This function looks up the argument string in the element database and returns the associated 298 K Gibbs Free energy of the element in its stable state.

Parameters
elemNameString. Only the first 3 characters are significant
Returns
value contains the Gibbs free energy for that element
Exceptions
CanteraErrorIf a match is not found, a CanteraError is thrown as well

Definition at line 561 of file PDSS_HKFT.cpp.

◆ convertDGFormation()

void convertDGFormation ( )
private

Translate a Gibbs free energy of formation value to a NIST-based Chemical potential.

Internally, this function is used to translate the input value, m_deltaG_formation_tr_pr, to the internally stored value, m_Mu0_tr_pr.

Definition at line 575 of file PDSS_HKFT.cpp.

Member Data Documentation

◆ m_tp

VPStandardStateTP* m_tp
private

Parent VPStandardStateTP (ThermoPhase) object.

Definition at line 100 of file PDSS_HKFT.h.

◆ m_spindex

size_t m_spindex
private

Index of this species within the parent phase.

Definition at line 101 of file PDSS_HKFT.h.

◆ m_waterSS

PDSS_Water* m_waterSS = nullptr
private

Water standard state calculator.

derived from the equation of state for water. This object doesn't own the object. Just a shallow pointer.

Definition at line 213 of file PDSS_HKFT.h.

◆ m_units

UnitSystem m_units
private

Definition at line 215 of file PDSS_HKFT.h.

◆ m_densWaterSS

double m_densWaterSS = -1.0
mutableprivate

density of standard-state water. internal temporary variable

Definition at line 218 of file PDSS_HKFT.h.

◆ m_waterProps

unique_ptr<WaterProps> m_waterProps
private

Pointer to the water property calculator.

Definition at line 221 of file PDSS_HKFT.h.

◆ m_deltaG_formation_tr_pr

double m_deltaG_formation_tr_pr = NAN
private

Input value of deltaG of Formation at Tr and Pr (cal gmol-1)

Tr = 298.15 Pr = 1 atm

This is the delta G for the formation reaction of the ion from elements in their stable state at Tr, Pr.

Definition at line 230 of file PDSS_HKFT.h.

◆ m_deltaH_formation_tr_pr

double m_deltaH_formation_tr_pr = NAN
private

Input value of deltaH of Formation at Tr and Pr (cal gmol-1)

Tr = 298.15 Pr = 1 atm

This is the delta H for the formation reaction of the ion from elements in their stable state at Tr, Pr.

Definition at line 239 of file PDSS_HKFT.h.

◆ m_Mu0_tr_pr

double m_Mu0_tr_pr = 0.0
private

Value of the Absolute Gibbs Free Energy NIST scale at T_r and P_r.

This is the NIST scale value of Gibbs free energy at T_r = 298.15 and P_r = 1 atm.

J kmol-1

Definition at line 248 of file PDSS_HKFT.h.

◆ m_Entrop_tr_pr

double m_Entrop_tr_pr = NAN
private

Input value of S_j at Tr and Pr (cal gmol-1 K-1)

Tr = 298.15 Pr = 1 atm

Definition at line 254 of file PDSS_HKFT.h.

◆ m_a1

double m_a1 = 0.0
private

Input a1 coefficient (cal gmol-1 bar-1)

Definition at line 257 of file PDSS_HKFT.h.

◆ m_a2

double m_a2 = 0.0
private

Input a2 coefficient (cal gmol-1)

Definition at line 260 of file PDSS_HKFT.h.

◆ m_a3

double m_a3 = 0.0
private

Input a3 coefficient (cal K gmol-1 bar-1)

Definition at line 263 of file PDSS_HKFT.h.

◆ m_a4

double m_a4 = 0.0
private

Input a4 coefficient (cal K gmol-1)

Definition at line 266 of file PDSS_HKFT.h.

◆ m_c1

double m_c1 = 0.0
private

Input c1 coefficient (cal gmol-1 K-1)

Definition at line 269 of file PDSS_HKFT.h.

◆ m_c2

double m_c2 = 0.0
private

Input c2 coefficient (cal K gmol-1)

Definition at line 272 of file PDSS_HKFT.h.

◆ m_omega_pr_tr

double m_omega_pr_tr = 0.0
private

Input omega_pr_tr coefficient(cal gmol-1)

Definition at line 275 of file PDSS_HKFT.h.

◆ m_Y_pr_tr

double m_Y_pr_tr = 0.0
private

y = dZdT = 1/(esp*esp) desp/dT at 298.15 and 1 bar

Definition at line 278 of file PDSS_HKFT.h.

◆ m_Z_pr_tr

double m_Z_pr_tr = 0.0
private

Z = -1 / relEpsilon at 298.15 and 1 bar.

Definition at line 281 of file PDSS_HKFT.h.

◆ m_presR_bar

double m_presR_bar = OneAtm * 1.0E-5
private

Reference pressure is 1 atm in units of bar= 1.0132.

Definition at line 284 of file PDSS_HKFT.h.

◆ m_domega_jdT_prtr

double m_domega_jdT_prtr = 0.0
private

small value that is not quite zero

Definition at line 287 of file PDSS_HKFT.h.

◆ m_charge_j

double m_charge_j = 0.0
private

Charge of the ion.

Definition at line 290 of file PDSS_HKFT.h.

◆ s_InputInconsistencyErrorExit

int s_InputInconsistencyErrorExit = 1
staticprivate

Static variable determining error exiting.

If true, then will error exit if there is an inconsistency in DG0, DH0, and DS0. If not, then will rewrite DH0 to be consistent with the other two.

Definition at line 297 of file PDSS_HKFT.h.


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