Cantera  2.5.1
Public Member Functions | List of all members

Class DebyeHuckel represents a dilute liquid electrolyte phase which obeys the Debye Huckel formulation for nonideality. More...

#include <DebyeHuckel.h>

Inheritance diagram for DebyeHuckel:
[legend]
Collaboration diagram for DebyeHuckel:
[legend]

Public Member Functions

 DebyeHuckel ()
 Default Constructor. More...
 
 DebyeHuckel (const std::string &inputFile, const std::string &id="")
 Full constructor for creating the phase. More...
 
 DebyeHuckel (XML_Node &phaseRef, const std::string &id="")
 Full constructor for creating the phase. More...
 
Utilities
virtual std::string type () const
 String indicating the thermodynamic model implemented. More...
 
Molar Thermodynamic Properties of the Solution
virtual doublereal enthalpy_mole () const
 Molar enthalpy. Units: J/kmol. More...
 
virtual doublereal entropy_mole () const
 Molar entropy. Units: J/kmol/K. More...
 
virtual doublereal gibbs_mole () const
 Molar Gibbs function. Units: J/kmol. More...
 
virtual doublereal cp_mole () const
 Molar heat capacity at constant pressure. Units: J/kmol/K. More...
 
Activities, Standard States, and Activity Concentrations

The activity \(a_k\) of a species in solution is related to the chemical potential by

\[ \mu_k = \mu_k^0(T) + \hat R T \log a_k. \]

The quantity \(\mu_k^0(T,P)\) is the chemical potential at unit activity, which depends only on temperature and the pressure. Activity is assumed to be molality-based here.

virtual void getActivityConcentrations (doublereal *c) const
 This method returns an array of generalized concentrations. More...
 
virtual doublereal standardConcentration (size_t k=0) const
 Return the standard concentration for the kth species. More...
 
virtual void getActivities (doublereal *ac) const
 Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration. More...
 
virtual void getMolalityActivityCoefficients (doublereal *acMolality) const
 Get the array of non-dimensional molality-based activity coefficients at the current solution temperature, pressure, and solution concentration. More...
 
- Public Member Functions inherited from MolalityVPSSTP
 MolalityVPSSTP ()
 Default Constructor. More...
 
virtual std::string phaseOfMatter () const
 String indicating the mechanical phase of the matter in this Phase. More...
 
void setpHScale (const int pHscaleType)
 Set the pH scale, which determines the scale for single-ion activity coefficients. More...
 
int pHScale () const
 Reports the pH scale, which determines the scale for single-ion activity coefficients. More...
 
void setMoleFSolventMin (doublereal xmolSolventMIN)
 Sets the minimum mole fraction in the molality formulation. More...
 
doublereal moleFSolventMin () const
 Returns the minimum mole fraction in the molality formulation. More...
 
void calcMolalities () const
 Calculates the molality of all species and stores the result internally. More...
 
void getMolalities (doublereal *const molal) const
 This function will return the molalities of the species. More...
 
void setMolalities (const doublereal *const molal)
 Set the molalities of the solutes in a phase. More...
 
void setMolalitiesByName (const compositionMap &xMap)
 Set the molalities of a phase. More...
 
void setMolalitiesByName (const std::string &name)
 Set the molalities of a phase. More...
 
int activityConvention () const
 We set the convention to molality here. More...
 
virtual void getActivityCoefficients (doublereal *ac) const
 Get the array of non-dimensional activity coefficients at the current solution temperature, pressure, and solution concentration. More...
 
virtual double osmoticCoefficient () const
 Calculate the osmotic coefficient. More...
 
virtual void setStateFromXML (const XML_Node &state)
 Set equation of state parameter values from XML entries. More...
 
void setState_TPM (doublereal t, doublereal p, const doublereal *const molalities)
 Set the temperature (K), pressure (Pa), and molalities (gmol kg-1) of the solutes. More...
 
void setState_TPM (doublereal t, doublereal p, const compositionMap &m)
 Set the temperature (K), pressure (Pa), and molalities. More...
 
void setState_TPM (doublereal t, doublereal p, const std::string &m)
 Set the temperature (K), pressure (Pa), and molalities. More...
 
virtual void setState (const AnyMap &state)
 Set the state using an AnyMap containing any combination of properties supported by the thermodynamic model. More...
 
virtual void getdlnActCoeffdlnN (const size_t ld, doublereal *const dlnActCoeffdlnN)
 Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers. More...
 
virtual std::string report (bool show_thermo=true, doublereal threshold=1e-14) const
 returns a summary of the state of the phase as a string More...
 
- Public Member Functions inherited from VPStandardStateTP
 VPStandardStateTP ()
 Constructor. More...
 
virtual bool isCompressible () const
 Return whether phase represents a compressible substance. More...
 
virtual int standardStateConvention () const
 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...
 
virtual void getdlnActCoeffdlnN_diag (doublereal *dlnActCoeffdlnN_diag) const
 Get the array of log species mole number derivatives of the log activity coefficients. More...
 
virtual void getChemPotentials_RT (doublereal *mu) const
 Get the array of non-dimensional species chemical potentials. More...
 
void installPDSS (size_t k, std::unique_ptr< PDSS > &&pdss)
 Install a PDSS object for species k More...
 
PDSSprovidePDSS (size_t k)
 
const PDSSprovidePDSS (size_t k) const
 
virtual bool addSpecies (shared_ptr< Species > spec)
 Add a Species to this Phase. More...
 
virtual void getStandardChemPotentials (doublereal *mu) const
 Get the array of chemical potentials at unit activity for the species at their standard states at the current T and P of the solution. More...
 
virtual void getEnthalpy_RT (doublereal *hrt) const
 Get the nondimensional Enthalpy functions for the species at their standard states at the current T and P of the solution. More...
 
virtual void getEntropy_R (doublereal *sr) const
 Get the array of nondimensional Entropy functions for the standard state species at the current T and P of the solution. More...
 
virtual void getGibbs_RT (doublereal *grt) const
 Get the nondimensional Gibbs functions for the species in their standard states at the current T and P of the solution. More...
 
virtual void getPureGibbs (doublereal *gpure) const
 Get the Gibbs functions for the standard state of the species at the current T and P of the solution. More...
 
virtual void getIntEnergy_RT (doublereal *urt) const
 Returns the vector of nondimensional Internal Energies of the standard state species at the current T and P of the solution. More...
 
virtual void getCp_R (doublereal *cpr) const
 Get the nondimensional Heat Capacities at constant pressure for the species standard states at the current T and P of the solution. More...
 
virtual void getStandardVolumes (doublereal *vol) const
 Get the molar volumes of the species standard states at the current T and P of the solution. More...
 
virtual const vector_fpgetStandardVolumes () const
 
virtual void setTemperature (const doublereal temp)
 Set the temperature of the phase. More...
 
virtual void setPressure (doublereal p)
 Set the internally stored pressure (Pa) at constant temperature and composition. More...
 
virtual void setState_TP (doublereal T, doublereal pres)
 Set the temperature and pressure at the same time. More...
 
virtual doublereal pressure () const
 Returns the current pressure of the phase. More...
 
virtual void updateStandardStateThermo () const
 Updates the standard state thermodynamic functions at the current T and P of the solution. More...
 
virtual double minTemp (size_t k=npos) const
 Minimum temperature for which the thermodynamic data for the species or phase are valid. More...
 
virtual double maxTemp (size_t k=npos) const
 Maximum temperature for which the thermodynamic data for the species are valid. More...
 
virtual void getEnthalpy_RT_ref (doublereal *hrt) const
 
virtual void getGibbs_RT_ref (doublereal *grt) 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...
 
virtual void getGibbs_ref (doublereal *g) const
 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...
 
virtual void getEntropy_R_ref (doublereal *er) const
 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...
 
virtual void getCp_R_ref (doublereal *cprt) const
 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...
 
virtual void getStandardVolumes_ref (doublereal *vol) const
 Get the molar volumes of the species reference states at the current T and P_ref of the solution. More...
 
- Public Member Functions inherited from ThermoPhase
 ThermoPhase ()
 Constructor. More...
 
virtual doublereal refPressure () const
 Returns the reference pressure in Pa. More...
 
doublereal 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 doublereal 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...
 
bool chargeNeutralityNecessary () const
 Returns the chargeNeutralityNecessity boolean. More...
 
virtual doublereal intEnergy_mole () const
 Molar internal energy. Units: J/kmol. More...
 
virtual doublereal cv_mole () const
 Molar heat capacity at constant volume. Units: J/kmol/K. More...
 
virtual doublereal isothermalCompressibility () const
 Returns the isothermal compressibility. Units: 1/Pa. More...
 
virtual doublereal thermalExpansionCoeff () const
 Return the volumetric thermal expansion coefficient. Units: 1/K. More...
 
void setElectricPotential (doublereal v)
 Set the electric potential of this phase (V). More...
 
doublereal electricPotential () const
 Returns the electric potential of this phase (V). More...
 
virtual Units standardConcentrationUnits () const
 Returns the units of the "standard concentration" for this phase. More...
 
virtual doublereal logStandardConc (size_t k=0) const
 Natural logarithm of the standard concentration of the kth species. More...
 
virtual void getLnActivityCoefficients (doublereal *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 (doublereal *mu) const
 Get the species electrochemical potentials. More...
 
virtual void getPartialMolarIntEnergies (doublereal *ubar) const
 Return an array of partial molar internal energies for the species in the mixture. More...
 
virtual void getIntEnergy_RT_ref (doublereal *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...
 
doublereal enthalpy_mass () const
 Specific enthalpy. Units: J/kg. More...
 
doublereal intEnergy_mass () const
 Specific internal energy. Units: J/kg. More...
 
doublereal entropy_mass () const
 Specific entropy. Units: J/kg/K. More...
 
doublereal gibbs_mass () const
 Specific Gibbs function. Units: J/kg. More...
 
doublereal cp_mass () const
 Specific heat at constant pressure. Units: J/kg/K. More...
 
doublereal cv_mass () const
 Specific heat at constant volume. Units: J/kg/K. More...
 
doublereal RT () const
 Return the Gas Constant multiplied by the current temperature. More...
 
virtual void setState_TPX (doublereal t, doublereal p, const doublereal *x)
 Set the temperature (K), pressure (Pa), and mole fractions. More...
 
virtual void setState_TPX (doublereal t, doublereal p, const compositionMap &x)
 Set the temperature (K), pressure (Pa), and mole fractions. More...
 
virtual void setState_TPX (doublereal t, doublereal p, const std::string &x)
 Set the temperature (K), pressure (Pa), and mole fractions. More...
 
virtual void setState_TPY (doublereal t, doublereal p, const doublereal *y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. More...
 
virtual void setState_TPY (doublereal t, doublereal p, const compositionMap &y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. More...
 
virtual void setState_TPY (doublereal t, doublereal p, const std::string &y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. More...
 
virtual void setState_PX (doublereal p, doublereal *x)
 Set the pressure (Pa) and mole fractions. More...
 
virtual void setState_PY (doublereal p, doublereal *y)
 Set the internally stored pressure (Pa) and mass fractions. 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_RP (doublereal rho, doublereal p)
 Set the density (kg/m**3) and pressure (Pa) at constant composition. More...
 
virtual void setState_RPX (doublereal rho, doublereal p, const doublereal *x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions. More...
 
virtual void setState_RPX (doublereal rho, doublereal p, const compositionMap &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions. More...
 
virtual void setState_RPX (doublereal rho, doublereal p, const std::string &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions. More...
 
virtual void setState_RPY (doublereal rho, doublereal p, const doublereal *y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions. More...
 
virtual void setState_RPY (doublereal rho, doublereal p, const compositionMap &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions. More...
 
virtual void setState_RPY (doublereal rho, doublereal p, const std::string &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions. 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 std::string &fuelComp, const std::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 compositionMap &fuelComp, const compositionMap &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 std::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 std::string &fuelComp, const std::string &oxComp, ThermoBasis basis=ThermoBasis::molar, const std::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 compositionMap &fuelComp, const compositionMap &oxComp, ThermoBasis basis=ThermoBasis::molar, const std::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 std::string &fuelComp, const std::string &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the equivalence ratio. More...
 
void setEquivalenceRatio (double phi, const compositionMap &fuelComp, const compositionMap &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 std::string &fuelComp, const std::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 compositionMap &fuelComp, const compositionMap &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
 Compute the equivalence ratio for the current mixture from available oxygen and required oxygen. More...
 
void equilibrate (const std::string &XY, const std::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 doublereal *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 doublereal critTemperature () const
 Critical temperature (K). More...
 
virtual doublereal critPressure () const
 Critical pressure (Pa). More...
 
virtual doublereal critVolume () const
 Critical volume (m3/kmol). More...
 
virtual doublereal critCompressibility () const
 Critical compressibility (unitless). More...
 
virtual doublereal critDensity () const
 Critical density (kg/m3). More...
 
virtual doublereal satTemperature (doublereal p) const
 Return the saturation temperature given the pressure. More...
 
virtual doublereal satPressure (doublereal t)
 Return the saturation pressure given the temperature. More...
 
virtual doublereal vaporFraction () const
 Return the fraction of vapor at the current conditions. More...
 
virtual void setState_Tsat (doublereal t, doublereal x)
 Set the state to a saturated system at a particular temperature. More...
 
virtual void setState_Psat (doublereal p, doublereal 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...
 
virtual void modifySpecies (size_t k, shared_ptr< Species > spec)
 Modify the thermodynamic data associated with a species. More...
 
void saveSpeciesData (const size_t k, const XML_Node *const data)
 Store a reference pointer to the XML tree containing the species data for this phase. More...
 
const std::vector< const XML_Node * > & speciesData () const
 Return a pointer to the vector of XML nodes containing the species data for this phase. More...
 
virtual MultiSpeciesThermospeciesThermo (int k=-1)
 Return a changeable reference to the calculation manager for species reference-state thermodynamic properties. More...
 
virtual const MultiSpeciesThermospeciesThermo (int k=-1) const
 
virtual void initThermoFile (const std::string &inputFile, const std::string &id)
 
virtual void setParameters (int n, doublereal *const c)
 Set the equation of state parameters. More...
 
virtual void getParameters (int &n, doublereal *const c) const
 Get the equation of state parameters in a vector. More...
 
virtual void setParameters (const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap())
 Set equation of state parameters from an AnyMap phase description. More...
 
const AnyMapinput () const
 Access input data associated with the phase description. More...
 
AnyMapinput ()
 
virtual void setParametersFromXML (const XML_Node &eosdata)
 Set equation of state parameter values from XML entries. More...
 
virtual void getdlnActCoeffds (const doublereal dTds, const doublereal *const dXds, doublereal *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 (doublereal *dlnActCoeffdlnX_diag) const
 Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only. More...
 
virtual void getdlnActCoeffdlnN_numderiv (const size_t ld, doublereal *const dlnActCoeffdlnN)
 
virtual void reportCSV (std::ofstream &csvFile) const
 returns a summary of the state of the phase to a comma separated file. 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 std::string &fuelComp, const std::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 compositionMap &fuelComp, const compositionMap &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. More...
 
- Public Member Functions inherited from Phase
 Phase ()
 Default constructor. More...
 
 Phase (const Phase &)=delete
 
Phaseoperator= (const Phase &)=delete
 
XML_Nodexml () const
 Returns a const reference to the XML_Node that describes the phase. More...
 
void setXMLdata (XML_Node &xmlPhase)
 Stores the XML tree information for the current phase. More...
 
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 std::map< std::string, size_t > nativeState () const
 Return a map of properties defining the native state of a substance. More...
 
virtual std::vector< std::string > fullStates () const
 Return a vector containing full states defining a phase. More...
 
virtual std::vector< std::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_fp &state) const
 Save the current internal state of the phase. More...
 
virtual void saveState (size_t lenstate, doublereal *state) const
 Write to array 'state' the current internal state. More...
 
void restoreState (const vector_fp &state)
 Restore a state saved on a previous call to saveState. More...
 
virtual void restoreState (size_t lenstate, const doublereal *state)
 Restore the state of the phase from a previously saved state vector. More...
 
doublereal molecularWeight (size_t k) const
 Molecular weight of species k. More...
 
void getMolecularWeights (vector_fp &weights) const
 Copy the vector of molecular weights into vector weights. More...
 
void getMolecularWeights (doublereal *weights) const
 Copy the vector of molecular weights into array weights. More...
 
const vector_fpmolecularWeights () 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 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, e.g. More...
 
virtual void setRoot (std::shared_ptr< Solution > root)
 Set root Solution holding all phase information. More...
 
vector_fp getCompositionFromMap (const compositionMap &comp) const
 Converts a compositionMap 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...
 
std::string id () const
 Return the string id for the phase. More...
 
void setID (const std::string &id)
 Set the string id for the phase. More...
 
std::string name () const
 Return the name of the phase. More...
 
void setName (const std::string &nm)
 Sets the string name for the phase. More...
 
std::string elementName (size_t m) const
 Name of the element with index m. More...
 
size_t elementIndex (const std::string &name) const
 Return the index of element named 'name'. More...
 
const std::vector< std::string > & elementNames () const
 Return a read-only reference to the vector of element names. More...
 
doublereal atomicWeight (size_t m) const
 Atomic weight of element m. More...
 
doublereal 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_fpatomicWeights () 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...
 
doublereal nAtoms (size_t k, size_t m) const
 Number of atoms of element m in species k. More...
 
void getAtoms (size_t k, double *atomArray) const
 Get a vector containing the atomic composition of species k. More...
 
size_t speciesIndex (const std::string &name) const
 Returns the index of a species named 'name' within the Phase object. More...
 
std::string speciesName (size_t k) const
 Name of the species with index k. More...
 
std::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. More...
 
const std::vector< std::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 compositionMap &xMap)
 Set the species mole fractions by name. More...
 
void setMoleFractionsByName (const std::string &x)
 Set the mole fractions of a group of species by name. More...
 
void setMassFractionsByName (const compositionMap &yMap)
 Set the species mass fractions by name. More...
 
void setMassFractionsByName (const std::string &x)
 Set the species mass fractions by name. More...
 
void setState_TRX (doublereal t, doublereal dens, const doublereal *x)
 Set the internally stored temperature (K), density, and mole fractions. More...
 
void setState_TRX (doublereal t, doublereal dens, const compositionMap &x)
 Set the internally stored temperature (K), density, and mole fractions. More...
 
void setState_TRY (doublereal t, doublereal dens, const doublereal *y)
 Set the internally stored temperature (K), density, and mass fractions. More...
 
void setState_TRY (doublereal t, doublereal dens, const compositionMap &y)
 Set the internally stored temperature (K), density, and mass fractions. More...
 
void setState_TNX (doublereal t, doublereal n, const doublereal *x)
 Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions. More...
 
void setState_TR (doublereal t, doublereal rho)
 Set the internally stored temperature (K) and density (kg/m^3) More...
 
void setState_TX (doublereal t, doublereal *x)
 Set the internally stored temperature (K) and mole fractions. More...
 
void setState_TY (doublereal t, doublereal *y)
 Set the internally stored temperature (K) and mass fractions. More...
 
void setState_RX (doublereal rho, doublereal *x)
 Set the density (kg/m^3) and mole fractions. More...
 
void setState_RY (doublereal rho, doublereal *y)
 Set the density (kg/m^3) and mass fractions. More...
 
compositionMap 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 std::string &name) const
 Return the mole fraction of a single species. More...
 
compositionMap 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 std::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...
 
void getConcentrations (double *const c) const
 Get the species concentrations (kmol/m^3). More...
 
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...
 
doublereal elementalMassFraction (const size_t m) const
 Elemental mass fraction of element m. More...
 
doublereal elementalMoleFraction (const size_t m) const
 Elemental mole fraction of element m. More...
 
const double * moleFractdivMMW () const
 Returns a const pointer to the start of the moleFraction/MW array. More...
 
doublereal 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...
 
doublereal 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...
 
doublereal temperature () const
 Temperature (K). More...
 
virtual double density () const
 Density (kg/m^3). More...
 
double molarDensity () const
 Molar density (kmol/m^3). More...
 
double molarVolume () const
 Molar volume (m^3/kmol). More...
 
doublereal mean_X (const doublereal *const Q) const
 Evaluate the mole-fraction-weighted mean of an array Q. More...
 
doublereal mean_X (const vector_fp &Q) const
 Evaluate the mole-fraction-weighted mean of an array Q. More...
 
doublereal meanMolecularWeight () const
 The mean molecular weight. Units: (kg/kmol) More...
 
doublereal sum_xlogx () const
 Evaluate \( \sum_k X_k \log X_k \). More...
 
size_t addElement (const std::string &symbol, doublereal weight=-12345.0, int atomicNumber=0, doublereal entropy298=ENTROPY298_UNKNOWN, int elem_type=CT_ELEM_TYPE_ABSPOS)
 Add an element. More...
 
void addSpeciesAlias (const std::string &name, const std::string &alias)
 Add a species alias (i.e. More...
 
virtual std::vector< std::string > findIsomers (const compositionMap &compMap) const
 Return a vector with isomers names matching a given composition map. More...
 
virtual std::vector< std::string > findIsomers (const std::string &comp) const
 Return a vector with isomers names matching a given composition string. More...
 
shared_ptr< Speciesspecies (const std::string &name) const
 Return the Species object for the named species. More...
 
shared_ptr< Speciesspecies (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...
 

Partial Molar Properties of the Solution

int m_formDH
 form of the Debye-Huckel parameterization used in the model. More...
 
vector_int m_electrolyteSpeciesType
 Vector containing the electrolyte species type. More...
 
vector_fp m_Aionic
 a_k = Size of the ionic species in the DH formulation. units = meters More...
 
double m_IionicMolality
 Current value of the ionic strength on the molality scale. More...
 
double m_maxIionicStrength
 Maximum value of the ionic strength allowed in the calculation of the activity coefficients. More...
 
double m_IionicMolalityStoich
 Stoichiometric ionic strength on the molality scale. More...
 
double m_A_Debye
 Current value of the Debye Constant, A_Debye. More...
 
double m_B_Debye
 Current value of the constant that appears in the denominator. More...
 
vector_fp m_B_Dot
 Array of B_Dot values. More...
 
PDSS_Waterm_waterSS
 Pointer to the Water standard state object. More...
 
double m_densWaterSS
 Storage for the density of water's standard state. More...
 
std::unique_ptr< WaterPropsm_waterProps
 Pointer to the water property calculator. More...
 
vector_fp m_tmpV
 vector of size m_kk, used as a temporary holding area. More...
 
vector_fp m_speciesCharge_Stoich
 Stoichiometric species charge -> This is for calculations of the ionic strength which ignore ion-ion pairing into neutral molecules. More...
 
Array2D m_Beta_ij
 Array of 2D data used in the DHFORM_BETAIJ formulation Beta_ij.value(i,j) is the coefficient of the jth species for the specification of the chemical potential of the ith species. More...
 
vector_fp m_lnActCoeffMolal
 Logarithm of the activity coefficients on the molality scale. More...
 
vector_fp m_dlnActCoeffMolaldT
 Derivative of log act coeff wrt T. More...
 
vector_fp m_d2lnActCoeffMolaldT2
 2nd Derivative of log act coeff wrt T More...
 
vector_fp m_dlnActCoeffMolaldP
 Derivative of log act coeff wrt P. More...
 
bool m_useHelgesonFixedForm
 If true, then the fixed for of Helgeson's activity for water is used instead of the rigorous form obtained from Gibbs-Duhem relation. More...
 
int m_form_A_Debye
 Form of the constant outside the Debye-Huckel term called A. More...
 
virtual void getChemPotentials (doublereal *mu) const
 Get the species chemical potentials. Units: J/kmol. More...
 
virtual void getPartialMolarEnthalpies (doublereal *hbar) const
 Returns an array of partial molar enthalpies for the species in the mixture. More...
 
virtual void getPartialMolarEntropies (doublereal *sbar) const
 Returns an array of partial molar entropies of the species in the solution. More...
 
virtual void getPartialMolarCp (doublereal *cpbar) const
 Return an array of partial molar heat capacities for the species in the mixture. More...
 
virtual void getPartialMolarVolumes (doublereal *vbar) const
 Return an array of partial molar volumes for the species in the mixture. More...
 
virtual bool addSpecies (shared_ptr< Species > spec)
 
virtual void initThermo ()
 
virtual void initThermoXML (XML_Node &phaseNode, const std::string &id)
 Import and initialize a ThermoPhase object using an XML tree. More...
 
virtual double A_Debye_TP (double temperature=-1.0, double pressure=-1.0) const
 Return the Debye Huckel constant as a function of temperature and pressure (Units = sqrt(kg/gmol)) More...
 
virtual double dA_DebyedT_TP (double temperature=-1.0, double pressure=-1.0) const
 Value of the derivative of the Debye Huckel constant with respect to temperature. More...
 
virtual double d2A_DebyedT2_TP (double temperature=-1.0, double pressure=-1.0) const
 Value of the 2nd derivative of the Debye Huckel constant with respect to temperature as a function of temperature and pressure. More...
 
virtual double dA_DebyedP_TP (double temperature=-1.0, double pressure=-1.0) const
 Value of the derivative of the Debye Huckel constant with respect to pressure, as a function of temperature and pressure. More...
 
double AionicRadius (int k=0) const
 Reports the ionic radius of the kth species. More...
 
void setDebyeHuckelModel (const std::string &form)
 Set the DebyeHuckel parameterization form. More...
 
int formDH () const
 Returns the form of the Debye-Huckel parameterization used. More...
 
void setA_Debye (double A)
 Set the A_Debye parameter. More...
 
void setB_Debye (double B)
 
void setB_dot (double bdot)
 
void setMaxIonicStrength (double Imax)
 
void useHelgesonFixedForm (bool mode=true)
 
void setDefaultIonicRadius (double value)
 Set the default ionic radius [m] for each species. More...
 
void setBeta (const std::string &sp1, const std::string &sp2, double value)
 Set the value for the beta interaction between species sp1 and sp2. More...
 
Array2Dget_Beta_ij ()
 Returns a reference to M_Beta_ij. More...
 
static double _nonpolarActCoeff (double IionicMolality)
 Static function that implements the non-polar species salt-out modifications. More...
 
double _osmoticCoeffHelgesonFixedForm () const
 Formula for the osmotic coefficient that occurs in the GWB. More...
 
double _lnactivityWaterHelgesonFixedForm () const
 Formula for the log of the water activity that occurs in the GWB. More...
 
void s_update_lnMolalityActCoeff () const
 Calculate the log activity coefficients. More...
 
void s_update_dlnMolalityActCoeff_dT () const
 Calculation of temperature derivative of activity coefficient. More...
 
void s_update_d2lnMolalityActCoeff_dT2 () const
 Calculate the temperature 2nd derivative of the activity coefficient. More...
 
void s_update_dlnMolalityActCoeff_dP () const
 Calculate the pressure derivative of the activity coefficient. More...
 

Mechanical Equation of State Properties

In this equation of state implementation, the density is a function only of the mole fractions.

Therefore, it can't be an independent variable. Instead, the pressure is used as the independent variable. Functions which try to set the thermodynamic state by calling setDensity() will cause an exception to be thrown.

virtual void setDensity (const doublereal rho)
 Set the internally stored density (gm/m^3) of the phase. More...
 
virtual void setMolarDensity (const doublereal conc)
 Set the internally stored molar density (kmol/m^3) of the phase. More...
 
virtual void calcDensity ()
 Calculate the density of the mixture using the partial molar volumes and mole fractions as input. More...
 

Additional Inherited Members

- Protected Member Functions inherited from MolalityVPSSTP
virtual void getCsvReportData (std::vector< std::string > &names, std::vector< vector_fp > &data) const
 Fills names and data with the column names and species thermo properties to be included in the output of the reportCSV method. More...
 
virtual void getUnscaledMolalityActivityCoefficients (doublereal *acMolality) const
 Get the array of unscaled non-dimensional molality based activity coefficients at the current solution temperature, pressure, and solution concentration. More...
 
virtual void applyphScale (doublereal *acMolality) const
 Apply the current phScale to a set of activity Coefficients or activities. More...
 
- Protected Member Functions inherited from VPStandardStateTP
virtual void invalidateCache ()
 Invalidate any cached values which are normally updated only when a change in state is detected. More...
 
virtual void _updateStandardStateThermo () const
 Updates the standard state thermodynamic functions at the current T and P of the solution. More...
 
const vector_fpGibbs_RT_ref () const
 
- Protected Member Functions inherited from ThermoPhase
- Protected Member Functions inherited from Phase
void assertCompressible (const std::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...
 
virtual void compositionChanged ()
 Apply changes to the state which are needed after the composition changes. More...
 
- Protected Attributes inherited from MolalityVPSSTP
int m_pHScalingType
 Scaling to be used for output of single-ion species activity coefficients. More...
 
size_t m_indexCLM
 Index of the phScale species. More...
 
doublereal m_weightSolvent
 Molecular weight of the Solvent. More...
 
doublereal m_xmolSolventMIN
 
doublereal m_Mnaught
 This is the multiplication factor that goes inside log expressions involving the molalities of species. More...
 
vector_fp m_molalities
 Current value of the molalities of the species in the phase. More...
 
- Protected Attributes inherited from VPStandardStateTP
doublereal m_Pcurrent
 Current value of the pressure - state variable. More...
 
double m_minTemp
 The minimum temperature at which data for all species is valid. More...
 
double m_maxTemp
 The maximum temperature at which data for all species is valid. More...
 
doublereal m_Tlast_ss
 The last temperature at which the standard state thermodynamic properties were calculated at. More...
 
doublereal m_Plast_ss
 The last pressure at which the Standard State thermodynamic properties were calculated at. More...
 
std::vector< std::unique_ptr< PDSS > > m_PDSS_storage
 Storage for the PDSS objects for the species. More...
 
vector_fp m_h0_RT
 Vector containing the species reference enthalpies at T = m_tlast and P = p_ref. More...
 
vector_fp m_cp0_R
 Vector containing the species reference constant pressure heat capacities at T = m_tlast and P = p_ref. More...
 
vector_fp m_g0_RT
 Vector containing the species reference Gibbs functions at T = m_tlast and P = p_ref. More...
 
vector_fp m_s0_R
 Vector containing the species reference entropies at T = m_tlast and P = p_ref. More...
 
vector_fp m_V0
 Vector containing the species reference molar volumes. More...
 
vector_fp m_hss_RT
 Vector containing the species Standard State enthalpies at T = m_tlast and P = m_plast. More...
 
vector_fp m_cpss_R
 Vector containing the species Standard State constant pressure heat capacities at T = m_tlast and P = m_plast. More...
 
vector_fp m_gss_RT
 Vector containing the species Standard State Gibbs functions at T = m_tlast and P = m_plast. More...
 
vector_fp m_sss_R
 Vector containing the species Standard State entropies at T = m_tlast and P = m_plast. More...
 
vector_fp m_Vss
 Vector containing the species standard state volumes at T = m_tlast and P = m_plast. 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...
 
std::vector< const XML_Node * > m_speciesData
 Vector of pointers to the species databases. More...
 
doublereal m_phi
 Stored value of the electric potential for this phase. Units are Volts. More...
 
bool m_chargeNeutralityNecessary
 Boolean indicating whether a charge neutrality condition is a necessity. More...
 
int m_ssConvention
 Contains the standard state convention. More...
 
doublereal m_tlast
 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
 Number of species in the phase. More...
 
size_t m_ndim
 Dimensionality of the phase. More...
 
vector_fp m_speciesComp
 Atomic composition of the species. More...
 
vector_fp m_speciesCharge
 Vector of species charges. length m_kk. More...
 
std::map< std::string, shared_ptr< Species > > m_species
 
UndefElement::behavior m_undefinedElementBehavior
 Flag determining behavior when adding species with an undefined element. More...
 
bool m_caseSensitiveSpecies
 Flag determining whether case sensitive species names are enforced. More...
 

Detailed Description

Class DebyeHuckel represents a dilute liquid electrolyte phase which obeys the Debye Huckel formulation for nonideality.

The concentrations of the ionic species are assumed to obey the electroneutrality condition.

Specification of Species Standard State Properties

The standard states are on the unit molality basis. Therefore, in the documentation below, the normal \( o \) superscript is replaced with the \( \triangle \) symbol. The reference state symbol is now \( \triangle, ref \).

It is assumed that the reference state thermodynamics may be obtained by a pointer to a populated species thermodynamic property manager class (see ThermoPhase::m_spthermo). How to relate pressure changes to the reference state thermodynamics is resolved at this level.

For an incompressible, stoichiometric substance, the molar internal energy is independent of pressure. Since the thermodynamic properties are specified by giving the standard-state enthalpy, the term \( P_0 \hat v\) is subtracted from the specified molar enthalpy to compute the molar internal energy. The entropy is assumed to be independent of the pressure.

The enthalpy function is given by the following relation.

  \f[
    h^\triangle_k(T,P) = h^{\triangle,ref}_k(T)
    + \tilde v \left( P - P_{ref} \right)
  \f]

For an incompressible, stoichiometric substance, the molar internal energy is independent of pressure. Since the thermodynamic properties are specified by giving the standard-state enthalpy, the term \( P_{ref} \tilde v\) is subtracted from the specified reference molar enthalpy to compute the molar internal energy.

  \f[
       u^\triangle_k(T,P) = h^{\triangle,ref}_k(T) - P_{ref} \tilde v
  \f]

The standard state heat capacity and entropy are independent of pressure. The standard state Gibbs free energy is obtained from the enthalpy and entropy functions.

The current model assumes that an incompressible molar volume for all solutes. The molar volume for the water solvent, however, is obtained from a pure water equation of state, waterSS. Therefore, the water standard state varies with both T and P. It is an error to request standard state water properties at a T and P where the water phase is not a stable phase, i.e., beyond its spinodal curve.

Specification of Solution Thermodynamic Properties

Chemical potentials of the solutes, \( \mu_k \), and the solvent, \( \mu_o \), which are based on the molality form, have the following general format:

\[ \mu_k = \mu^{\triangle}_k(T,P) + R T ln(\gamma_k^{\triangle} \frac{m_k}{m^\triangle}) \]

\[ \mu_o = \mu^o_o(T,P) + RT ln(a_o) \]

where \( \gamma_k^{\triangle} \) is the molality based activity coefficient for species \(k\).

Individual activity coefficients of ions can not be independently measured. Instead, only binary pairs forming electroneutral solutions can be measured.

Ionic Strength

Most of the parameterizations within the model use the ionic strength as a key variable. The ionic strength, \( I\) is defined as follows

\[ I = \frac{1}{2} \sum_k{m_k z_k^2} \]

\( m_k \) is the molality of the kth species. \( z_k \) is the charge of the kth species. Note, the ionic strength is a defined units quantity. The molality has defined units of gmol kg-1, and therefore the ionic strength has units of sqrt( gmol kg-1).

In some instances, from some authors, a different formulation is used for the ionic strength in the equations below. The different formulation is due to the possibility of the existence of weak acids and how association wrt to the weak acid equilibrium relation affects the calculation of the activity coefficients via the assumed value of the ionic strength.

If we are to assume that the association reaction doesn't have an effect on the ionic strength, then we will want to consider the associated weak acid as in effect being fully dissociated, when we calculate an effective value for the ionic strength. We will call this calculated value, the stoichiometric ionic strength, \( I_s \), putting a subscript s to denote it from the more straightforward calculation of \( I \).

\[ I_s = \frac{1}{2} \sum_k{m_k^s z_k^2} \]

Here, \( m_k^s \) is the value of the molalities calculated assuming that all weak acid-base pairs are in their fully dissociated states. This calculation may be simplified by considering that the weakly associated acid may be made up of two charged species, k1 and k2, each with their own charges, obeying the following relationship:

\[ z_k = z_{k1} + z_{k2} \]

Then, we may only need to specify one charge value, say, \( z_{k1}\), the cation charge number, in order to get both numbers, since we have already specified \( z_k \) in the definition of original species. Then, the stoichiometric ionic strength may be calculated via the following formula.

\[ I_s = \frac{1}{2} \left(\sum_{k,ions}{m_k z_k^2}+ \sum_{k,weak_assoc}(m_k z_{k1}^2 + m_k z_{k2}^2) \right) \]

The specification of which species are weakly associated acids is made in the input file via the stoichIsMods XML block, where the charge for k1 is also specified. An example is given below:

<stoichIsMods>
NaCl(aq):-1.0
</stoichIsMods>

Because we need the concept of a weakly associated acid in order to calculate \( I_s \) we need to catalog all species in the phase. This is done using the following categories:

Polar and non-polar neutral species are differentiated, because some additions to the activity coefficient expressions distinguish between these two types of solutes. This is the so-called salt-out effect.

The type of species is specified in the electrolyteSpeciesType XML block. Note, this is not considered a part of the specification of the standard state for the species, at this time. Therefore, this information is put under the activityCoefficient XML block. An example is given below

<electrolyteSpeciesType>
H2L(L):solvent
H+:chargedSpecies
NaOH(aq):weakAcidAssociated
NaCl(aq):strongAcidAssociated
NH3(aq):polarNeutral
O2(aq):nonpolarNeutral
</electrolyteSpeciesType>

Much of the species electrolyte type information is inferred from other information in the input file. For example, as species which is charged is given the "chargedSpecies" default category. A neutral solute species is put into the "nonpolarNeutral" category by default.

The specification of solute activity coefficients depends on the model assumed for the Debye-Huckel term. The model is set by the internal parameter m_formDH. We will now describe each category in its own section.

Debye-Huckel Dilute Limit

DHFORM_DILUTE_LIMIT = 0

This form assumes a dilute limit to DH, and is mainly for informational purposes:

\[ \ln(\gamma_k^\triangle) = - z_k^2 A_{Debye} \sqrt{I} \]

where \( I\) is the ionic strength

\[ I = \frac{1}{2} \sum_k{m_k z_k^2} \]

The activity for the solvent water, \( a_o \), is not independent and must be determined from the Gibbs-Duhem relation.

\[ \ln(a_o) = \frac{X_o - 1.0}{X_o} + \frac{ 2 A_{Debye} \tilde{M}_o}{3} (I)^{3/2} \]

Bdot Formulation

DHFORM_BDOT_AK = 1

This form assumes Bethke's format for the Debye Huckel activity coefficient:

\[ \ln(\gamma_k^\triangle) = -z_k^2 \frac{A_{Debye} \sqrt{I}}{ 1 + B_{Debye} a_k \sqrt{I}} + \log(10) B^{dot}_k I \]

Note, this particular form where \( a_k \) can differ in multielectrolyte solutions has problems with respect to a Gibbs-Duhem analysis. However, we include it here because there is a lot of data fit to it.

The activity for the solvent water, \( a_o \), is not independent and must be determined from the Gibbs-Duhem relation. Here, we use:

\[ \ln(a_o) = \frac{X_o - 1.0}{X_o} + \frac{ 2 A_{Debye} \tilde{M}_o}{3} (I)^{1/2} \left[ \sum_k{\frac{1}{2} m_k z_k^2 \sigma( B_{Debye} a_k \sqrt{I} ) } \right] - \frac{\log(10)}{2} \tilde{M}_o I \sum_k{ B^{dot}_k m_k} \]

where

\[ \sigma (y) = \frac{3}{y^3} \left[ (1+y) - 2 \ln(1 + y) - \frac{1}{1+y} \right] \]

Additionally, Helgeson's formulation for the water activity is offered as an alternative.

Bdot Formulation with Uniform Size Parameter in the Denominator

DHFORM_BDOT_AUNIFORM = 2

This form assumes Bethke's format for the Debye-Huckel activity coefficient

\[ \ln(\gamma_k^\triangle) = -z_k^2 \frac{A_{Debye} \sqrt{I}}{ 1 + B_{Debye} a \sqrt{I}} + \log(10) B^{dot}_k I \]

The value of a is determined at the beginning of the calculation, and not changed.

\[ \ln(a_o) = \frac{X_o - 1.0}{X_o} + \frac{ 2 A_{Debye} \tilde{M}_o}{3} (I)^{3/2} \sigma( B_{Debye} a \sqrt{I} ) - \frac{\log(10)}{2} \tilde{M}_o I \sum_k{ B^{dot}_k m_k} \]

Beta_IJ formulation

DHFORM_BETAIJ        = 3

This form assumes a linear expansion in a virial coefficient form. It is used extensively in the book by Newmann, "Electrochemistry Systems", and is the beginning of more complex treatments for stronger electrolytes, fom Pitzer and from Harvey, Moller, and Weire.

\[ \ln(\gamma_k^\triangle) = -z_k^2 \frac{A_{Debye} \sqrt{I}}{ 1 + B_{Debye} a \sqrt{I}} + 2 \sum_j \beta_{j,k} m_j \]

In the current treatment the binary interaction coefficients, \( \beta_{j,k}\), are independent of temperature and pressure.

\[ \ln(a_o) = \frac{X_o - 1.0}{X_o} + \frac{ 2 A_{Debye} \tilde{M}_o}{3} (I)^{3/2} \sigma( B_{Debye} a \sqrt{I} ) - \tilde{M}_o \sum_j \sum_k \beta_{j,k} m_j m_k \]

In this formulation the ionic radius, \( a \), is a constant. This must be supplied to the model, in an ionicRadius XML block.

The \( \beta_{j,k} \) parameters are binary interaction parameters. They are supplied to the object in an DHBetaMatrix XML block. There are in principle \( N (N-1) /2 \) different, symmetric interaction parameters, where \( N \) are the number of solute species in the mechanism. An example is given below.

An example activityCoefficients XML block for this formulation is supplied below

<activityCoefficients model="Beta_ij">
<!-- A_Debye units = sqrt(kg/gmol) -->
<A_Debye> 1.172576 </A_Debye>
<!-- B_Debye units = sqrt(kg/gmol)/m -->
<B_Debye> 3.28640E9 </B_Debye>
<ionicRadius default="3.042843" units="Angstroms">
</ionicRadius>
<DHBetaMatrix>
H+:Cl-:0.27
Na+:Cl-:0.15
Na+:OH-:0.06
</DHBetaMatrix>
<stoichIsMods>
NaCl(aq):-1.0
</stoichIsMods>
<electrolyteSpeciesType>
H+:chargedSpecies
NaCl(aq):weakAcidAssociated
</electrolyteSpeciesType>
</activityCoefficients>

Pitzer Beta_IJ formulation

DHFORM_PITZER_BETAIJ  = 4

This form assumes an activity coefficient formulation consistent with a truncated form of Pitzer's formulation. Pitzer's formulation is equivalent to the formulations above in the dilute limit, where rigorous theory may be applied.

\[ \ln(\gamma_k^\triangle) = -z_k^2 \frac{A_{Debye}}{3} \frac{\sqrt{I}}{ 1 + B_{Debye} a \sqrt{I}} -2 z_k^2 \frac{A_{Debye}}{3} \frac{\ln(1 + B_{Debye} a \sqrt{I})}{ B_{Debye} a} + 2 \sum_j \beta_{j,k} m_j \]

\[ \ln(a_o) = \frac{X_o - 1.0}{X_o} + \frac{ 2 A_{Debye} \tilde{M}_o}{3} \frac{(I)^{3/2} }{1 + B_{Debye} a \sqrt{I} } - \tilde{M}_o \sum_j \sum_k \beta_{j,k} m_j m_k \]

Specification of the Debye Huckel Constants

In the equations above, the formulas for \( A_{Debye} \) and \( B_{Debye} \) are needed. The DebyeHuckel object uses two methods for specifying these quantities. The default method is to assume that \( A_{Debye} \) is a constant, given in the initialization process, and stored in the member double, m_A_Debye. Optionally, a full water treatment may be employed that makes \( A_{Debye} \) a full function of T and P.

\[ A_{Debye} = \frac{F e B_{Debye}}{8 \pi \epsilon R T} {\left( C_o \tilde{M}_o \right)}^{1/2} \]

where

\[ B_{Debye} = \frac{F} {{(\frac{\epsilon R T}{2})}^{1/2}} \]

Therefore:

\[ A_{Debye} = \frac{1}{8 \pi} {\left(\frac{2 N_a \rho_o}{1000}\right)}^{1/2} {\left(\frac{N_a e^2}{\epsilon R T }\right)}^{3/2} \]

where

Nominal value at 298 K and 1 atm = 1.172576 (kg/gmol)^(1/2) based on:

An example of a fixed value implementation is given below.

<activityCoefficients model="Beta_ij">
<!-- A_Debye units = sqrt(kg/gmol) -->
<A_Debye> 1.172576 </A_Debye>
<!-- B_Debye units = sqrt(kg/gmol)/m -->
<B_Debye> 3.28640E9 </B_Debye>
</activityCoefficients>

An example of a variable value implementation is given below.

<activityCoefficients model="Beta_ij">
<A_Debye model="water" />
<!-- B_Debye units = sqrt(kg/gmol)/m -->
<B_Debye> 3.28640E9 </B_Debye>
</activityCoefficients>

Currently, \( B_{Debye} \) is a constant in the model, specified either by a default water value, or through the input file. This may have to be looked at, in the future.

%Application within Kinetics Managers

For the time being, we have set the standard concentration for all species in this phase equal to the default concentration of the solvent at 298 K and 1 atm. This means that the kinetics operator essentially works on an activities basis, with units specified as if it were on a concentration basis.

For example, a bulk-phase binary reaction between liquid species j and k, producing a new liquid species l would have the following equation for its rate of progress variable, \( R^1 \), which has units of kmol m-3 s-1.

\[ R^1 = k^1 C_j^a C_k^a = k^1 (C_o a_j) (C_o a_k) \]

where

\[ C_j^a = C_o a_j \quad and \quad C_k^a = C_o a_k \]

\( C_j^a \) is the activity concentration of species j, and \( C_k^a \) is the activity concentration of species k. \( C_o \) is the concentration of water at 298 K and 1 atm. \( a_j \) is the activity of species j at the current temperature and pressure and concentration of the liquid phase. \(k^1 \) has units of m3 kmol-1 s-1.

The reverse rate constant can then be obtained from the law of microscopic reversibility and the equilibrium expression for the system.

\[ \frac{a_j a_k}{ a_l} = K^{o,1} = \exp(\frac{\mu^o_l - \mu^o_j - \mu^o_k}{R T} ) \]

\( K^{o,1} \) is the dimensionless form of the equilibrium constant.

\[ R^{-1} = k^{-1} C_l^a = k^{-1} (C_o a_l) \]

where

\[ k^{-1} = k^1 K^{o,1} C_o \]

\(k^{-1} \) has units of s-1.

Note, this treatment may be modified in the future, as events dictate.

Instantiation of the Class

The constructor for this phase is NOT located in the default ThermoFactory for Cantera. However, a new DebyeHuckel object may be created by the following code snippets:

DebyeHuckel *DH = new DebyeHuckel("DH_NaCl.xml", "NaCl_electrolyte");
DebyeHuckel()
Default Constructor.
Definition: DebyeHuckel.cpp:29

or

XML_Node *xm = get_XML_NameID("phase", "DH_NaCl.xml#NaCl_electrolyte", 0);
DebyeHuckel *dh = new DebyeHuckel(*xm);
XML_Node * get_XML_NameID(const std::string &nameTarget, const std::string &file_ID, XML_Node *root)
This routine will locate an XML node in either the input XML tree or in another input file specified ...
Definition: global.cpp:232

or by the following call to importPhase():

XML_Node *xm = get_XML_NameID("phase", "DH_NaCl.xml#NaCl_electrolyte", 0);
DebyeHuckel dhphase;
importPhase(*xm, &dhphase);
void importPhase(XML_Node &phase, ThermoPhase *th)
Import a phase information into an empty ThermoPhase object.

XML Example

The phase model name for this is called StoichSubstance. It must be supplied as the model attribute of the thermo XML element entry. Within the phase XML block, the density of the phase must be specified. An example of an XML file this phase is given below.

<phase id="NaCl_electrolyte" dim="3">
<speciesArray datasrc="#species_waterSolution">
H2O(L) Na+ Cl- H+ OH- NaCl(aq) NaOH(aq)
</speciesArray>
<state>
<temperature units="K"> 300 </temperature>
<pressure units="Pa">101325.0</pressure>
<soluteMolalities>
Na+:3.0
Cl-:3.0
H+:1.0499E-8
OH-:1.3765E-6
NaCl(aq):0.98492
NaOH(aq):3.8836E-6
</soluteMolalities>
</state>
<!-- thermo model identifies the inherited class
from ThermoPhase that will handle the thermodynamics.
-->
<thermo model="DebyeHuckel">
<standardConc model="solvent_volume" />
<activityCoefficients model="Beta_ij">
<!-- A_Debye units = sqrt(kg/gmol) -->
<A_Debye> 1.172576 </A_Debye>
<!-- B_Debye units = sqrt(kg/gmol)/m -->
<B_Debye> 3.28640E9 </B_Debye>
<ionicRadius default="3.042843" units="Angstroms">
</ionicRadius>
<DHBetaMatrix>
H+:Cl-:0.27
Na+:Cl-:0.15
Na+:OH-:0.06
</DHBetaMatrix>
<stoichIsMods>
NaCl(aq):-1.0
</stoichIsMods>
<electrolyteSpeciesType>
H+:chargedSpecies
NaCl(aq):weakAcidAssociated
</electrolyteSpeciesType>
</activityCoefficients>
<solvent> H2O(L) </solvent>
</thermo>
<elementArray datasrc="elements.xml"> O H Na Cl </elementArray>
</phase>
doublereal temperature() const
Temperature (K).
Definition: Phase.h:667
ThermoPhase()
Constructor.
Definition: ThermoPhase.cpp:27
virtual doublereal pressure() const
Returns the current pressure of the phase.

Definition at line 558 of file DebyeHuckel.h.

Constructor & Destructor Documentation

◆ DebyeHuckel() [1/3]

Default Constructor.

Definition at line 29 of file DebyeHuckel.cpp.

◆ DebyeHuckel() [2/3]

DebyeHuckel ( const std::string &  inputFile,
const std::string &  id = "" 
)

Full constructor for creating the phase.

Parameters
inputFileFile name containing the definition of the phase
idid attribute containing the name of the phase.

Definition at line 43 of file DebyeHuckel.cpp.

◆ DebyeHuckel() [3/3]

DebyeHuckel ( XML_Node phaseRef,
const std::string &  id = "" 
)

Full constructor for creating the phase.

Parameters
phaseRefXML phase node containing the description of the phase
idid attribute containing the name of the phase.
Deprecated:
The XML input format is deprecated and will be removed in Cantera 3.0.

Definition at line 59 of file DebyeHuckel.cpp.

Member Function Documentation

◆ type()

virtual std::string type ( ) const
inlinevirtual

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 ThermoPhase.

Definition at line 586 of file DebyeHuckel.h.

◆ enthalpy_mole()

doublereal enthalpy_mole ( ) const
virtual

Molar enthalpy. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 80 of file DebyeHuckel.cpp.

References DebyeHuckel::getPartialMolarEnthalpies(), DebyeHuckel::m_tmpV, and Phase::mean_X().

◆ entropy_mole()

doublereal entropy_mole ( ) const
virtual

Molar entropy. Units: J/kmol/K.

For an ideal, constant partial molar volume solution mixture with pure species phases which exhibit zero volume expansivity:

\[ \hat s(T, P, X_k) = \sum_k X_k \hat s^0_k(T) - \hat R \sum_k X_k log(X_k) \]

The reference-state pure-species entropies \( \hat s^0_k(T,p_{ref}) \) are computed by the species thermodynamic property manager. The pure species entropies are independent of temperature since the volume expansivities are equal to zero.

See also
MultiSpeciesThermo

Reimplemented from ThermoPhase.

Definition at line 86 of file DebyeHuckel.cpp.

References DebyeHuckel::getPartialMolarEntropies(), DebyeHuckel::m_tmpV, and Phase::mean_X().

◆ gibbs_mole()

doublereal gibbs_mole ( ) const
virtual

Molar Gibbs function. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 92 of file DebyeHuckel.cpp.

References DebyeHuckel::getChemPotentials(), DebyeHuckel::m_tmpV, and Phase::mean_X().

◆ cp_mole()

doublereal cp_mole ( ) const
virtual

Molar heat capacity at constant pressure. Units: J/kmol/K.

Reimplemented from ThermoPhase.

Definition at line 98 of file DebyeHuckel.cpp.

References DebyeHuckel::getPartialMolarCp(), DebyeHuckel::m_tmpV, and Phase::mean_X().

◆ calcDensity()

void calcDensity ( )
protectedvirtual

Calculate the density of the mixture using the partial molar volumes and mole fractions as input.

The formula for this is

\[ \rho = \frac{\sum_k{X_k W_k}}{\sum_k{X_k V_k}} \]

where \(X_k\) are the mole fractions, \(W_k\) are the molecular weights, and \(V_k\) are the pure species molar volumes.

Note, the basis behind this formula is that in an ideal solution the partial molar volumes are equal to the pure species molar volumes. We have additionally specified in this class that the pure species molar volumes are independent of temperature and pressure.

NOTE: This function is not a member of the ThermoPhase base class.

Reimplemented from VPStandardStateTP.

Definition at line 106 of file DebyeHuckel.cpp.

References Phase::assignDensity(), PDSS_Water::density(), DebyeHuckel::getPartialMolarVolumes(), DebyeHuckel::m_densWaterSS, DebyeHuckel::m_tmpV, DebyeHuckel::m_waterSS, Phase::mean_X(), and Phase::meanMolecularWeight().

◆ setDensity()

void setDensity ( const doublereal  rho)
virtual

Set the internally stored density (gm/m^3) of the phase.

Overridden setDensity() function is necessary because the density is not an independent variable.

This function will now throw an error condition

May have to adjust the strategy here to make the eos for these materials slightly compressible, in order to create a condition where the density is a function of the pressure.

This function will now throw an error condition if the input isn't exactly equal to the current density.

Todo:
Now have a compressible ss equation for liquid water. Therefore, this phase is compressible. May still want to change the independent variable however.
Parameters
rhoInput density (kg/m^3).
Deprecated:
Functionality merged with base function after Cantera 2.5. (superseded by isCompressible check in Phase::setDensity)

Reimplemented from Phase.

Definition at line 118 of file DebyeHuckel.cpp.

References Phase::density(), and Cantera::warn_deprecated().

◆ setMolarDensity()

void setMolarDensity ( const doublereal  conc)
virtual

Set the internally stored molar density (kmol/m^3) of the phase.

Overridden setMolarDensity() function is necessary because the density is not an independent variable.

This function will now throw an error condition if the input isn't exactly equal to the current molar density.

Parameters
concInput molar density (kmol/m^3).
Deprecated:
Functionality merged with base function after Cantera 2.5. (superseded by isCompressible check in Phase::setDensity)

Reimplemented from Phase.

Definition at line 130 of file DebyeHuckel.cpp.

References Phase::molarDensity(), and Cantera::warn_deprecated().

◆ getActivityConcentrations()

void getActivityConcentrations ( doublereal *  c) const
virtual

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 MolalityVPSSTP.

Definition at line 144 of file DebyeHuckel.cpp.

References DebyeHuckel::getActivities(), Phase::m_kk, and DebyeHuckel::standardConcentration().

◆ standardConcentration()

doublereal standardConcentration ( size_t  k = 0) const
virtual

Return the standard concentration for the kth species.

The standard concentration \( C^0_k \) used to normalize the activity (i.e., generalized) concentration in kinetics calculations.

For the time being, we will use the concentration of pure solvent for the the standard concentration of all species. This has the effect of making reaction rates based on the molality of species proportional to the molality of the species.

Parameters
kOptional parameter indicating the species. The default is to assume this refers to species 0.
Returns
the standard Concentration in units of m^3/kmol

Reimplemented from MolalityVPSSTP.

Definition at line 153 of file DebyeHuckel.cpp.

References PDSS::molarVolume().

Referenced by DebyeHuckel::getActivityConcentrations().

◆ getActivities()

void getActivities ( doublereal *  ac) const
virtual

Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration.

(note solvent activity coefficient is on molar scale).

Parameters
acOutput vector of activities. Length: m_kk.

Reimplemented from MolalityVPSSTP.

Definition at line 159 of file DebyeHuckel.cpp.

References VPStandardStateTP::_updateStandardStateThermo(), Phase::m_kk, DebyeHuckel::m_lnActCoeffMolal, MolalityVPSSTP::m_molalities, Phase::moleFraction(), and DebyeHuckel::s_update_lnMolalityActCoeff().

Referenced by DebyeHuckel::getActivityConcentrations().

◆ getMolalityActivityCoefficients()

void getMolalityActivityCoefficients ( doublereal *  acMolality) const
virtual

Get the array of non-dimensional molality-based activity coefficients at the current solution temperature, pressure, and solution concentration.

note solvent is on molar scale. The solvent molar based activity coefficient is returned.

Note, most of the work is done in an internal private routine

Parameters
acMolalityVector of Molality-based activity coefficients Length: m_kk

Reimplemented from MolalityVPSSTP.

Definition at line 173 of file DebyeHuckel.cpp.

References VPStandardStateTP::_updateStandardStateThermo(), DebyeHuckel::A_Debye_TP(), Phase::m_kk, DebyeHuckel::m_lnActCoeffMolal, and DebyeHuckel::s_update_lnMolalityActCoeff().

◆ getChemPotentials()

void getChemPotentials ( doublereal *  mu) const
virtual

Get the species chemical potentials. Units: J/kmol.

This function returns a vector of chemical potentials of the species in solution.

\[ \mu_k = \mu^{\triangle}_k(T,P) + R T ln(\gamma_k^{\triangle} m_k) \]

Parameters
muOutput vector of species chemical potentials. Length: m_kk. Units: J/kmol

Reimplemented from ThermoPhase.

Definition at line 186 of file DebyeHuckel.cpp.

References VPStandardStateTP::getStandardChemPotentials(), Phase::m_kk, DebyeHuckel::m_lnActCoeffMolal, MolalityVPSSTP::m_molalities, Phase::moleFraction(), ThermoPhase::RT(), DebyeHuckel::s_update_lnMolalityActCoeff(), and Cantera::SmallNumber.

Referenced by DebyeHuckel::gibbs_mole().

◆ getPartialMolarEnthalpies()

void getPartialMolarEnthalpies ( doublereal *  hbar) const
virtual

Returns an array of partial molar enthalpies for the species in the mixture.

Units (J/kmol)

For this phase, the partial molar enthalpies are equal to the standard state enthalpies modified by the derivative of the molality-based activity coefficient wrt temperature

\[ \bar h_k(T,P) = h^{\triangle}_k(T,P) - R T^2 \frac{d \ln(\gamma_k^\triangle)}{dT} \]

The solvent partial molar enthalpy is equal to

\[ \bar h_o(T,P) = h^{o}_o(T,P) - R T^2 \frac{d \ln(a_o}{dT} \]

The temperature dependence of the activity coefficients currently only occurs through the temperature dependence of the Debye constant.

Parameters
hbarOutput vector of species partial molar enthalpies. Length: m_kk. units are J/kmol.

Reimplemented from ThermoPhase.

Definition at line 206 of file DebyeHuckel.cpp.

References DebyeHuckel::dA_DebyedT_TP(), VPStandardStateTP::getEnthalpy_RT(), DebyeHuckel::m_dlnActCoeffMolaldT, Phase::m_kk, ThermoPhase::RT(), DebyeHuckel::s_update_dlnMolalityActCoeff_dT(), DebyeHuckel::s_update_lnMolalityActCoeff(), and Phase::temperature().

Referenced by DebyeHuckel::enthalpy_mole().

◆ getPartialMolarEntropies()

void getPartialMolarEntropies ( doublereal *  sbar) const
virtual

Returns an array of partial molar entropies of the species in the solution.

Units: J/kmol/K. Maxwell's equations provide an insight in how to calculate this (p.215 Smith and Van Ness)

\[ \frac{d\mu_i}{dT} = -\bar{s}_i \]

For this phase, the partial molar entropies are equal to the SS species entropies plus the ideal solution contribution:

\[ \bar s_k(T,P) = \hat s^0_k(T) - R log(M0 * molality[k]) \]

\[ \bar s_{solvent}(T,P) = \hat s^0_{solvent}(T) - R ((xmolSolvent - 1.0) / xmolSolvent) \]

The reference-state pure-species entropies, \( \hat s^0_k(T) \), at the reference pressure, \( P_{ref} \), are computed by the species thermodynamic property manager. They are polynomial functions of temperature.

See also
MultiSpeciesThermo
Parameters
sbarOutput vector of species partial molar entropies. Length = m_kk. units are J/kmol/K.

Reimplemented from ThermoPhase.

Definition at line 231 of file DebyeHuckel.cpp.

References DebyeHuckel::dA_DebyedT_TP(), Cantera::GasConstant, VPStandardStateTP::getEntropy_R(), DebyeHuckel::m_dlnActCoeffMolaldT, Phase::m_kk, DebyeHuckel::m_lnActCoeffMolal, MolalityVPSSTP::m_molalities, Phase::moleFraction(), ThermoPhase::RT(), DebyeHuckel::s_update_dlnMolalityActCoeff_dT(), DebyeHuckel::s_update_lnMolalityActCoeff(), and Cantera::SmallNumber.

Referenced by DebyeHuckel::entropy_mole().

◆ getPartialMolarCp()

void getPartialMolarCp ( doublereal *  cpbar) const
virtual

Return an array of partial molar heat capacities for the species in the mixture.

Units: J/kmol/K

Parameters
cpbarOutput vector of species partial molar heat capacities at constant pressure. Length = m_kk. units are J/kmol/K.

Reimplemented from ThermoPhase.

Definition at line 281 of file DebyeHuckel.cpp.

References DebyeHuckel::dA_DebyedT_TP(), Cantera::GasConstant, VPStandardStateTP::getCp_R(), DebyeHuckel::m_d2lnActCoeffMolaldT2, DebyeHuckel::m_dlnActCoeffMolaldT, Phase::m_kk, ThermoPhase::RT(), DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2(), DebyeHuckel::s_update_dlnMolalityActCoeff_dT(), DebyeHuckel::s_update_lnMolalityActCoeff(), and Phase::temperature().

Referenced by DebyeHuckel::cp_mole().

◆ getPartialMolarVolumes()

void getPartialMolarVolumes ( doublereal *  vbar) const
virtual

Return an array of partial molar volumes for the species in the mixture.

Units: m^3/kmol.

For this solution, the partial molar volumes are normally equal to the constant species molar volumes, except when the activity coefficients depend on pressure.

The general relation is

  vbar_i = d(chemPot_i)/dP at const T, n
         = V0_i + d(Gex)/dP)_T,M
         = V0_i + RT d(lnActCoeffi)dP _T,M
Parameters
vbarOutput vector of species partial molar volumes. Length = m_kk. units are m^3/kmol.

Reimplemented from ThermoPhase.

Definition at line 269 of file DebyeHuckel.cpp.

References VPStandardStateTP::getStandardVolumes(), DebyeHuckel::m_dlnActCoeffMolaldP, Phase::m_kk, ThermoPhase::RT(), DebyeHuckel::s_update_dlnMolalityActCoeff_dP(), and DebyeHuckel::s_update_lnMolalityActCoeff().

Referenced by DebyeHuckel::calcDensity().

◆ addSpecies()

bool addSpecies ( shared_ptr< Species spec)
virtual

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

Reimplemented from MolalityVPSSTP.

Definition at line 738 of file DebyeHuckel.cpp.

References VPStandardStateTP::addSpecies(), Cantera::cEST_solvent, AnyMap::convert(), Cantera::interp_est(), DebyeHuckel::m_Aionic, DebyeHuckel::m_B_Dot, DebyeHuckel::m_d2lnActCoeffMolaldT2, DebyeHuckel::m_dlnActCoeffMolaldP, DebyeHuckel::m_dlnActCoeffMolaldT, DebyeHuckel::m_electrolyteSpeciesType, DebyeHuckel::m_lnActCoeffMolal, DebyeHuckel::m_speciesCharge_Stoich, and DebyeHuckel::m_tmpV.

◆ initThermo()

void initThermo ( )
virtual

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

Reimplemented from MolalityVPSSTP.

Definition at line 573 of file DebyeHuckel.cpp.

References AnyMap::hasKey(), MolalityVPSSTP::initThermo(), DebyeHuckel::m_form_A_Debye, ThermoPhase::m_input, DebyeHuckel::m_waterSS, DebyeHuckel::setA_Debye(), DebyeHuckel::setBeta(), DebyeHuckel::setDebyeHuckelModel(), DebyeHuckel::setDefaultIonicRadius(), and Phase::species().

◆ initThermoXML()

void initThermoXML ( XML_Node phaseNode,
const std::string &  id 
)
virtual

Import and initialize a ThermoPhase object using an XML tree.

Here we read extra information about the XML description of a phase. Regular information about elements and species and their reference state thermodynamic information have already been read at this point. For example, we do not need to call this function for ideal gas equations of state. This function is called from importPhase() after the elements and the species are initialized with default ideal solution level data.

The default implementation in ThermoPhase calls the virtual function initThermo() and then sets the "state" of the phase by looking for an XML element named "state", and then interpreting its contents by calling the virtual function setStateFromXML().

Parameters
phaseNodeThis object must be the phase node of a complete XML tree description of the phase, including all of the species data. In other words while "phase" must point to an XML phase object, it must have sibling nodes "speciesData" that describe the species in the phase.
idID of the phase. If nonnull, a check is done to see if phaseNode is pointing to the phase with the correct id.
Deprecated:
The XML input format is deprecated and will be removed in Cantera 3.0.

Reimplemented from ThermoPhase.

Definition at line 412 of file DebyeHuckel.cpp.

References XML_Node::attrib(), Cantera::caseInsensitiveEquals(), XML_Node::child(), XML_Node::findByName(), Cantera::fpValue(), Cantera::getFloat(), XML_Node::hasAttrib(), XML_Node::hasChild(), XML_Node::id(), DebyeHuckel::m_formDH, DebyeHuckel::setA_Debye(), DebyeHuckel::setDebyeHuckelModel(), DebyeHuckel::setDefaultIonicRadius(), and Cantera::toSI().

◆ A_Debye_TP()

double A_Debye_TP ( double  temperature = -1.0,
double  pressure = -1.0 
) const
virtual

Return the Debye Huckel constant as a function of temperature and pressure (Units = sqrt(kg/gmol))

The default is to assume that it is constant, given in the initialization process, and stored in the member double, m_A_Debye. Optionally, a full water treatment may be employed that makes \( A_{Debye} \) a full function of T and P.

\[ A_{Debye} = \frac{F e B_{Debye}}{8 \pi \epsilon R T} {\left( C_o \tilde{M}_o \right)}^{1/2} \]

where

\[ B_{Debye} = \frac{F} {{(\frac{\epsilon R T}{2})}^{1/2}} \]

Therefore:

\[ A_{Debye} = \frac{1}{8 \pi} {\left(\frac{2 N_a \rho_o}{1000}\right)}^{1/2} {\left(\frac{N_a e^2}{\epsilon R T }\right)}^{3/2} \]

where

  • Units = sqrt(kg/gmol)
  • \( N_a \) is Avogadro's number
  • \( \rho_w \) is the density of water
  • \( e \) is the electronic charge
  • \( \epsilon = K \epsilon_o \) is the permittivity of water
  • \( K \) is the dielectric constant of water,
  • \( \epsilon_o \) is the permittivity of free space.
  • \( \rho_o \) is the density of the solvent in its standard state.

Nominal value at 298 K and 1 atm = 1.172576 (kg/gmol)^(1/2) based on:

  • \( \epsilon / \epsilon_0 \) = 78.54 (water at 25C)
  • T = 298.15 K
  • B_Debye = 3.28640E9 (kg/gmol)^(1/2)/m
Parameters
temperatureTemperature in kelvin. Defaults to -1, in which case the temperature of the phase is assumed.
pressurePressure (Pa). Defaults to -1, in which case the pressure of the phase is assumed.

Definition at line 631 of file DebyeHuckel.cpp.

Referenced by DebyeHuckel::getMolalityActivityCoefficients(), and DebyeHuckel::s_update_lnMolalityActCoeff().

◆ dA_DebyedT_TP()

double dA_DebyedT_TP ( double  temperature = -1.0,
double  pressure = -1.0 
) const
virtual

Value of the derivative of the Debye Huckel constant with respect to temperature.

This is a function of temperature and pressure. See A_Debye_TP() for a definition of \( A_{Debye} \).

Units = sqrt(kg/gmol) K-1

Parameters
temperatureTemperature in kelvin. Defaults to -1, in which case the temperature of the phase is assumed.
pressurePressure (Pa). Defaults to -1, in which case the pressure of the phase is assumed.

Definition at line 657 of file DebyeHuckel.cpp.

Referenced by DebyeHuckel::getPartialMolarCp(), DebyeHuckel::getPartialMolarEnthalpies(), DebyeHuckel::getPartialMolarEntropies(), DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2(), and DebyeHuckel::s_update_dlnMolalityActCoeff_dT().

◆ d2A_DebyedT2_TP()

double d2A_DebyedT2_TP ( double  temperature = -1.0,
double  pressure = -1.0 
) const
virtual

Value of the 2nd derivative of the Debye Huckel constant with respect to temperature as a function of temperature and pressure.

This is a function of temperature and pressure. See A_Debye_TP() for a definition of \( A_{Debye} \).

Units = sqrt(kg/gmol) K-2

Parameters
temperatureTemperature in kelvin. Defaults to -1, in which case the temperature of the phase is assumed.
pressurePressure (Pa). Defaults to -1, in which case the pressure of the phase is assumed.

Definition at line 681 of file DebyeHuckel.cpp.

Referenced by DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2().

◆ dA_DebyedP_TP()

double dA_DebyedP_TP ( double  temperature = -1.0,
double  pressure = -1.0 
) const
virtual

Value of the derivative of the Debye Huckel constant with respect to pressure, as a function of temperature and pressure.

This is a function of temperature and pressure. See A_Debye_TP() for a definition of \( A_{Debye} \).

Units = sqrt(kg/gmol) Pa-1

Parameters
temperatureTemperature in kelvin. Defaults to -1, in which case the temperature of the phase is assumed.
pressurePressure (Pa). Defaults to -1, in which case the pressure of the phase is assumed.

Definition at line 705 of file DebyeHuckel.cpp.

Referenced by DebyeHuckel::s_update_dlnMolalityActCoeff_dP().

◆ AionicRadius()

double AionicRadius ( int  k = 0) const

Reports the ionic radius of the kth species.

Parameters
kspecies index.

Definition at line 731 of file DebyeHuckel.cpp.

References DebyeHuckel::m_Aionic.

◆ setDebyeHuckelModel()

void setDebyeHuckelModel ( const std::string &  form)

Set the DebyeHuckel parameterization form.

Must be one of 'dilute-limit', 'B-dot-with-variable-a', 'B-dot-with-common-a', 'beta_ij', or 'Pitzer-with-beta_ij'.

Definition at line 337 of file DebyeHuckel.cpp.

References Cantera::caseInsensitiveEquals(), and DebyeHuckel::m_formDH.

Referenced by DebyeHuckel::initThermo(), and DebyeHuckel::initThermoXML().

◆ formDH()

int formDH ( ) const
inline

Returns the form of the Debye-Huckel parameterization used.

Definition at line 928 of file DebyeHuckel.h.

References DebyeHuckel::m_formDH.

◆ setA_Debye()

void setA_Debye ( double  A)

Set the A_Debye parameter.

If a negative value is provided, enables calculation of A_Debye using the detailed water equation of state.

Definition at line 361 of file DebyeHuckel.cpp.

Referenced by DebyeHuckel::initThermo(), and DebyeHuckel::initThermoXML().

◆ setDefaultIonicRadius()

void setDefaultIonicRadius ( double  value)

Set the default ionic radius [m] for each species.

Definition at line 388 of file DebyeHuckel.cpp.

References DebyeHuckel::m_Aionic, and Phase::m_kk.

Referenced by DebyeHuckel::initThermo(), and DebyeHuckel::initThermoXML().

◆ setBeta()

void setBeta ( const std::string &  sp1,
const std::string &  sp2,
double  value 
)

Set the value for the beta interaction between species sp1 and sp2.

Definition at line 397 of file DebyeHuckel.cpp.

References DebyeHuckel::m_Beta_ij, Cantera::npos, and Phase::speciesIndex().

Referenced by DebyeHuckel::initThermo().

◆ get_Beta_ij()

Array2D& get_Beta_ij ( )
inline

Returns a reference to M_Beta_ij.

Definition at line 948 of file DebyeHuckel.h.

References DebyeHuckel::m_Beta_ij.

◆ _nonpolarActCoeff()

double _nonpolarActCoeff ( double  IionicMolality)
staticprivate

Static function that implements the non-polar species salt-out modifications.

Returns the calculated activity coefficients.

Parameters
IionicMolalityValue of the ionic molality (sqrt(gmol/kg))

Definition at line 785 of file DebyeHuckel.cpp.

◆ _osmoticCoeffHelgesonFixedForm()

double _osmoticCoeffHelgesonFixedForm ( ) const
private

Formula for the osmotic coefficient that occurs in the GWB.

It is originally from Helgeson for a variable NaCl brine. It's to be used with extreme caution.

Definition at line 799 of file DebyeHuckel.cpp.

References DebyeHuckel::m_A_Debye, and DebyeHuckel::m_IionicMolalityStoich.

Referenced by DebyeHuckel::_lnactivityWaterHelgesonFixedForm().

◆ _lnactivityWaterHelgesonFixedForm()

double _lnactivityWaterHelgesonFixedForm ( ) const
private

Formula for the log of the water activity that occurs in the GWB.

It is originally from Helgeson for a variable NaCl brine. It's to be used with extreme caution.

Definition at line 818 of file DebyeHuckel.cpp.

References DebyeHuckel::_osmoticCoeffHelgesonFixedForm(), MolalityVPSSTP::calcMolalities(), Phase::m_kk, DebyeHuckel::m_maxIionicStrength, MolalityVPSSTP::m_Mnaught, and MolalityVPSSTP::m_molalities.

◆ s_update_lnMolalityActCoeff()

void s_update_lnMolalityActCoeff ( ) const
private

◆ s_update_dlnMolalityActCoeff_dT()

void s_update_dlnMolalityActCoeff_dT ( ) const
private

Calculation of temperature derivative of activity coefficient.

Using internally stored values, this function calculates the temperature derivative of the logarithm of the activity coefficient for all species in the mechanism.

We assume that the activity coefficients are current in this routine. The solvent activity coefficient is on the molality scale. Its derivative is too.

Definition at line 1056 of file DebyeHuckel.cpp.

References DebyeHuckel::dA_DebyedT_TP(), DebyeHuckel::m_B_Debye, DebyeHuckel::m_dlnActCoeffMolaldT, DebyeHuckel::m_formDH, DebyeHuckel::m_IionicMolality, Phase::m_kk, and Phase::moleFraction().

Referenced by DebyeHuckel::getPartialMolarCp(), DebyeHuckel::getPartialMolarEnthalpies(), and DebyeHuckel::getPartialMolarEntropies().

◆ s_update_d2lnMolalityActCoeff_dT2()

void s_update_d2lnMolalityActCoeff_dT2 ( ) const
private

Calculate the temperature 2nd derivative of the activity coefficient.

Using internally stored values, this function calculates the temperature 2nd derivative of the logarithm of the activity coefficient for all species in the mechanism.

We assume that the activity coefficients are current in this routine. Solvent activity coefficient is on the molality scale. Its derivatives are too.

Definition at line 1166 of file DebyeHuckel.cpp.

References DebyeHuckel::d2A_DebyedT2_TP(), DebyeHuckel::dA_DebyedT_TP(), DebyeHuckel::m_B_Debye, DebyeHuckel::m_d2lnActCoeffMolaldT2, DebyeHuckel::m_formDH, DebyeHuckel::m_IionicMolality, Phase::m_kk, and Phase::moleFraction().

Referenced by DebyeHuckel::getPartialMolarCp().

◆ s_update_dlnMolalityActCoeff_dP()

void s_update_dlnMolalityActCoeff_dP ( ) const
private

Calculate the pressure derivative of the activity coefficient.

Using internally stored values, this function calculates the pressure derivative of the logarithm of the activity coefficient for all species in the mechanism.

We assume that the activity coefficients, molalities, and A_Debye are current. Solvent activity coefficient is on the molality scale. Its derivatives are too.

Definition at line 1272 of file DebyeHuckel.cpp.

References DebyeHuckel::dA_DebyedP_TP(), DebyeHuckel::m_B_Debye, DebyeHuckel::m_dlnActCoeffMolaldP, DebyeHuckel::m_formDH, DebyeHuckel::m_IionicMolality, Phase::m_kk, and Phase::moleFraction().

Referenced by DebyeHuckel::getPartialMolarVolumes().

Member Data Documentation

◆ m_formDH

int m_formDH
protected

form of the Debye-Huckel parameterization used in the model.

The options are described at the top of this document, and in the general documentation. The list is repeated here:

DHFORM_DILUTE_LIMIT = 0 (default) DHFORM_BDOT_AK = 1 DHFORM_BDOT_AUNIFORM = 2 DHFORM_BETAIJ = 3 DHFORM_PITZER_BETAIJ = 4

Definition at line 990 of file DebyeHuckel.h.

Referenced by DebyeHuckel::formDH(), DebyeHuckel::initThermoXML(), DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2(), DebyeHuckel::s_update_dlnMolalityActCoeff_dP(), DebyeHuckel::s_update_dlnMolalityActCoeff_dT(), DebyeHuckel::s_update_lnMolalityActCoeff(), and DebyeHuckel::setDebyeHuckelModel().

◆ m_electrolyteSpeciesType

vector_int m_electrolyteSpeciesType
protected

Vector containing the electrolyte species type.

The possible types are:

  • solvent
  • Charged Species
  • weakAcidAssociated
  • strongAcidAssociated
  • polarNeutral
  • nonpolarNeutral

Definition at line 1003 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies().

◆ m_Aionic

vector_fp m_Aionic
protected

a_k = Size of the ionic species in the DH formulation. units = meters

Definition at line 1006 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies(), DebyeHuckel::AionicRadius(), and DebyeHuckel::setDefaultIonicRadius().

◆ m_IionicMolality

double m_IionicMolality
mutableprotected

◆ m_maxIionicStrength

double m_maxIionicStrength
protected

Maximum value of the ionic strength allowed in the calculation of the activity coefficients.

Definition at line 1013 of file DebyeHuckel.h.

Referenced by DebyeHuckel::_lnactivityWaterHelgesonFixedForm(), and DebyeHuckel::s_update_lnMolalityActCoeff().

◆ m_useHelgesonFixedForm

bool m_useHelgesonFixedForm

If true, then the fixed for of Helgeson's activity for water is used instead of the rigorous form obtained from Gibbs-Duhem relation.

This should be used with caution, and is really only included as a validation exercise.

Definition at line 1020 of file DebyeHuckel.h.

◆ m_IionicMolalityStoich

double m_IionicMolalityStoich
mutableprotected

Stoichiometric ionic strength on the molality scale.

Definition at line 1023 of file DebyeHuckel.h.

Referenced by DebyeHuckel::_osmoticCoeffHelgesonFixedForm(), and DebyeHuckel::s_update_lnMolalityActCoeff().

◆ m_form_A_Debye

int m_form_A_Debye

Form of the constant outside the Debye-Huckel term called A.

It's normally a function of temperature and pressure. However, it can be set from the input file in order to aid in numerical comparisons. Acceptable forms:

  A_DEBYE_CONST  0
  A_DEBYE_WATER  1

The A_DEBYE_WATER form may be used for water solvents with needs to cover varying temperatures and pressures. Note, the dielectric constant of water is a relatively strong function of T, and its variability must be accounted for,

Definition at line 1042 of file DebyeHuckel.h.

Referenced by DebyeHuckel::initThermo().

◆ m_A_Debye

double m_A_Debye
mutableprotected

Current value of the Debye Constant, A_Debye.

A_Debye -> this expression appears on the top of the ln actCoeff term in the general Debye-Huckel expression It depends on temperature and pressure.

A_Debye = (F e B_Debye) / (8 Pi epsilon R T)

Units = sqrt(kg/gmol)

Nominal value(298K, atm) = 1.172576 sqrt(kg/gmol) based on: epsilon/epsilon_0 = 78.54 (water at 25C) T = 298.15 K B_Debye = 3.28640E9 sqrt(kg/gmol)/m

note in Pitzer's nomenclature, A_phi = A_Debye/3.0

Definition at line 1064 of file DebyeHuckel.h.

Referenced by DebyeHuckel::_osmoticCoeffHelgesonFixedForm(), and DebyeHuckel::s_update_lnMolalityActCoeff().

◆ m_B_Debye

double m_B_Debye
protected

Current value of the constant that appears in the denominator.

B_Debye -> this expression appears on the bottom of the ln actCoeff term in the general Debye-Huckel expression It depends on temperature

B_Bebye = F / sqrt( epsilon R T / 2 )

Units = sqrt(kg/gmol) / m

Nominal value = 3.28640E9 sqrt(kg/gmol) / m based on: epsilon/epsilon_0 = 78.54 (water at 25C) T = 298.15 K

Definition at line 1082 of file DebyeHuckel.h.

Referenced by DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2(), DebyeHuckel::s_update_dlnMolalityActCoeff_dP(), DebyeHuckel::s_update_dlnMolalityActCoeff_dT(), and DebyeHuckel::s_update_lnMolalityActCoeff().

◆ m_B_Dot

vector_fp m_B_Dot
protected

Array of B_Dot values.

This expression is an extension of the Debye-Huckel expression used in some formulations to extend DH to higher molalities. B_dot is specific to the major ionic pair.

Definition at line 1090 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies().

◆ m_waterSS

PDSS_Water* m_waterSS
protected

Pointer to the Water standard state object.

derived from the equation of state for water.

Definition at line 1096 of file DebyeHuckel.h.

Referenced by DebyeHuckel::calcDensity(), and DebyeHuckel::initThermo().

◆ m_densWaterSS

double m_densWaterSS
protected

Storage for the density of water's standard state.

Density depends on temperature and pressure.

Definition at line 1102 of file DebyeHuckel.h.

Referenced by DebyeHuckel::calcDensity().

◆ m_waterProps

std::unique_ptr<WaterProps> m_waterProps
protected

Pointer to the water property calculator.

Definition at line 1105 of file DebyeHuckel.h.

◆ m_tmpV

vector_fp m_tmpV
mutableprotected

vector of size m_kk, used as a temporary holding area.

Definition at line 1108 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies(), DebyeHuckel::calcDensity(), DebyeHuckel::cp_mole(), DebyeHuckel::enthalpy_mole(), DebyeHuckel::entropy_mole(), and DebyeHuckel::gibbs_mole().

◆ m_speciesCharge_Stoich

vector_fp m_speciesCharge_Stoich
protected

Stoichiometric species charge -> This is for calculations of the ionic strength which ignore ion-ion pairing into neutral molecules.

The Stoichiometric species charge is the charge of one of the ion that would occur if the species broke into two charged ion pairs. NaCl -> m_speciesCharge_Stoich = -1; HSO4- -> H+ + SO42- = -2 -> The other charge is calculated. For species that aren't ion pairs, it's equal to the m_speciesCharge[] value.

Definition at line 1122 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies(), and DebyeHuckel::s_update_lnMolalityActCoeff().

◆ m_Beta_ij

Array2D m_Beta_ij
protected

Array of 2D data used in the DHFORM_BETAIJ formulation Beta_ij.value(i,j) is the coefficient of the jth species for the specification of the chemical potential of the ith species.

Definition at line 1130 of file DebyeHuckel.h.

Referenced by DebyeHuckel::get_Beta_ij(), and DebyeHuckel::setBeta().

◆ m_lnActCoeffMolal

vector_fp m_lnActCoeffMolal
mutableprotected

Logarithm of the activity coefficients on the molality scale.

mutable because we change this if the composition or temperature or pressure changes.

Definition at line 1137 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies(), DebyeHuckel::getActivities(), DebyeHuckel::getChemPotentials(), DebyeHuckel::getMolalityActivityCoefficients(), and DebyeHuckel::getPartialMolarEntropies().

◆ m_dlnActCoeffMolaldT

vector_fp m_dlnActCoeffMolaldT
mutableprotected

◆ m_d2lnActCoeffMolaldT2

vector_fp m_d2lnActCoeffMolaldT2
mutableprotected

2nd Derivative of log act coeff wrt T

Definition at line 1143 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies(), DebyeHuckel::getPartialMolarCp(), and DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2().

◆ m_dlnActCoeffMolaldP

vector_fp m_dlnActCoeffMolaldP
mutableprotected

Derivative of log act coeff wrt P.

Definition at line 1146 of file DebyeHuckel.h.

Referenced by DebyeHuckel::addSpecies(), DebyeHuckel::getPartialMolarVolumes(), and DebyeHuckel::s_update_dlnMolalityActCoeff_dP().


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