32 for (
size_t k = 0; k <
m_kk; k++) {
47 for (
size_t k = 0; k <
m_kk; k++) {
56 vector<double> hbar(
m_kk);
58 for (
size_t i = 0; i <
m_kk; i++) {
67 vector<double> sbar(
m_kk);
69 for (
size_t i = 0; i <
m_kk; i++) {
78 vector<double> cpbar(
m_kk);
80 for (
size_t i = 0; i <
m_kk; i++) {
97 for (
size_t k = 0; k <
m_kk; k++) {
105 for (
size_t k = 0; k <
m_kk; k++) {
114 for (
size_t k = 0; k <
m_kk; k++) {
130 for (
size_t k = 0; k <
m_kk; k++) {
135 for (
size_t k = 0; k <
m_kk; k++) {
144 for (
size_t iK = 0; iK <
m_kk; iK++) {
152 for (
const auto& item :
m_input[
"interactions"].asVector<
AnyMap>()) {
153 auto&
species = item[
"species"].asVector<
string>(2);
154 vector<double> h_excess = item.convertVector(
"excess-enthalpy",
"J/kmol");
155 vector<double> s_excess = item.convertVector(
"excess-entropy",
"J/kmol/K");
157 h_excess.data(), h_excess.size(),
158 s_excess.data(), s_excess.size());
168 vector<AnyMap> interactions;
171 interaction[
"species"] = vector<string>{
175 while (h.size() > 1 && h.back() == 0) {
178 while (s.size() > 1 && s.back() == 0) {
181 interaction[
"excess-enthalpy"].setQuantity(std::move(h),
"J/kmol");
182 interaction[
"excess-entropy"].setQuantity(std::move(s),
"J/kmol/K");
183 interactions.push_back(std::move(interaction));
185 phaseNode[
"interactions"] = std::move(interactions);
202 for (
size_t i = 0; i <
m_HE_m_ij.size(); i++) {
207 double deltaX = XA - XB;
208 const vector<double>& he_vec =
m_HE_m_ij[i];
209 const vector<double>& se_vec =
m_SE_m_ij[i];
211 double polyMm1 = 1.0;
215 for (
size_t m = 0; m < he_vec.size(); m++) {
216 double A_ge = (he_vec[m] - T * se_vec[m]) / (
GasConstant * T);
218 sum2 += A_ge * (m + 1) * poly;
221 sumMm1 += (A_ge * polyMm1 * m);
225 double oneMXA = 1.0 - XA;
226 double oneMXB = 1.0 - XB;
227 for (
size_t k = 0; k <
m_kk; k++) {
230 }
else if (iB == k) {
245 for (
size_t i = 0; i <
m_HE_m_ij.size(); i++) {
250 double deltaX = XA - XB;
253 const vector<double>& he_vec =
m_HE_m_ij[i];
255 double polyMm1 = 1.0;
257 for (
size_t m = 0; m < he_vec.size(); m++) {
260 sum2 += h_e * (m + 1) * poly;
263 sumMm1 += (h_e * polyMm1 * m);
267 double oneMXA = 1.0 - XA;
268 double oneMXB = 1.0 - XB;
269 for (
size_t k = 0; k <
m_kk; k++) {
272 }
else if (iB == k) {
280 for (
size_t k = 0; k <
m_kk; k++) {
288 for (
size_t k = 0; k <
m_kk; k++) {
296 for (
size_t k = 0; k <
m_kk; k++) {
306 for (
size_t i = 0; i <
m_HE_m_ij.size(); i++) {
311 double deltaX = XA - XB;
314 const vector<double>& he_vec =
m_HE_m_ij[i];
315 const vector<double>& se_vec =
m_SE_m_ij[i];
317 double polyMm1 = 1.0;
318 double polyMm2 = 1.0;
320 for (
size_t m = 0; m < he_vec.size(); m++) {
321 double A_ge = (he_vec[m] - T * se_vec[m]) / (
GasConstant * T);;
325 sumMm1 += (A_ge * polyMm1 * m);
329 sumMm2 += (A_ge * polyMm2 * m * (m - 1.0));
334 for (
size_t k = 0; k <
m_kk; k++) {
337 XA * (- (1-XA+XB) * sum + 2*(1.0 - XA) * XB * sumMm1
338 + sumMm1 * (XB * (1 - 2*XA + XB) - XA * (1 - XA + 2*XB))
339 + 2 * XA * XB * sumMm2 * (1.0 - XA + XB));
340 }
else if (iB == k) {
342 XB * (- (1-XB+XA) * sum - 2*(1.0 - XB) * XA * sumMm1
343 + sumMm1 * (XA * (2*XB - XA - 1) - XB * (-2*XA + XB - 1))
344 - 2 * XA * XB * sumMm2 * (-XA - 1 + XB));
355 for (
size_t i = 0; i <
m_HE_m_ij.size(); i++) {
360 double deltaX = XA - XB;
363 const vector<double>& he_vec =
m_HE_m_ij[i];
364 const vector<double>& se_vec =
m_SE_m_ij[i];
366 double polyMm1 = 1.0;
367 double polyMm2 = 1.0;
369 double sum2Mm1 = 0.0;
371 for (
size_t m = 0; m < he_vec.size(); m++) {
372 double A_ge = he_vec[m] - T * se_vec[m];
374 sum2 += A_ge * (m + 1) * poly;
377 sumMm1 += (A_ge * polyMm1 * m);
378 sum2Mm1 += (A_ge * polyMm1 * m * (1.0 + m));
382 sumMm2 += (A_ge * polyMm2 * m * (m - 1.0));
387 for (
size_t k = 0; k <
m_kk; k++) {
390 + XB * sumMm1 * (1.0 - 2.0 * XA + XB)
391 + XA * XB * sumMm2 * (1.0 - XA + XB));
394 + XA * sumMm1 * (1.0 + 2.0 * XB - XA)
395 - XA * XB * sumMm2 * (1.0 - XA + XB));
396 }
else if (iB == k) {
398 + XB * sumMm1 * (1.0 - 2.0 * XA + XB)
399 + XA * XB * sumMm2 * (1.0 - XA + XB));
402 + XA * sumMm1 * (XB - XA - (1.0 - XB))
403 - XA * XB * sumMm2 * (-XA - (1.0 - XB)));
413 double* dlnActCoeffds)
const
417 for (
size_t k = 0; k <
m_kk; k++) {
419 for (
size_t j = 0; j <
m_kk; j++) {
428 for (
size_t j = 0; j <
m_kk; j++) {
430 for (
size_t k = 0; k <
m_kk; k++) {
439 for (
size_t k = 0; k <
m_kk; k++) {
448 for (
size_t k = 0; k <
m_kk; k++) {
449 for (
size_t m = 0; m <
m_kk; m++) {
450 dlnActCoeffdlnN[ld * k + m] = data[
m_kk * k + m];
456 const string& speciesA,
const string& speciesB,
457 const double* excess_enthalpy,
size_t n_enthalpy,
458 const double* excess_entropy,
size_t n_entropy)
463 throw CanteraError(
"RedlichKisterVPSSTP::addBinaryInteraction",
464 "Species '{}' not present in phase", speciesA);
465 }
else if (kB ==
npos) {
466 throw CanteraError(
"RedlichKisterVPSSTP::addBinaryInteraction",
467 "Species '{}' not present in phase", speciesB);
470 throw CanteraError(
"RedlichKisterVPSSTP::addBinaryInteraction",
471 "Species '{}' should be neutral", speciesA);
472 }
else if (
charge(kB) != 0) {
473 throw CanteraError(
"RedlichKisterVPSSTP::addBinaryInteraction",
474 "Species '{}' should be neutral", speciesB);
479 m_HE_m_ij.emplace_back(excess_enthalpy, excess_enthalpy + n_enthalpy);
480 m_SE_m_ij.emplace_back(excess_entropy, excess_entropy + n_entropy);
481 size_t N = max(n_enthalpy, n_entropy);
(see Thermodynamic Properties and class RedlichKisterVPSSTP).
Headers for the factory class that can create known ThermoPhase objects (see Thermodynamic Properties...
A map of string keys to values whose type can vary at runtime.
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
void zero()
Set all of the entries to zero.
virtual void resize(size_t n, size_t m, double v=0.0)
Resize the array, and fill the new entries with 'v'.
Base class for exceptions thrown by Cantera classes.
vector< double > d2lnActCoeffdT2_Scaled_
Storage for the current derivative values of the gradients with respect to temperature of the log of ...
Array2D dlnActCoeffdlnN_
Storage for the current derivative values of the gradients with respect to logarithm of the species m...
vector< double > lnActCoeff_Scaled_
Storage for the current values of the activity coefficients of the species.
vector< double > dlnActCoeffdlnX_diag_
Storage for the current derivative values of the gradients with respect to logarithm of the mole frac...
vector< double > moleFractions_
Storage for the current values of the mole fractions of the species.
vector< double > dlnActCoeffdT_Scaled_
Storage for the current derivative values of the gradients with respect to temperature of the log of ...
size_t m_kk
Number of species in the phase.
double temperature() const
Temperature (K).
string speciesName(size_t k) const
Name of the species with index k.
size_t speciesIndex(const string &name) const
Returns the index of a species named 'name' within the Phase object.
shared_ptr< Species > species(const string &name) const
Return the Species object for the named species.
double charge(size_t k) const
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the...
void getdlnActCoeffds(const double dTds, const double *const dXds, double *dlnActCoeffds) const override
Get the change in activity coefficients wrt changes in state (temp, mole fraction,...
void s_update_dlnActCoeff_dX_() const
Internal routine that calculates the derivative of the activity coefficients wrt the mole fractions.
double enthalpy_mole() const override
Molar enthalpy. Units: J/kmol.
void getPartialMolarEnthalpies(double *hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
void getChemPotentials(double *mu) const override
Get the species chemical potentials. Units: J/kmol.
Array2D dlnActCoeff_dX_
Two dimensional array of derivatives of activity coefficients wrt mole fractions.
vector< vector< double > > m_SE_m_ij
Entropy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
void getParameters(AnyMap &phaseNode) const override
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
void getPartialMolarVolumes(double *vbar) const override
Return an array of partial molar volumes for the species in the mixture.
vector< size_t > m_pSpecies_A_ij
vector of species indices representing species A in the interaction
double cv_mole() const override
Molar heat capacity at constant volume. Units: J/kmol/K.
void s_update_dlnActCoeff_dT() const
Update the derivative of the log of the activity coefficients wrt T.
vector< size_t > m_pSpecies_B_ij
vector of species indices representing species B in the interaction
void addBinaryInteraction(const string &speciesA, const string &speciesB, const double *excess_enthalpy, size_t n_enthalpy, const double *excess_entropy, size_t n_entropy)
Add a binary species interaction with the specified parameters.
vector< vector< double > > m_HE_m_ij
Enthalpy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
double entropy_mole() const override
Molar entropy. Units: J/kmol/K.
void getdlnActCoeffdT(double *dlnActCoeffdT) const override
Get the array of temperature derivatives of the log activity coefficients.
double cp_mole() const override
Molar heat capacity at constant pressure. Units: J/kmol/K.
void getPartialMolarCp(double *cpbar) const override
Returns an array of partial molar heat capacities for the species in the mixture.
void initLengths()
Initialize lengths of local variables after all species have been identified.
void getLnActivityCoefficients(double *lnac) const override
Get the array of non-dimensional molar-based ln activity coefficients at the current solution tempera...
void s_update_dlnActCoeff_dlnX_diag() const
Internal routine that calculates the total derivative of the activity coefficients with respect to th...
void s_update_lnActCoeff() const
Update the activity coefficients.
void getdlnActCoeffdlnX_diag(double *dlnActCoeffdlnX_diag) const override
Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component o...
RedlichKisterVPSSTP(const string &inputFile="", const string &id="")
Construct a RedlichKisterVPSSTP object from an input file.
void getdlnActCoeffdlnN_diag(double *dlnActCoeffdlnN_diag) const override
Get the array of log species mole number derivatives of the log activity coefficients.
void getd2lnActCoeffdT2(double *d2lnActCoeffdT2) const
Get the array of temperature second derivatives of the log activity coefficients.
void getPartialMolarEntropies(double *sbar) const override
Returns an array of partial molar entropies for the species in the mixture.
void getdlnActCoeffdlnN(const size_t ld, double *const dlnActCoeffdlnN) override
Get the array of derivatives of the log activity coefficients with respect to the log of the species ...
virtual void getParameters(AnyMap &phaseNode) const
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
double RT() const
Return the Gas Constant multiplied by the current temperature.
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
void initThermoFile(const string &inputFile, const string &id)
Initialize a ThermoPhase object using an input file.
AnyMap m_input
Data supplied via setParameters.
void getEntropy_R(double *sr) const override
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
void getStandardChemPotentials(double *mu) const override
Get the array of chemical potentials at unit activity for the species at their standard states at the...
void getCp_R(double *cpr) const override
Get the nondimensional Heat Capacities at constant pressure for the species standard states at the cu...
void getEnthalpy_RT(double *hrt) const override
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
void getStandardVolumes(double *vol) const override
Get the molar volumes of the species standard states at the current T and P of the solution.
const double GasConstant
Universal Gas Constant [J/kmol/K].
Namespace for the Cantera kernel.
const size_t npos
index returned by functions to indicate "no position"
const double SmallNumber
smallest number to compare to zero.
Contains declarations for string manipulation functions within Cantera.