Cantera
2.1.2
|
This is the main structure used to hold the internal data used in vcs_solve_TP(), and to solve TP systems. More...
#include <vcs_solve.h>
Public Member Functions | |
void | vcs_initSizes (const size_t nspecies0, const size_t nelements, const size_t nphase0) |
Initialize the sizes within the VCS_SOLVE object. More... | |
int | vcs (VCS_PROB *vprob, int ifunc, int ipr, int ip1, int maxit) |
Solve an equilibrium problem. More... | |
int | vcs_solve_TP (int print_lvl, int printDetails, int maxit) |
Main routine that solves for equilibrium at constant T and P using a variant of the VCS method. More... | |
int | vcs_PS (VCS_PROB *vprob, int iph, int printLvl, double &feStable) |
void | vcs_reinsert_deleted (size_t kspec) |
int | vcs_basopt (const bool doJustComponents, double aw[], double sa[], double sm[], double ss[], double test, bool *const usedZeroedSpecies) |
Choose the optimum species basis for the calculations. More... | |
size_t | vcs_basisOptMax (const double *const molNum, const size_t j, const size_t n) |
Choose a species to test for the next component. More... | |
int | vcs_species_type (const size_t kspec) const |
Evaluate the species category for the indicated species. More... | |
bool | vcs_evaluate_speciesType () |
This routine evaluates the species type for all species. More... | |
void | vcs_chemPotPhase (const int stateCalc, const size_t iph, const double *const molNum, double *const ac, double *const mu_i, const bool do_deleted=false) |
We calculate the dimensionless chemical potentials of all species in a single phase. More... | |
void | vcs_dfe (const int stateCalc, const int ll, const size_t lbot, const size_t ltop) |
Calculate the dimensionless chemical potentials of all species or of certain groups of species, at a fixed temperature and pressure. More... | |
void | vcs_printSpeciesChemPot (const int stateCalc) const |
Print out a table of chemical potentials. More... | |
void | vcs_updateVP (const int stateCalc) |
This routine uploads the state of the system into all of the vcs_VolumePhase objects in the current problem. More... | |
bool | vcs_popPhasePossible (const size_t iphasePop) const |
Utility function that evaluates whether a phase can be popped into existence. More... | |
int | vcs_phasePopDeterminePossibleList () |
Determine the list of problems that need to be checked to see if there are any phases pops. More... | |
size_t | vcs_popPhaseID (std::vector< size_t > &phasePopPhaseIDs) |
Decision as to whether a phase pops back into existence. More... | |
int | vcs_popPhaseRxnStepSizes (const size_t iphasePop) |
Calculates the deltas of the reactions due to phases popping into existence. More... | |
size_t | vcs_RxnStepSizes (int &forceComponentCalc, size_t &kSpecial) |
Calculates formation reaction step sizes. More... | |
double | vcs_tmoles () |
Calculates the total number of moles of species in all phases. More... | |
void | check_tmoles () const |
void | vcs_deltag (const int l, const bool doDeleted, const int vcsState, const bool alterZeroedPhases=true) |
This subroutine calculates reaction free energy changes for all noncomponent formation reactions. More... | |
void | vcs_printDeltaG (const int stateCalc) |
void | vcs_deltag_Phase (const size_t iphase, const bool doDeleted, const int stateCalc, const bool alterZeroedPhases=true) |
Calculate deltag of formation for all species in a single phase. More... | |
void | vcs_switch_pos (const bool ifunc, const size_t k1, const size_t k2) |
Swaps the indices for all of the global data for two species, k1 and k2. More... | |
double | vcs_birthGuess (const int kspec) |
Birth guess returns the number of moles of a species that is coming back to life. More... | |
int | vcs_solve_phaseStability (const int iphase, int ifunc, double &funcval, int print_lvl) |
Routine that independently determines whether a phase should be popped under the current conditions. More... | |
double | vcs_phaseStabilityTest (const size_t iph) |
Main program to test whether a deleted phase should be brought back into existence. More... | |
int | vcs_TP (int ipr, int ip1, int maxit, double T, double pres) |
Solve an equilibrium problem at a particular fixed temperature and pressure. More... | |
int | vcs_evalSS_TP (int ipr, int ip1, double Temp, double pres) |
void | vcs_fePrep_TP () |
Initialize the chemical potential of single species phases. More... | |
double | vcs_VolTotal (const double tkelvin, const double pres, const double w[], double volPM[]) |
Calculation of the total volume and the partial molar volumes. More... | |
int | vcs_prep_oneTime (int printLvl) |
This routine is mostly concerned with changing the private data to be consistent with what's needed for solution. More... | |
int | vcs_prep () |
Prepare the object for solution. More... | |
bool | vcs_wellPosed (VCS_PROB *vprob) |
In this routine, we check for things that will cause the algorithm to fail. More... | |
int | vcs_elem_rearrange (double *const aw, double *const sa, double *const sm, double *const ss) |
Rearrange the constraint equations represented by the Formula Matrix so that the operational ones are in the front. More... | |
void | vcs_switch_elem_pos (size_t ipos, size_t jpos) |
Swaps the indices for all of the global data for two elements, ipos and jpos. More... | |
int | vcs_rxn_adj_cg (void) |
Calculates reaction adjustments using a full Hessian approximation. More... | |
double | vcs_Hessian_diag_adj (size_t irxn, double hessianDiag_Ideal) |
Calculates the diagonal contribution to the Hessian due to the dependence of the activity coefficients on the mole numbers. More... | |
double | vcs_Hessian_actCoeff_diag (size_t irxn) |
Calculates the diagonal contribution to the Hessian due to the dependence of the activity coefficients on the mole numbers. More... | |
void | vcs_CalcLnActCoeffJac (const double *const moleSpeciesVCS) |
Recalculate all of the activity coefficients in all of the phases based on input mole numbers. More... | |
double | vcs_line_search (const size_t irxn, const double dx_orig, char *const ANOTE) |
A line search algorithm is carried out on one reaction. More... | |
int | vcs_report (int iconv) |
Print out a report on the state of the equilibrium problem to standard output. More... | |
int | vcs_rearrange () |
Switch all species data back to the original order. More... | |
double | vcs_nondim_Farad (int mu_units, double TKelvin) const |
Returns the multiplier for electric charge terms. More... | |
double | vcs_nondimMult_TP (int mu_units, double TKelvin) const |
Returns the multiplier for the nondimensionalization of the equations. More... | |
void | vcs_nondim_TP () |
Nondimensionalize the problem data. More... | |
void | vcs_redim_TP () |
Redimensionalize the problem data. More... | |
void | vcs_printChemPotUnits (int unitsFormat) const |
Print the string representing the Chemical potential units. More... | |
void | vcs_elab () |
Computes the current elemental abundances vector. More... | |
bool | vcs_elabcheck (int ibound) |
void | vcs_elabPhase (size_t iphase, double *const elemAbundPhase) |
int | vcs_elcorr (double aa[], double x[]) |
int | vcs_inest_TP () |
Create an initial estimate of the solution to the thermodynamic equilibrium problem. More... | |
int | vcs_setMolesLinProg () |
Estimate the initial mole numbers by constrained linear programming. More... | |
double | vcs_Total_Gibbs (double *w, double *fe, double *tPhMoles) |
Calculate the total dimensionless Gibbs free energy. More... | |
double | vcs_GibbsPhase (size_t iphase, const double *const w, const double *const fe) |
Calculate the total dimensionless Gibbs free energy of a single phase. More... | |
int | vcs_prob_update (VCS_PROB *pub) |
Transfer the results of the equilibrium calculation back to VCS_PROB. More... | |
int | vcs_prob_specifyFully (const VCS_PROB *pub) |
Fully specify the problem to be solved using VCS_PROB. More... | |
int | vcs_prob_specify (const VCS_PROB *pub) |
Specify the problem to be solved using VCS_PROB, incrementally. More... | |
int | vcs_rank (const double *awtmp, size_t numSpecies, const double *matrix, size_t numElemConstraints, std::vector< size_t > &compRes, std::vector< size_t > &elemComp, int *const usedZeroedSpecies) const |
Calculate the rank of a matrix and return the rows and columns that will generate an independent basis for that rank. More... | |
Public Attributes | |
size_t | NSPECIES0 |
value of the number of species used to malloc data structures More... | |
size_t | NPHASE0 |
value of the number of phases used to malloc data structures More... | |
size_t | m_numSpeciesTot |
Total number of species in the problems. More... | |
size_t | m_numElemConstraints |
Number of element constraints in the problem. More... | |
size_t | m_numComponents |
Number of components calculated for the problem. More... | |
size_t | m_numRxnTot |
Total number of non-component species in the problem. More... | |
size_t | m_numSpeciesRdc |
Current number of species in the problems. More... | |
size_t | m_numRxnRdc |
Current number of non-component species in the problem. More... | |
size_t | m_numRxnMinorZeroed |
Number of active species which are currently either treated as minor species. More... | |
size_t | m_numPhases |
Number of Phases in the problem. More... | |
DoubleStarStar | m_formulaMatrix |
Formula matrix for the problem. More... | |
DoubleStarStar | m_stoichCoeffRxnMatrix |
Stoichiometric coefficient matrix for the reaction mechanism expressed in Reduced Canonical Form. More... | |
std::vector< double > | m_scSize |
Absolute size of the stoichiometric coefficients. More... | |
std::vector< double > | m_spSize |
total size of the species More... | |
std::vector< double > | m_SSfeSpecies |
Standard state chemical potentials for species K at the current temperature and pressure. More... | |
std::vector< double > | m_feSpecies_old |
Free energy vector from the start of the current iteration. More... | |
std::vector< double > | m_feSpecies_new |
Dimensionless new free energy for all the species in the mechanism at the new tentatite T, P, and mole numbers. More... | |
int | m_doEstimateEquil |
Setting for whether to do an initial estimate. More... | |
std::vector< double > | m_molNumSpecies_old |
Total moles of the species. More... | |
std::vector< int > | m_speciesUnknownType |
Specifies the species unknown type. More... | |
DoubleStarStar | m_deltaMolNumPhase |
Change in the number of moles of phase, iphase, due to the noncomponent formation reaction, irxn, for species, k: More... | |
IntStarStar | m_phaseParticipation |
This is 1 if the phase, iphase, participates in the formation reaction irxn, and zero otherwise. More... | |
std::vector< double > | m_phasePhi |
electric potential of the iph phase More... | |
std::vector< double > | m_molNumSpecies_new |
Tentative value of the mole number vector. More... | |
std::vector< double > | m_deltaGRxn_new |
Delta G(irxn) for the noncomponent species in the mechanism. More... | |
std::vector< double > | m_deltaGRxn_old |
Last deltag[irxn] from the previous step. More... | |
std::vector< double > | m_deltaGRxn_Deficient |
Last deltag[irxn] from the previous step with additions for possible births of zeroed phases. More... | |
std::vector< double > | m_deltaGRxn_tmp |
Temporary vector of Rxn DeltaG's. More... | |
std::vector< double > | m_deltaMolNumSpecies |
Reaction Adjustments for each species during the current step. More... | |
std::vector< double > | m_elemAbundances |
Element abundances vector. More... | |
std::vector< double > | m_elemAbundancesGoal |
Element abundances vector Goals. More... | |
double | m_totalMolNum |
Total number of kmoles in all phases. More... | |
std::vector< double > | m_tPhaseMoles_old |
Total kmols of species in each phase. More... | |
std::vector< double > | m_tPhaseMoles_new |
total kmols of species in each phase in the tentative soln vector More... | |
std::vector< double > | m_TmpPhase |
Temporary vector of length NPhase. More... | |
std::vector< double > | m_TmpPhase2 |
Temporary vector of length NPhase. More... | |
std::vector< double > | m_deltaPhaseMoles |
Change in the total moles in each phase. More... | |
double | m_temperature |
Temperature (Kelvin) More... | |
double | m_pressurePA |
Pressure (units are determined by m_VCS_UnitsFormat. More... | |
std::vector< double > | TPhInertMoles |
Total kmoles of inert to add to each phase. More... | |
double | m_tolmaj |
Tolerance requirement for major species. More... | |
double | m_tolmin |
Tolerance requirements for minor species. More... | |
double | m_tolmaj2 |
Below this, major species aren't refined any more. More... | |
double | m_tolmin2 |
Below this, minor species aren't refined any more. More... | |
std::vector< size_t > | m_speciesMapIndex |
Index vector that keeps track of the species vector rearrangement. More... | |
std::vector< size_t > | m_speciesLocalPhaseIndex |
Index that keeps track of the index of the species within the local phase. More... | |
std::vector< size_t > | m_elementMapIndex |
Index vector that keeps track of the rearrangement of the elements. More... | |
std::vector< size_t > | m_indexRxnToSpecies |
Mapping between the species index for noncomponent species and the full species index. More... | |
std::vector< int > | m_speciesStatus |
Major -Minor status vector for the species in the problem. More... | |
std::vector< size_t > | m_phaseID |
Mapping from the species number to the phase number. More... | |
std::vector< char > | m_SSPhase |
Boolean indicating whether a species belongs to a single-species phase. More... | |
std::vector< std::string > | m_speciesName |
Species string name for the kth species. More... | |
std::vector< std::string > | m_elementName |
Vector of strings containing the element names. More... | |
std::vector< int > | m_elType |
Type of the element constraint. More... | |
std::vector< int > | m_elementActive |
Specifies whether an element constraint is active. More... | |
std::vector< vcs_VolPhase * > | m_VolPhaseList |
Array of Phase Structures. Length = number of phases. More... | |
std::string | m_title |
String containing the title of the run. More... | |
char | m_unitsState |
This specifies the current state of units for the Gibbs free energy properties in the program. More... | |
double | m_totalMoleScale |
Multiplier for the mole numbers within the nondimensionless formulation. More... | |
std::vector< int > | m_actConventionSpecies |
specifies the activity convention of the phase containing the species More... | |
std::vector< int > | m_phaseActConvention |
specifies the activity convention of the phase. More... | |
std::vector< double > | m_lnMnaughtSpecies |
specifies the ln(Mnaught) used to calculate the chemical potentials More... | |
std::vector< double > | m_actCoeffSpecies_new |
Molar-based Activity Coefficients for Species. More... | |
std::vector< double > | m_actCoeffSpecies_old |
Molar-based Activity Coefficients for Species based on old mole numbers. More... | |
DoubleStarStar | m_np_dLnActCoeffdMolNum |
Change in the log of the activity coefficient with respect to the mole number multiplied by the phase mole number. More... | |
std::vector< double > | m_wtSpecies |
Molecular weight of each species. More... | |
std::vector< double > | m_chargeSpecies |
Charge of each species. Length = number of species. More... | |
std::vector< std::vector < size_t > > | phasePopProblemLists_ |
std::vector< VCS_SPECIES_THERMO * > | m_speciesThermoList |
Vector of pointers to thermostructures which identify the model and parameters for evaluating the thermodynamic functions for that particular species. More... | |
int | m_useActCoeffJac |
Choice of Hessians. More... | |
double | m_totalVol |
Total volume of all phases. Units are m^3. More... | |
std::vector< double > | m_PMVolumeSpecies |
Partial molar volumes of the species. More... | |
double | m_Faraday_dim |
dimensionless value of Faraday's constant, F / RT (1/volt) More... | |
VCS_COUNTERS * | m_VCount |
Timing and iteration counters for the vcs object. More... | |
int | m_debug_print_lvl |
Debug printing lvl. More... | |
int | m_timing_print_lvl |
printing level of timing information More... | |
int | m_VCS_UnitsFormat |
Units for the chemical potential data. More... | |
Private Member Functions | |
int | vcs_zero_species (const size_t kspec) |
Zero out the concentration of a species. More... | |
int | vcs_delete_species (const size_t kspec) |
Change a single species from active to inactive status. More... | |
bool | vcs_delete_multiphase (const size_t iph) |
This routine handles the bookkeeping involved with the deletion of multiphase phases from the problem. More... | |
int | delta_species (const size_t kspec, double *const delta_ptr) |
Change the concentration of a species by delta moles. More... | |
size_t | vcs_add_all_deleted () |
Provide an estimate for the deleted species in phases that are not zeroed out. More... | |
int | vcs_recheck_deleted () |
Recheck deleted species in multispecies phases. More... | |
bool | recheck_deleted_phase (const int iphase) |
Recheck deletion condition for multispecies phases. More... | |
double | vcs_minor_alt_calc (size_t kspec, size_t irxn, bool *do_delete, char *ANOTE) const |
Minor species alternative calculation. More... | |
bool | vcs_globStepDamp () |
This routine optimizes the minimization of the total gibbs free energy by making sure the slope of the following functional stays negative: More... | |
void | vcs_switch2D (double *const *const Jac, const size_t k1, const size_t k2) const |
Switch rows and columns of a square matrix. More... | |
double | l2normdg (double dg[]) const |
Calculate the norm of a deltaGibbs free energy vector. More... | |
void | prneav () const |
Print out and check the elemental abundance vector. More... | |
void | checkDelta1 (double *const ds, double *const delTPhMoles, int kspec) |
void | vcs_inest (double *const aw, double *const sa, double *const sm, double *const ss, double test) |
Estimate equilibrium compositions. More... | |
void | vcs_SSPhase (void) |
Calculate the status of single species phases. More... | |
double | deltaG_Recalc_Rxn (const int stateCalc, const size_t irxn, const double *const molNum, double *const ac, double *const mu_i) |
This function recalculates the deltaG for reaction, irxn. More... | |
void | vcs_delete_memory () |
Delete memory that isn't just resizable STL containers. More... | |
void | vcs_counters_init (int ifunc) |
Initialize the internal counters. More... | |
void | vcs_TCounters_report (int timing_print_lvl=1) |
Create a report on the plog file containing timing and its information. More... | |
void | vcs_setFlagsVolPhases (const bool upToDate, const int stateCalc) |
void | vcs_setFlagsVolPhase (const size_t iph, const bool upToDate, const int stateCalc) |
void | vcs_updateMolNumVolPhases (const int stateCalc) |
Update all underlying vcs_VolPhase objects. More... | |
Friends | |
class | vcs_phaseStabilitySolve |
This is the main structure used to hold the internal data used in vcs_solve_TP(), and to solve TP systems.
The indices of information in this structure may change when the species basis changes or when phases pop in and out of existence. Both of these operations change the species ordering.
Definition at line 50 of file vcs_solve.h.
void vcs_initSizes | ( | const size_t | nspecies0, |
const size_t | nelements, | ||
const size_t | nphase0 | ||
) |
Initialize the sizes within the VCS_SOLVE object.
This resizes all of the internal arrays within the object. This routine operates in two modes. If all of the parameters are the same as it currently exists in the object, nothing is done by this routine; a quick exit is carried out and all of the data in the object persists.
If any of the parameters are different than currently exists in the object, then all of the data in the object must be redone. It may not be zeroed, but it must be redone.
nspecies0 | Number of species within the object |
nelements | Number of element constraints within the problem |
nphase0 | Number of phases defined within the problem. |
Definition at line 60 of file vcs_solve.cpp.
References plogf, VCS_ELEM_TYPE_ABSPOS, VCS_SPECIES_TYPE_MOLNUM, and VCSnonideal::vcs_timing_print_lvl.
int vcs | ( | VCS_PROB * | vprob, |
int | ifunc, | ||
int | ipr, | ||
int | ip1, | ||
int | maxit | ||
) |
Solve an equilibrium problem.
This is the main interface routine to the equilibrium solver
Input:
vprob | Object containing the equilibrium Problem statement |
ifunc | Determines the operation to be done: Valid values: 0 -> Solve a new problem by initializing structures first. An initial estimate may or may not have been already determined. This is indicated in the VCS_PROB structure. 1 -> The problem has already been initialized and set up. We call this routine to resolve it using the problem statement and solution estimate contained in the VCS_PROB structure. 2 -> Don't solve a problem. Destroy all the private structures. |
ipr | Printing of results ipr = 1 -> Print problem statement and final results to standard output 0 -> don't report on anything |
ip1 | Printing of intermediate results IP1 = 1 -> Print intermediate results. |
maxit | Maximum number of iterations for the algorithm |
Output:
Definition at line 247 of file vcs_solve.cpp.
References VCS_PROB::ne, VCS_PROB::NPhase, VCS_PROB::nspecies, plogf, VCS_PROB::PresPA, clockWC::secondsWC(), VCS_PROB::T, and VCS_SUCCESS.
Referenced by vcs_MultiPhaseEquil::equilibrate_TP().
int vcs_solve_TP | ( | int | print_lvl, |
int | printDetails, | ||
int | maxit | ||
) |
Main routine that solves for equilibrium at constant T and P using a variant of the VCS method.
This is the main routine that solves for equilibrium at constant T and P using a variant of the VCS method. Nonideal phases can be accommodated as well.
Any number of single-species phases and multi-species phases can be handled by the present version.
print_lvl | 1 -> Print results to standard output; 0 -> don't report on anything |
printDetails | 1 -> Print intermediate results. |
maxit | Maximum number of iterations for the algorithm |
Definition at line 60 of file vcs_solve_TP.cpp.
References vcs_VolPhase::exists(), vcs_VolPhase::m_singleSpecies, vcs_VolPhase::molefraction(), Cantera::npos, vcs_VolPhase::nSpecies(), vcs_VolPhase::PhaseName, plogendl, plogf, clockWC::secondsWC(), VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_DELETE_PHASE_CUTOFF, VCS_ELEM_TYPE_ABSPOS, VCS_PHASE_EXIST_NO, VCSnonideal::vcs_print_line(), VCS_RELDELETE_SPECIES_CUTOFF, VCS_SPECIES_INTERFACIALVOLTAGE, VCS_SPECIES_MAJOR, VCS_SPECIES_MINOR, VCS_SPECIES_STOICHZERO, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_TYPE_MOLNUM, VCS_SPECIES_ZEROEDMS, VCS_SPECIES_ZEROEDSS, VCS_STATECALC_NEW, VCS_STATECALC_OLD, and VCS_SUCCESS.
Referenced by VCS_SOLVE::vcs_TP().
void vcs_reinsert_deleted | ( | size_t | kspec | ) |
We make decisions on the initial mole number, and major-minor status here. We also fix up the total moles in a phase.
irxn = id of the noncomponent species formation reaction for the species to be added in.
The algorithm proceeds to implement these decisions in the previous position of the species. Then, vcs_switch_pos is called to move the species into the last active species slot, incrementing the number of active species at the same time.
This routine is responsible for the global data manipulation only.
Definition at line 2277 of file vcs_solve_TP.cpp.
References vcs_VolPhase::exists(), plogf, vcs_VolPhase::setExistence(), vcs_VolPhase::setMolesFromVCSCheck(), VCS_DATA_PTR, VCS_PHASE_EXIST_NO, VCS_PHASE_EXIST_YES, VCS_RELDELETE_SPECIES_CUTOFF, VCS_SPECIES_DELETED, VCS_SPECIES_MAJOR, VCS_SPECIES_MINOR, and VCS_STATECALC_OLD.
Referenced by VCS_SOLVE::vcs_elcorr().
int vcs_basopt | ( | const bool | doJustComponents, |
double | aw[], | ||
double | sa[], | ||
double | sm[], | ||
double | ss[], | ||
double | test, | ||
bool *const | usedZeroedSpecies | ||
) |
Choose the optimum species basis for the calculations.
This is done by choosing the species with the largest mole fraction not currently a linear combination of the previous components. Then, calculate the stoichiometric coefficient matrix for that basis.
Rearranges the solution data to put the component data at the front of the species list.
Then, calculates m_stoichCoeffRxnMatrix[irxn][jcomp] the formation reactions for all noncomponent species in the mechanism. Also calculates DNG(I) and DNL(I), the net mole change for each formation reaction. Also, initializes IR(I) to the default state.
[in] | doJustComponents | If true, the m_stoichCoeffRxnMatrix[][] and m_deltaMolNumPhase[] are not calculated. |
[in] | aw | Vector of mole fractions which will be used to construct an optimal basis from. |
[in] | sa | Gram-Schmidt orthog work space (nc in length) sa[j] |
[in] | ss | Gram-Schmidt orthog work space (nc in length) ss[j] |
[in] | sm | QR matrix work space (nc*ne in length) sm[i+j*ne] |
[in] | test | This is a small negative number dependent upon whether an estimate is supplied or not. |
[out] | usedZeroedSpecies | If true, then a species with a zero concentration was used as a component. The problem may be converged. Or, the problem may have a range space error and may not have a proper solution. |
Definition at line 2888 of file vcs_solve_TP.cpp.
References Cantera::npos, plogendl, plogf, clockWC::secondsWC(), VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_ELEM_TYPE_ABSPOS, VCSnonideal::vcs_print_stringTrunc(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SUCCESS, and VCSnonideal::vcsUtil_gaussj().
Referenced by VCS_SOLVE::vcs_inest(), and VCS_SOLVE::vcs_prep_oneTime().
size_t vcs_basisOptMax | ( | const double *const | molNum, |
const size_t | j, | ||
const size_t | n | ||
) |
Choose a species to test for the next component.
We make the choice based on testing (molNum[i] * spSize[i]) for its maximum value. Preference for single species phases is also made.
molNum | Mole number vector |
j | index into molNum[] that indicates where the search will start from Previous successful components are swapped into the front of molNum[]. |
n | Length of molNum[] |
Definition at line 3468 of file vcs_solve_TP.cpp.
int vcs_species_type | ( | const size_t | kspec | ) | const |
Evaluate the species category for the indicated species.
All evaluations are done using the "old" version of the solution.
kspec | Species to be evaluated |
Definition at line 3523 of file vcs_solve_TP.cpp.
References vcs_VolPhase::exists(), plogf, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_ELEM_TYPE_ABSPOS, VCS_PHASE_EXIST_YES, VCS_PHASE_EXIST_ZEROEDPHASE, VCS_SPECIES_ACTIVEBUTZERO, VCS_SPECIES_INTERFACIALVOLTAGE, VCS_SPECIES_MAJOR, VCS_SPECIES_MINOR, VCS_SPECIES_STOICHZERO, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_ZEROEDMS, VCS_SPECIES_ZEROEDPHASE, and VCS_SPECIES_ZEROEDSS.
bool vcs_evaluate_speciesType | ( | ) |
This routine evaluates the species type for all species.
Definition at line 4276 of file vcs_solve_TP.cpp.
References plogendl, plogf, VCS_SPECIES_ACTIVEBUTZERO, VCS_SPECIES_COMPONENT, VCS_SPECIES_DELETED, VCS_SPECIES_INTERFACIALVOLTAGE, VCS_SPECIES_MAJOR, VCS_SPECIES_MINOR, VCS_SPECIES_STOICHZERO, VCS_SPECIES_ZEROEDMS, VCS_SPECIES_ZEROEDPHASE, VCS_SPECIES_ZEROEDSS, and VCSnonideal::vcs_speciesType_string().
void vcs_chemPotPhase | ( | const int | stateCalc, |
const size_t | iph, | ||
const double *const | molNum, | ||
double *const | ac, | ||
double *const | mu_i, | ||
const bool | do_deleted = false |
||
) |
We calculate the dimensionless chemical potentials of all species in a single phase.
Note, for multispecies phases which are currently zeroed out, the chemical potential is filled out with the standard chemical potential.
For species in multispecies phases whose concentration is zero, we need to set the mole fraction to a very low value. Its chemical potential is then calculated using the VCS_DELETE_MINORSPECIES_CUTOFF concentration to keep numbers positive.
m_feSpecies(I) = m_SSfeSpecies(I) + ln(z(I)) - ln(m_tPhaseMoles[iph]) + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase];
(This is equivalent to the adding the log of the mole fraction onto the standard chemical potential. )
m_feSpecies(I) = m_SSfeSpecies(I) + ln(ActCoeff[I] * z(I)) - ln(m_tPhaseMoles[iph]) + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase];
( This is equivalent to the adding the log of the mole fraction multiplied by the activity coefficient onto the standard chemical potential. )
m_feSpecies(I) = m_SSfeSpecies(I) + ln(ActCoeff[I] * z(I)) - ln(m_tPhaseMoles[iph]) - ln(Mnaught * m_units) + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase];
Note: m_SSfeSpecies(I)
is the molality based standard state. However, ActCoeff[I]
is the molar based activity coefficient We have used the formulas:
ActCoeff_M[I] = ActCoeff[I] / Xmol[N]
where Xmol[N]
is the mole fraction of the solvent and ActCoeff_M[I]
is the molality based act coeff.
Note: This is equivalent to the "normal" molality formulation:
m_feSpecies(I) = m_SSfeSpecies(I) + ln(ActCoeff_M[I] * m(I)) + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase]
where m[I]
is the molality of the ith solute
m[I] = Xmol[I] / ( Xmol[N] * Mnaught * m_units)
z(I)/tPhMoles_ptr[iph] = Xmol[i]
is the mole fraction of i in the phase.
NOTE: As per the discussion in vcs_dfe(), for small species where the mole fraction is small:
z(i) < VCS_DELETE_MINORSPECIES_CUTOFF
The chemical potential is calculated as:
m_feSpecies(I) = m_SSfeSpecies(I) + ln(ActCoeff[i](VCS_DELETE_MINORSPECIES_CUTOFF))
[in] | iph | Phase to be calculated |
[in] | molNum(i) | Number of moles of species i (VCS species order) |
[out] | ac | Activity coefficients for species in phase (VCS species order) |
[out] | mu_i | Dimensionless chemical potentials for phase species (VCS species order) |
Definition at line 3722 of file vcs_solve_TP.cpp.
References vcs_VolPhase::electricPotential(), vcs_VolPhase::nSpecies(), plogf, vcs_VolPhase::sendToVCS_ActCoeff(), vcs_VolPhase::setMolesFromVCS(), vcs_VolPhase::spGlobalIndexVCS(), VCS_DELETE_MINORSPECIES_CUTOFF, VCS_SPECIES_DELETED, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
Referenced by VCS_SOLVE::deltaG_Recalc_Rxn().
void vcs_dfe | ( | const int | stateCalc, |
const int | ll, | ||
const size_t | lbot, | ||
const size_t | ltop | ||
) |
Calculate the dimensionless chemical potentials of all species or of certain groups of species, at a fixed temperature and pressure.
We calculate the dimensionless chemical potentials of all species or certain groups of species here, at a fixed temperature and pressure, for the input mole vector z[] in the parameter list. Nondimensionalization is achieved by division by RT.
Note, for multispecies phases which are currently zeroed out, the chemical potential is filled out with the standard chemical potential.
For species in multispecies phases whose concentration is zero, we need to set the mole fraction to a very low value. Its chemical potential is then calculated using the VCS_DELETE_MINORSPECIES_CUTOFF concentration to keep numbers positive.
For formulas, see vcs_chemPotPhase().
As per the discussion above, for small species where the mole fraction
z(i) < VCS_DELETE_MINORSPECIES_CUTOFF
The chemical potential is calculated as:
m_feSpecies(I)(I) = m_SSfeSpecies(I) + ln(ActCoeff[i](VCS_DELETE_MINORSPECIES_CUTOFF))
Species in the following categories are treated as "small species"
For species in multispecies phases which are currently not active, the treatment is different. These species are in the following species categories:
For these species, the ln( ActCoeff[I] X[I])
term is dropped altogether. The following equation is used:
m_feSpecies(I) = m_SSfeSpecies(I) + Charge[I] * Faraday_dim * phasePhi[iphase];
These species have species types of VCS_SPECIES_TYPE_INTERFACIALVOLTAGE The chemical potentials for these "species" refer to electrons in metal electrodes. They have the following formula
m_feSpecies(I) = m_SSfeSpecies(I) - F z[I] / RT
F
is Faraday's constant.R
= gas constantT
= temperatureV
= potential of the interface = phi_electrode - phi_solutionFor these species, the solution vector unknown, z[I], is V, the phase voltage, in volts.
ll | Determine which group of species gets updated
|
lbot | Restricts the calculation of the chemical potential to the species between LBOT <= i < LTOP. Usually LBOT and LTOP will be equal to 0 and MR, respectively. |
ltop | Top value of the loops |
stateCalc | Determines whether z is old or new or tentative:
|
Also needed: ff : standard state chemical potentials. These are the chemical potentials of the standard states at the same T and P as the solution. tg : Total Number of moles in the phase.
Definition at line 3788 of file vcs_solve_TP.cpp.
References vcs_VolPhase::electricPotential(), vcs_VolPhase::m_singleSpecies, plogendl, plogf, vcs_VolPhase::sendToVCS_ActCoeff(), vcs_VolPhase::updateFromVCS_MoleNumbers(), VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_DIMENSIONAL_G, VCSnonideal::vcs_doubleEqual(), VCS_SPECIES_MINOR, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_ZEROEDMS, VCS_SPECIES_ZEROEDPHASE, VCS_STATECALC_NEW, and VCS_STATECALC_OLD.
Referenced by VCS_SOLVE::vcs_inest(), and VCS_SOLVE::vcs_report().
void vcs_printSpeciesChemPot | ( | const int | stateCalc | ) | const |
Print out a table of chemical potentials.
stateCalc | Determines where to get the mole numbers from.
|
Definition at line 4062 of file vcs_solve_TP.cpp.
References vcs_VolPhase::electricPotential(), Cantera::GasConstant, Cantera::int2str(), VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCSnonideal::vcs_print_line(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_ZEROEDMS, VCS_SPECIES_ZEROEDPHASE, VCS_SPECIES_ZEROEDSS, and VCS_STATECALC_NEW.
void vcs_updateVP | ( | const int | stateCalc | ) |
This routine uploads the state of the system into all of the vcs_VolumePhase objects in the current problem.
stateCalc | Determines where to get the mole numbers from.
|
Definition at line 4252 of file vcs_solve_TP.cpp.
References plogendl, plogf, vcs_VolPhase::setMolesFromVCSCheck(), VCS_DATA_PTR, VCS_STATECALC_NEW, and VCS_STATECALC_OLD.
bool vcs_popPhasePossible | ( | const size_t | iphasePop | ) | const |
Utility function that evaluates whether a phase can be popped into existence.
A phase can be popped iff the stoichiometric coefficients for the component species, whose concentrations will be lowered during the process, are positive by at least a small degree.
If one of the phase species is a zeroed component, then the phase can be popped if the component increases in mole number as the phase moles are increased.
iphasePop | id of the phase, which is currently zeroed, |
Definition at line 19 of file vcs_phaseStability.cpp.
References vcs_VolPhase::exists(), vcs_VolPhase::nSpecies(), vcs_VolPhase::spGlobalIndexVCS(), VCS_DELETE_ELEMENTABS_CUTOFF, and VCS_ELEM_TYPE_ABSPOS.
int vcs_phasePopDeterminePossibleList | ( | ) |
Determine the list of problems that need to be checked to see if there are any phases pops.
This routine evaluates and fills in #phasePopProblemLists_. Need to work in species that are zeroed by element constraints.
Definition at line 115 of file vcs_phaseStability.cpp.
References vcs_VolPhase::exists(), vcs_VolPhase::nSpecies(), vcs_VolPhase::spGlobalIndexVCS(), and VCS_ELEM_TYPE_ABSPOS.
size_t vcs_popPhaseID | ( | std::vector< size_t > & | phasePopPhaseIDs | ) |
Decision as to whether a phase pops back into existence.
phasePopPhaseIDs | Vector containing the phase ids of the phases that will be popped this step. |
Definition at line 242 of file vcs_phaseStability.cpp.
References vcs_VolPhase::exists(), vcs_VolPhase::m_singleSpecies, Cantera::npos, vcs_VolPhase::PhaseName, plogf, and vcs_VolPhase::spGlobalIndexVCS().
int vcs_popPhaseRxnStepSizes | ( | const size_t | iphasePop | ) |
Calculates the deltas of the reactions due to phases popping into existence.
Updates m_deltaMolNumSpecies[irxn] : reaction adjustments, where irxn refers to the irxn'th species formation reaction. This adjustment is for species irxn + M, where M is the number of components.
iphasePop | Phase id of the phase that will come into existence |
Definition at line 373 of file vcs_phaseStability.cpp.
References vcs_VolPhase::creationMoleNumbers(), DATA_PTR, vcs_VolPhase::exists(), vcs_VolPhase::m_singleSpecies, vcs_VolPhase::nSpecies(), vcs_VolPhase::PhaseName, plogf, vcs_VolPhase::spGlobalIndexVCS(), VCS_DELETE_PHASE_CUTOFF, VCS_ELEM_TYPE_ABSPOS, VCS_SPECIES_COMPONENT, VCS_SPECIES_MAJOR, and VCS_SPECIES_MINOR.
size_t vcs_RxnStepSizes | ( | int & | forceComponentCalc, |
size_t & | kSpecial | ||
) |
Calculates formation reaction step sizes.
This is equation 6.4-16, p. 143 in Smith and Missen.
m_deltaMolNumSpecies(irxn) : reaction adjustments, where irxn refers to the irxn'th species formation reaction. This adjustment is for species irxn + M, where M is the number of components.
Special branching occurs sometimes. This causes the component basis to be reevaluated
forceComponentCalc | integer flagging whether a component recalculation needs to be carried out. |
kSpecial | species number of phase being zeroed. |
Definition at line 20 of file vcs_rxnadj.cpp.
References vcs_VolPhase::exists(), VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_deltaGRxn_new, VCS_SOLVE::m_deltaMolNumPhase, VCS_SOLVE::m_deltaMolNumSpecies, VCS_SOLVE::m_indexRxnToSpecies, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numRxnRdc, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_phaseID, vcs_VolPhase::m_singleSpecies, VCS_SOLVE::m_speciesName, VCS_SOLVE::m_speciesStatus, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_SSPhase, VCS_SOLVE::m_stoichCoeffRxnMatrix, VCS_SOLVE::m_tolmaj2, VCS_SOLVE::m_totalMolNum, VCS_SOLVE::m_tPhaseMoles_old, VCS_SOLVE::m_useActCoeffJac, VCS_SOLVE::m_VolPhaseList, Cantera::npos, plogendl, plogf, VCS_SOLVE::vcs_CalcLnActCoeffJac(), VCS_DATA_PTR, VCS_DELETE_PHASE_CUTOFF, VCS_SOLVE::vcs_Hessian_diag_adj(), VCSnonideal::vcs_print_line(), VCS_SMALL_MULTIPHASE_SPECIES, VCS_SPECIES_MAJOR, VCS_SPECIES_STOICHZERO, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_ZEROEDPHASE, and VCSnonideal::vcs_speciesType_string().
double vcs_tmoles | ( | ) |
Calculates the total number of moles of species in all phases.
Also updates the variable m_totalMolNum and Reconciles Phase existence flags with total moles in each phase.
Definition at line 4198 of file vcs_solve_TP.cpp.
References vcs_VolPhase::setTotalMoles(), and VCS_SPECIES_TYPE_MOLNUM.
Referenced by VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().
void vcs_deltag | ( | const int | l, |
const bool | doDeleted, | ||
const int | vcsState, | ||
const bool | alterZeroedPhases = true |
||
) |
This subroutine calculates reaction free energy changes for all noncomponent formation reactions.
Formation reactions are reactions which create each noncomponent species from the component species. m_stoichCoeffRxnMatrix[irxn][jcomp] are the stoichiometric coefficients for these reactions. A stoichiometric coefficient of one is assumed for species irxn in this reaction.
l |
|
doDeleted | Do deleted species |
vcsState | Calculate deltaG corresponding to either old or new free energies |
alterZeroedPhases | boolean indicating whether we should add in a special section for zeroed phases. |
Note we special case one important issue. If the component has zero moles, then we do not allow deltaG < 0.0 for formation reactions which would lead to the loss of more of that same component. This dG < 0.0 condition feeds back into the algorithm in several places, and leads to a infinite loop in at least one case.
Definition at line 4380 of file vcs_solve_TP.cpp.
References Cantera::checkFinite(), vcs_VolPhase::m_singleSpecies, vcs_VolPhase::nSpecies(), plogf, vcs_VolPhase::spGlobalIndexVCS(), VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_SPECIES_MINOR, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_STATECALC_NEW, and VCS_STATECALC_OLD.
Referenced by VCS_SOLVE::vcs_inest().
void vcs_deltag_Phase | ( | const size_t | iphase, |
const bool | doDeleted, | ||
const int | stateCalc, | ||
const bool | alterZeroedPhases = true |
||
) |
Calculate deltag of formation for all species in a single phase.
Calculate deltag of formation for all species in a single phase. It is assumed that the fe[] is up to date for all species. However, if the phase is currently zeroed out, a subproblem is calculated to solve for AC[i] and pseudo-X[i] for that phase.
iphase | phase index of the phase to be calculated |
doDeleted | boolean indicating whether to do deleted species or not |
stateCalc | integer describing which set of free energies to use and where to stick the results. |
alterZeroedPhases | boolean indicating whether we should add in a special section for zeroed phases. |
NOTE: this is currently not used used anywhere. It may be in the future?
Definition at line 4725 of file vcs_solve_TP.cpp.
References vcs_VolPhase::m_singleSpecies, plogendl, plogf, vcs_VolPhase::spGlobalIndexVCS(), VCS_DATA_PTR, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_STATECALC_NEW, and VCS_STATECALC_OLD.
void vcs_switch_pos | ( | const bool | ifunc, |
const size_t | k1, | ||
const size_t | k2 | ||
) |
Swaps the indices for all of the global data for two species, k1 and k2.
ifunc | If true, switch the species data and the noncomponent reaction data. This must be called for a non-component species only. If false, switch the species data only. Typically, we use this option when determining the component species and at the end of the calculation, when we want to return unscrambled results. All rxn data will be out-of-date. |
k1 | First species index |
k2 | Second species index |
Definition at line 4879 of file vcs_solve_TP.cpp.
References plogf, vcs_VolPhase::setSpGlobalIndexVCS(), and vcs_VolPhase::spGlobalIndexVCS().
Referenced by VCS_SOLVE::vcs_rearrange().
double vcs_birthGuess | ( | const int | kspec | ) |
Birth guess returns the number of moles of a species that is coming back to life.
Birth guess returns the number of moles of a species that is coming back to life. Note, this routine is not applicable if the whole phase is coming back to life, not just one species in that phase.
Do a minor alt calculation. But, cap the mole numbers at 1.0E-15. For SS phases use VCS_DELETE_SPECIES_CUTOFF * 100.
The routine makes sure the guess doesn't reduce the concentration of a component by more than 1/3. Note this may mean that the vlaue coming back from this routine is zero or a very small number.
kspec | Species number that is coming back to life |
Definition at line 4984 of file vcs_solve_TP.cpp.
References plogf, VCS_DELETE_MINORSPECIES_CUTOFF, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
int vcs_solve_phaseStability | ( | const int | iphase, |
int | ifunc, | ||
double & | funcval, | ||
int | print_lvl | ||
) |
Routine that independently determines whether a phase should be popped under the current conditions.
Definition at line 159 of file vcs_solve_phaseStability.cpp.
References VCS_DATA_PTR, and VCS_STATECALC_OLD.
double vcs_phaseStabilityTest | ( | const size_t | iph | ) |
Main program to test whether a deleted phase should be brought back into existence.
iph | Phase id of the deleted phase |
Definition at line 569 of file vcs_phaseStability.cpp.
References vcs_VolPhase::creationMoleNumbers(), Cantera::npos, vcs_VolPhase::nSpecies(), plogf, vcs_VolPhase::sendToVCS_ActCoeff(), vcs_VolPhase::setCreationMoleNumbers(), vcs_VolPhase::setMoleFractionsState(), vcs_VolPhase::spGlobalIndexVCS(), VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCSnonideal::vcs_l2norm(), VCS_STATECALC_OLD, and VCS_STATECALC_PHASESTABILITY.
int vcs_TP | ( | int | ipr, |
int | ip1, | ||
int | maxit, | ||
double | T, | ||
double | pres | ||
) |
Solve an equilibrium problem at a particular fixed temperature and pressure.
The actual problem statement is assumed to be in the structure already. This is a wrapper around the solve_TP() function. In this wrapper, we nondimensionalize the system we calculate the standard state gibbs free energies of the species, and we decide whether to we need to use the initial guess algorithm.
ipr | = 1 -> Print results to standard output; 0 -> don't report on anything |
ip1 | = 1 -> Print intermediate results; 0 -> Dont print any intermediate results |
maxit | Maximum number of iterations for the algorithm |
T | Value of the Temperature (Kelvin) |
pres | Value of the Pressure (units given by m_VCS_UnitsFormat variable |
Definition at line 9 of file vcs_TP.cpp.
References VCS_SOLVE::m_doEstimateEquil, VCS_SOLVE::m_pressurePA, VCS_SOLVE::m_temperature, plogf, VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_fePrep_TP(), VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_solve_TP(), and VCS_SUCCESS.
int vcs_evalSS_TP | ( | int | ipr, |
int | ip1, | ||
double | Temp, | ||
double | pres | ||
) |
Evaluate the standard state free energies at the current temperature and pressure. Ideal gas pressure contribution is added in here.
ipr | 1 -> Print results to standard output; 0 -> don't report on anything |
ip1 | 1 -> Print intermediate results; 0 -> don't. |
Temp | Temperature (Kelvin) |
pres | Pressure (Pascal) |
Definition at line 61 of file vcs_TP.cpp.
References VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_pressurePA, VCS_SOLVE::m_SSfeSpecies, VCS_SOLVE::m_temperature, VCS_SOLVE::m_VCS_UnitsFormat, VCS_SOLVE::m_VolPhaseList, vcs_VolPhase::sendToVCS_GStar(), vcs_VolPhase::setState_TP(), VCS_DATA_PTR, and VCS_SUCCESS.
Referenced by VCS_SOLVE::vcs_prep_oneTime(), and VCS_SOLVE::vcs_TP().
void vcs_fePrep_TP | ( | void | ) |
Initialize the chemical potential of single species phases.
For single species phases, initialize the chemical potential with the value of the standard state chemical potential. This value doesn't change during the calculation
Definition at line 109 of file vcs_TP.cpp.
References VCS_SOLVE::m_feSpecies_new, VCS_SOLVE::m_feSpecies_old, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_SSfeSpecies, and VCS_SOLVE::m_SSPhase.
Referenced by VCS_SOLVE::vcs_TP().
double vcs_VolTotal | ( | const double | tkelvin, |
const double | pres, | ||
const double | w[], | ||
double | volPM[] | ||
) |
Calculation of the total volume and the partial molar volumes.
This function calculates the partial molar volume for all species, kspec, in the thermo problem at the temperature TKelvin and pressure, Pres, pres is in atm. And, it calculates the total volume of the combined system.
[in] | tkelvin | Temperature in kelvin() |
[in] | pres | Pressure in Pascal |
[in] | w | w[] is the vector containing the current mole numbers in units of kmol. |
[out] | volPM[] | For species in all phase, the entries are the partial molar volumes units of M**3 / kmol. |
Definition at line 990 of file vcs_solve.cpp.
References vcs_VolPhase::sendToVCS_VolPM(), vcs_VolPhase::setMolesFromVCS(), vcs_VolPhase::setState_TP(), and VCS_STATECALC_OLD.
Referenced by VCS_SOLVE::vcs_report().
int vcs_prep_oneTime | ( | int | printLvl | ) |
This routine is mostly concerned with changing the private data to be consistent with what's needed for solution.
It is called one time for each new problem structure definition.
This routine is always followed by vcs_prep(). Therefore, tasks that need to be done for every call to vcsc() should be placed in vcs_prep() and not in this routine.
The problem structure refers to:
Tasks:
printLvl | Print level of the routine |
Definition at line 65 of file vcs_prep.cpp.
References vcs_SpeciesProperties::FormulaMatrixCol, VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_doEstimateEquil, VCS_SOLVE::m_indexRxnToSpecies, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numRxnRdc, VCS_SOLVE::m_numRxnTot, VCS_SOLVE::m_numSpeciesRdc, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_phaseID, VCS_SOLVE::m_pressurePA, VCS_SOLVE::m_speciesLocalPhaseIndex, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_spSize, VCS_SOLVE::m_SSfeSpecies, VCS_SOLVE::m_temperature, VCS_SOLVE::m_VolPhaseList, plogf, vcs_VolPhase::speciesProperty(), VCS_SOLVE::vcs_basopt(), VCS_DATA_PTR, VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_evalSS_TP(), VCS_SPECIES_TYPE_MOLNUM, VCS_SOLVE::vcs_SSPhase(), and VCS_SUCCESS.
int vcs_prep | ( | ) |
Prepare the object for solution.
This routine is mostly concerned with changing the private data to be consistent with that needed for solution. It is called for every invocation of the vcs_solve() except for the cleanup invocation.
Tasks:
Definition at line 219 of file vcs_prep.cpp.
References VCS_SOLVE::m_deltaMolNumPhase, VCS_SOLVE::m_deltaPhaseMoles, VCS_SOLVE::m_feSpecies_new, VCS_SOLVE::m_feSpecies_old, VCS_SOLVE::m_molNumSpecies_new, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_phaseParticipation, VCS_SOLVE::m_tPhaseMoles_new, VCS_DATA_PTR, VCS_SUCCESS, and VCS_SOLVE::vcs_tmoles().
bool vcs_wellPosed | ( | VCS_PROB * | vprob | ) |
In this routine, we check for things that will cause the algorithm to fail.
We check to see if the problem is well posed. If it is not, we return false and print out error conditions.
Current there is one condition. If all the element abundances are zero, the algorithm will fail.
vprob | VCS_PROB pointer to the definition of the equilibrium problem |
Definition at line 238 of file vcs_prep.cpp.
References VCS_PROB::gai, VCS_PROB::ne, and plogf.
int vcs_elem_rearrange | ( | double *const | aw, |
double *const | sa, | ||
double *const | sm, | ||
double *const | ss | ||
) |
Rearrange the constraint equations represented by the Formula Matrix so that the operational ones are in the front.
This subroutine handles the rearrangement of the constraint equations represented by the Formula Matrix. Rearrangement is only necessary when the number of components is less than the number of elements. For this case, some constraints can never be satisfied exactly, because the range space represented by the Formula Matrix of the components can't span the extra space. These constraints, which are out of the range space of the component Formula matrix entries, are migrated to the back of the Formula matrix.
A prototypical example is an extra element column in FormulaMatrix[], which is identically zero. For example, let's say that argon is has an element column in FormulaMatrix[], but no species in the mechanism actually contains argon. Then, nc < ne. Also, without perturbation of FormulaMatrix[] vcs_basopt[] would produce a zero pivot because the matrix would be singular (unless the argon element column was already the last column of FormulaMatrix[].
This routine borrows heavily from vcs_basopt's algorithm. It finds nc constraints which span the range space of the Component Formula matrix, and assigns them as the first nc components in the formula matrix. This guarantees that vcs_basopt[] has a nonsingular matrix to invert.
Other Variables
aw | Mole fraction work space (ne in length) |
sa | Gram-Schmidt orthog work space (ne in length) |
sm | QR matrix work space (ne*ne in length) |
ss | Gram-Schmidt orthog work space (ne in length) |
Definition at line 19 of file vcs_elem_rearrange.cpp.
References VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_elementActive, VCS_SOLVE::m_elementName, VCS_SOLVE::m_formulaMatrix, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numElemConstraints, Cantera::npos, plogendl, plogf, VCS_SUCCESS, and VCS_SOLVE::vcs_switch_elem_pos().
Referenced by VCS_SOLVE::vcs_prep_oneTime().
void vcs_switch_elem_pos | ( | size_t | ipos, |
size_t | jpos | ||
) |
Swaps the indices for all of the global data for two elements, ipos and jpos.
This function knows all of the element information with VCS_SOLVE, and can therefore switch element positions
ipos | first global element index |
jpos | second global element index |
Definition at line 179 of file vcs_elem_rearrange.cpp.
References vcs_VolPhase::elemGlobalIndex(), VCS_SOLVE::m_elemAbundances, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_elementActive, VCS_SOLVE::m_elementMapIndex, VCS_SOLVE::m_elementName, VCS_SOLVE::m_elType, VCS_SOLVE::m_formulaMatrix, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_VolPhaseList, vcs_VolPhase::nElemConstraints(), plogendl, plogf, and vcs_VolPhase::setElemGlobalIndex().
Referenced by VCS_SOLVE::vcs_elem_rearrange().
int vcs_rxn_adj_cg | ( | void | ) |
Calculates reaction adjustments using a full Hessian approximation.
This does what equation 6.4-16, p. 143 in Smith and Missen is supposed to do. However, a full matrix is formed and then solved via a conjugate gradient algorithm. No preconditioning is done.
If special branching is warranted, then the program bails out.
DS(I) : reaction adjustment, where I refers to the Ith species Special branching occurs sometimes. This causes the component basis to be reevaluated return = 0 : normal return 1 : A single species phase species has been zeroed out in this routine. The species is a noncomponent 2 : Same as one but, the zeroed species is a component.
Special attention is taken to flag cases where the direction of the update is contrary to the steepest descent rule. This is an important attribute of the regular vcs algorithm. We don't want to violate this.
NOTE: currently this routine is not used.
Definition at line 383 of file vcs_rxnadj.cpp.
References VCS_SOLVE::m_deltaGRxn_new, VCS_SOLVE::m_deltaMolNumPhase, VCS_SOLVE::m_deltaMolNumSpecies, VCS_SOLVE::m_indexRxnToSpecies, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numRxnMinorZeroed, VCS_SOLVE::m_numRxnRdc, VCS_SOLVE::m_phaseID, VCS_SOLVE::m_speciesName, VCS_SOLVE::m_speciesStatus, VCS_SOLVE::m_SSPhase, VCS_SOLVE::m_stoichCoeffRxnMatrix, VCS_SOLVE::m_tolmaj2, VCS_SOLVE::m_tPhaseMoles_old, VCS_SOLVE::m_VolPhaseList, plogf, VCS_SPECIES_MAJOR, and VCS_SPECIES_MINOR.
double vcs_Hessian_diag_adj | ( | size_t | irxn, |
double | hessianDiag_Ideal | ||
) |
Calculates the diagonal contribution to the Hessian due to the dependence of the activity coefficients on the mole numbers.
(See framemaker notes, Eqn. 20 - VCS Equations document)
We allow the diagonal to be increased positively to any degree. We allow the diagonal to be decreased to 1/3 of the ideal solution value, but no more -> it must remain positive.
NOTE: currently this routine is not used
Definition at line 593 of file vcs_rxnadj.cpp.
References plogf, and VCS_SOLVE::vcs_Hessian_actCoeff_diag().
Referenced by VCS_SOLVE::vcs_RxnStepSizes().
double vcs_Hessian_actCoeff_diag | ( | size_t | irxn | ) |
Calculates the diagonal contribution to the Hessian due to the dependence of the activity coefficients on the mole numbers.
(See framemaker notes, Eqn. 20 - VCS Equations document)
NOTE: currently this routine is not used
Definition at line 611 of file vcs_rxnadj.cpp.
References VCS_SOLVE::m_indexRxnToSpecies, VCS_SOLVE::m_np_dLnActCoeffdMolNum, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_phaseID, VCS_SOLVE::m_SSPhase, VCS_SOLVE::m_stoichCoeffRxnMatrix, and VCS_SOLVE::m_tPhaseMoles_old.
Referenced by VCS_SOLVE::vcs_Hessian_diag_adj().
void vcs_CalcLnActCoeffJac | ( | const double *const | moleSpeciesVCS | ) |
Recalculate all of the activity coefficients in all of the phases based on input mole numbers.
Definition at line 650 of file vcs_rxnadj.cpp.
References DoubleStarStar::baseDataAddr(), vcs_VolPhase::isIdealSoln(), VCS_SOLVE::m_np_dLnActCoeffdMolNum, VCS_SOLVE::m_numPhases, vcs_VolPhase::m_singleSpecies, VCS_SOLVE::m_VolPhaseList, vcs_VolPhase::setMolesFromVCS(), and VCS_STATECALC_OLD.
Referenced by VCS_SOLVE::vcs_RxnStepSizes().
double vcs_line_search | ( | const size_t | irxn, |
const double | dx_orig, | ||
char *const | ANOTE | ||
) |
A line search algorithm is carried out on one reaction.
In this routine we carry out a rough line search algorithm to make sure that the m_deltaGRxn_new doesn't switch signs prematurely.
irxn | Reaction number |
dx_orig | Original step length |
ANOTE | Output character string stating the conclusions of the line search |
Definition at line 694 of file vcs_rxnadj.cpp.
References VCS_DATA_PTR, VCS_STATECALC_NEW, and VCS_STATECALC_OLD.
int vcs_report | ( | int | iconv | ) |
Print out a report on the state of the equilibrium problem to standard output.
iconv | Indicator of convergence, to be printed out in the report:
|
Definition at line 29 of file vcs_report.cpp.
References VCS_COUNTERS::Basis_Opts, VCS_COUNTERS::Its, VCS_SOLVE::m_actCoeffSpecies_old, VCS_SOLVE::m_chargeSpecies, VCS_SOLVE::m_deltaGRxn_new, VCS_SOLVE::m_elemAbundances, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_elementActive, VCS_SOLVE::m_elementName, VCS_SOLVE::m_elType, VCS_SOLVE::m_Faraday_dim, VCS_SOLVE::m_feSpecies_old, VCS_SOLVE::m_indexRxnToSpecies, VCS_SOLVE::m_lnMnaughtSpecies, VCS_SOLVE::m_molNumSpecies_new, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numRxnTot, VCS_SOLVE::m_numSpeciesRdc, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_phaseID, VCS_SOLVE::m_phasePhi, VCS_SOLVE::m_PMVolumeSpecies, VCS_SOLVE::m_pressurePA, VCS_SOLVE::m_speciesName, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_SSfeSpecies, VCS_SOLVE::m_stoichCoeffRxnMatrix, VCS_SOLVE::m_temperature, VCS_SOLVE::m_timing_print_lvl, VCS_SOLVE::m_totalMoleScale, VCS_SOLVE::m_totalVol, VCS_SOLVE::m_tPhaseMoles_old, VCS_SOLVE::m_unitsState, VCS_SOLVE::m_VCount, VCS_SOLVE::m_VCS_UnitsFormat, VCS_SOLVE::m_VolPhaseList, vcs_VolPhase::PhaseName, plogf, VCS_COUNTERS::Time_basopt, VCS_COUNTERS::Time_vcs_TP, vcs_VolPhase::totalMoles(), VCS_SOLVE::TPhInertMoles, VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_SOLVE::vcs_dfe(), VCS_DIMENSIONAL_G, VCSnonideal::vcs_doubleEqual(), VCS_SOLVE::vcs_elabPhase(), VCS_SOLVE::vcs_GibbsPhase(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_nondimMult_TP(), VCSnonideal::vcs_optMax(), VCS_SOLVE::vcs_printChemPotUnits(), VCS_SOLVE::vcs_redim_TP(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_TYPE_MOLNUM, VCS_STATECALC_OLD, VCS_SUCCESS, VCS_SOLVE::vcs_tmoles(), VCS_SOLVE::vcs_Total_Gibbs(), and VCS_SOLVE::vcs_VolTotal().
int vcs_rearrange | ( | ) |
Switch all species data back to the original order.
This destroys the data based on reaction ordering.
Definition at line 16 of file vcs_rearrange.cpp.
References VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_speciesMapIndex, and VCS_SOLVE::vcs_switch_pos().
double vcs_nondim_Farad | ( | int | mu_units, |
double | TKelvin | ||
) | const |
Returns the multiplier for electric charge terms.
Definition at line 19 of file vcs_nondim.cpp.
References Cantera::Avogadro, Cantera::GasConstant, plogendl, and plogf.
Referenced by VCS_SOLVE::vcs_nondim_TP().
double vcs_nondimMult_TP | ( | int | mu_units, |
double | TKelvin | ||
) | const |
Returns the multiplier for the nondimensionalization of the equations.
This is basically equal to RT
mu_units | integer representing the dimensional units system |
TKelvin | double Temperature in Kelvin |
Definition at line 46 of file vcs_nondim.cpp.
References Cantera::GasConst_cal_mol_K, Cantera::GasConstant, plogendl, and plogf.
Referenced by VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().
void vcs_nondim_TP | ( | ) |
Nondimensionalize the problem data.
Nondimensionalize the free energies using the divisor, R * T
Essentially the internal data can either be in dimensional form or in nondimensional form. This routine switches the data from dimensional form into nondimensional form.
Definition at line 76 of file vcs_nondim.cpp.
References Cantera::fp2str(), VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_deltaGRxn_new, VCS_SOLVE::m_deltaGRxn_old, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_elType, VCS_SOLVE::m_Faraday_dim, VCS_SOLVE::m_feSpecies_old, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_SSfeSpecies, VCS_SOLVE::m_temperature, VCS_SOLVE::m_totalMoleScale, VCS_SOLVE::m_unitsState, VCS_SOLVE::m_VCS_UnitsFormat, VCS_SOLVE::m_VolPhaseList, plogendl, plogf, vcs_VolPhase::setTotalMolesInert(), VCS_SOLVE::TPhInertMoles, VCS_DIMENSIONAL_G, VCS_ELEM_TYPE_ABSPOS, VCS_SOLVE::vcs_nondim_Farad(), VCS_NONDIMENSIONAL_G, VCS_SOLVE::vcs_nondimMult_TP(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, and VCS_SOLVE::vcs_tmoles().
Referenced by VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_TP().
void vcs_redim_TP | ( | void | ) |
Redimensionalize the problem data.
Reddimensionalize the free energies using the multiplier R * T
Essentially the internal data can either be in dimensional form or in nondimensional form. This routine switches the data from nondimensional form into dimensional form.
Definition at line 167 of file vcs_nondim.cpp.
References VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_deltaGRxn_new, VCS_SOLVE::m_deltaGRxn_old, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_Faraday_dim, VCS_SOLVE::m_feSpecies_old, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_SSfeSpecies, VCS_SOLVE::m_temperature, VCS_SOLVE::m_totalMoleScale, VCS_SOLVE::m_unitsState, VCS_SOLVE::m_VCS_UnitsFormat, VCS_SOLVE::m_VolPhaseList, plogendl, plogf, vcs_VolPhase::setTotalMolesInert(), VCS_SOLVE::TPhInertMoles, VCS_DIMENSIONAL_G, VCS_SOLVE::vcs_nondimMult_TP(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, and VCS_SOLVE::vcs_tmoles().
Referenced by VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_TP().
void vcs_printChemPotUnits | ( | int | unitsFormat | ) | const |
Print the string representing the Chemical potential units.
This gets printed using plogf()
unitsFormat | Integer representing the units system |
Definition at line 214 of file vcs_nondim.cpp.
References plogf.
Referenced by VCS_SOLVE::vcs_report().
void vcs_elab | ( | ) |
Computes the current elemental abundances vector.
Computes the elemental abundances vector, m_elemAbundances[], and stores it back into the global structure
Definition at line 14 of file vcs_elem.cpp.
References VCS_SOLVE::m_elemAbundances, VCS_SOLVE::m_formulaMatrix, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_speciesUnknownType, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
Referenced by VCS_SOLVE::vcs_elcorr(), and VCS_SOLVE::vcs_inest_TP().
bool vcs_elabcheck | ( | int | ibound | ) |
Checks to see if the element abundances are in compliance. If they are, then TRUE is returned. If not, FALSE is returned. Note the number of constraints checked is usually equal to the number of components in the problem. This routine can check satisfaction of all of the constraints in the problem, which is equal to ne. However, the solver can't fix breakage of constraints above nc, because that nc is the range space by definition. Satisfaction of extra constraints would have had to occur in the problem specification.
The constraints should be broken up into 2 sections. If a constraint involves a formula matrix with positive and negative signs, and eaSet = 0.0, then you can't expect that the sum will be zero. There may be roundoff that inhibits this. However, if the formula matrix is all of one sign, then this requires that all species with nonzero entries in the formula matrix be identically zero. We put this into the logic below.
ibound | 1: Checks constraints up to the number of elements; 0: Checks constraints up to the number of components. |
Definition at line 26 of file vcs_elem.cpp.
References VCS_SOLVE::m_elemAbundances, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_elementActive, VCS_SOLVE::m_elType, VCS_SOLVE::m_formulaMatrix, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numSpeciesTot, Cantera::scale(), VCS_DELETE_MINORSPECIES_CUTOFF, VCS_ELEM_TYPE_ABSPOS, VCS_ELEM_TYPE_CHARGENEUTRALITY, and VCS_ELEM_TYPE_ELECTRONCHARGE.
Referenced by VCS_SOLVE::vcs_elcorr(), and VCS_SOLVE::vcs_inest_TP().
void vcs_elabPhase | ( | size_t | iphase, |
double *const | elemAbundPhase | ||
) |
Computes the elemental abundances vector for a single phase, elemAbundPhase[], and returns it through the argument list. The mole numbers of species are taken from the current value in m_molNumSpecies_old[].
Definition at line 95 of file vcs_elem.cpp.
References VCS_SOLVE::m_formulaMatrix, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_phaseID, VCS_SOLVE::m_speciesUnknownType, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
Referenced by VCS_SOLVE::vcs_report().
int vcs_elcorr | ( | double | aa[], |
double | x[] | ||
) |
This subroutine corrects for element abundances. At the end of the subroutine, the total moles in all phases are recalculated again, because we have changed the number of moles in this routine.
aa | temporary work vector, length ne*ne |
x | temporary work vector, length ne |
3 = The solution changed significantly, The element abundances are still bad and a component species got zeroed out.
Internal data to be worked on::
NOTES: This routine is turning out to be very problematic. There are lots of special cases and problems with zeroing out species.
Still need to check out when we do loops over nc vs. ne.
Definition at line 109 of file vcs_elem.cpp.
References VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_elemAbundances, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_elementName, VCS_SOLVE::m_elType, VCS_SOLVE::m_formulaMatrix, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numSpeciesRdc, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_speciesName, VCS_SOLVE::m_speciesStatus, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_SSPhase, Cantera::npos, plogf, VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabcheck(), VCS_ELEM_TYPE_ABSPOS, VCS_ELEM_TYPE_CHARGENEUTRALITY, VCS_ELEM_TYPE_ELECTRONCHARGE, VCS_SOLVE::vcs_reinsert_deleted(), VCS_SPECIES_ACTIVEBUTZERO, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_ZEROEDSS, VCS_SOLVE::vcs_tmoles(), and VCSnonideal::vcsUtil_mlequ().
Referenced by VCS_SOLVE::vcs_inest_TP().
int vcs_inest_TP | ( | ) |
Create an initial estimate of the solution to the thermodynamic equilibrium problem.
Definition at line 345 of file vcs_inest.cpp.
References VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_doEstimateEquil, VCS_SOLVE::m_feSpecies_new, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_tPhaseMoles_old, VCS_SOLVE::m_VCount, plogendl, plogf, clockWC::secondsWC(), VCS_COUNTERS::T_Calls_Inest, VCS_COUNTERS::T_Time_inest, VCS_DATA_PTR, VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_inest(), and VCS_SOLVE::vcs_Total_Gibbs().
Referenced by VCS_SOLVE::vcs_TP().
int vcs_setMolesLinProg | ( | ) |
Estimate the initial mole numbers by constrained linear programming.
This is done by running each reaction as far forward or backward as possible, subject to the constraint that all mole numbers remain non- negative. Reactions for which \( \Delta \mu^0 \) are positive are run in reverse, and ones for which it is negative are run in the forward direction. The end result is equivalent to solving the linear programming problem of minimizing the linear Gibbs function subject to the element and non-negativity constraints.
Definition at line 38 of file vcs_setMolesLinProg.cpp.
References plogf, VCS_DATA_PTR, VCS_SPECIES_INTERFACIALVOLTAGE, and VCS_SUCCESS.
Referenced by VCS_SOLVE::vcs_inest().
double vcs_Total_Gibbs | ( | double * | w, |
double * | fe, | ||
double * | tPhMoles | ||
) |
Calculate the total dimensionless Gibbs free energy.
Inert species are handled as if they had a standard free energy of zero. Note, for this algorithm this function should be monotonically decreasing.
Definition at line 18 of file vcs_Gibbs.cpp.
References vcs_VolPhase::m_gasPhase, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numSpeciesRdc, VCS_SOLVE::m_pressurePA, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_VolPhaseList, VCS_SOLVE::TPhInertMoles, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
Referenced by VCS_SOLVE::vcs_inest_TP(), and VCS_SOLVE::vcs_report().
double vcs_GibbsPhase | ( | size_t | iphase, |
const double *const | w, | ||
const double *const | fe | ||
) |
Calculate the total dimensionless Gibbs free energy of a single phase.
Inert species are handled as if they had a standard free energy of zero and if they obeyed ideal solution/gas theory.
iphase | ID of the phase |
w | Species mole number vector for all species |
fe | vector of partial molar free energies of all of the species |
Definition at line 43 of file vcs_Gibbs.cpp.
References vcs_VolPhase::m_gasPhase, VCS_SOLVE::m_numSpeciesRdc, VCS_SOLVE::m_phaseID, VCS_SOLVE::m_pressurePA, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_VolPhaseList, VCS_SOLVE::TPhInertMoles, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
Referenced by VCS_SOLVE::vcs_report().
int vcs_prob_update | ( | VCS_PROB * | pub | ) |
Transfer the results of the equilibrium calculation back to VCS_PROB.
The VCS_PROB structure is returned to the user.
pub | Pointer to VCS_PROB object that will get the results of the equilibrium calculation transfered to it. |
Definition at line 888 of file vcs_solve.cpp.
References vcs_VolPhase::electricPotential(), VCS_PROB::m_gibbsSpecies, VCS_PROB::m_Iterations, VCS_PROB::m_NumBasisOptimizations, VCS_PROB::mf, vcs_VolPhase::molefraction(), vcs_VolPhase::moleFractions(), VCS_PROB::NPhase, vcs_VolPhase::nSpecies(), vcs_VolPhase::phiVarIndex(), plogf, VCS_PROB::PresPA, vcs_VolPhase::setElectricPotential(), vcs_VolPhase::setMoleFractionsState(), vcs_VolPhase::setTotalMoles(), vcs_VolPhase::setTotalMolesInert(), VCS_PROB::SpeciesUnknownType, vcs_VolPhase::speciesUnknownType(), vcs_VolPhase::spGlobalIndexVCS(), VCS_PROB::T, vcs_VolPhase::totalMoles(), vcs_VolPhase::totalMolesInert(), VCS_DATA_PTR, VCSnonideal::vcs_doubleEqual(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_STATECALC_TMP, VCS_SUCCESS, VCS_PROB::Vol, VCS_PROB::VolPM, VCS_PROB::VPhaseList, and VCS_PROB::w.
int vcs_prob_specifyFully | ( | const VCS_PROB * | pub | ) |
Fully specify the problem to be solved using VCS_PROB.
Use the contents of the VCS_PROB to specify the contents of the private data, VCS_SOLVE.
pub | Pointer to VCS_PROB that will be used to initialize the current equilibrium problem |
Definition at line 381 of file vcs_solve.cpp.
References VCS_PROB::Charge, VCS_SPECIES_THERMO::duplMyselfAsVCS_SPECIES_THERMO(), VCS_PROB::ElActive, VCS_PROB::ElName, VCS_PROB::FormulaMatrix, VCS_PROB::gai, VCS_PROB::iest, VCS_PROB::m_elType, VCS_PROB::m_VCS_UnitsFormat, VCS_PROB::ne, VCS_PROB::NPhase, VCS_PROB::nspecies, vcs_VolPhase::nSpecies(), Cantera::OneAtm, vcs_VolPhase::p_activityConvention, VCS_PROB::PhaseID, vcs_VolPhase::PhaseName, plogf, VCS_PROB::PresPA, vcs_VolPhase::speciesProperty(), VCS_PROB::SpeciesThermo, VCS_PROB::SpeciesUnknownType, vcs_VolPhase::spGlobalIndexVCS(), VCS_PROB::SpName, VCS_PROB::T, VCS_PROB::tolmaj, VCS_PROB::tolmin, vcs_VolPhase::totalMolesInert(), VCS_DATA_PTR, VCS_PROB::vcs_debug_print_lvl, VCS_ELEM_TYPE_ABSPOS, VCS_ELEM_TYPE_CHARGENEUTRALITY, VCS_ELEM_TYPE_LATTICERATIO, VCS_SPECIES_MAJOR, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SUCCESS, VCS_PROB::Vol, VCS_PROB::VolPM, VCS_PROB::VPhaseList, VCS_PROB::w, and VCS_PROB::WtSpecies.
int vcs_prob_specify | ( | const VCS_PROB * | pub | ) |
Specify the problem to be solved using VCS_PROB, incrementally.
Use the contents of the VCS_PROB to specify the contents of the private data, VCS_SOLVE.
It's assumed we are solving the same problem.
pub | Pointer to VCS_PROB that will be used to initialize the current equilibrium problem |
Definition at line 763 of file vcs_solve.cpp.
References vcs_VolPhase::electricPotential(), VCS_PROB::gai, VCS_PROB::iest, vcs_VolPhase::m_eqnState, vcs_VolPhase::m_gasPhase, VCS_PROB::m_gibbsSpecies, vcs_VolPhase::m_singleSpecies, VCS_PROB::m_VCS_UnitsFormat, VCS_PROB::mf, vcs_VolPhase::nSpecies(), vcs_VolPhase::PhaseName, plogf, VCS_PROB::PresPA, vcs_VolPhase::setElectricPotential(), vcs_VolPhase::setExistence(), vcs_VolPhase::setTotalMolesInert(), VCS_PROB::T, VCS_PROB::tolmaj, VCS_PROB::tolmin, vcs_VolPhase::totalMolesInert(), VCS_SUCCESS, VCS_PROB::Vol, vcs_VolPhase::VP_ID_, VCS_PROB::VPhaseList, and VCS_PROB::w.
|
private |
Zero out the concentration of a species.
Make sure to conserveelements and keep track of the total moles in all phases.
kspec | Species index |
Definition at line 2174 of file vcs_solve_TP.cpp.
References plogendl, plogf, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
|
private |
Change a single species from active to inactive status.
Rearrange data when species is added or removed. The Lth species is moved to the back of the species vector. The back of the species vector is indicated by the value of MR, the current number of active species in the mechanism.
kspec | Species Index |
Definition at line 2199 of file vcs_solve_TP.cpp.
References vcs_VolPhase::exists(), plogendl, plogf, vcs_VolPhase::setMolesFromVCSCheck(), VCS_DATA_PTR, VCS_PHASE_EXIST_ALWAYS, VCS_SPECIES_DELETED, VCS_SPECIES_MAJOR, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, and VCS_STATECALC_OLD.
|
private |
This routine handles the bookkeeping involved with the deletion of multiphase phases from the problem.
When they are deleted, all of their species become active species, even though their mole numbers are set to zero. The routine does not make the decision to eliminate multiphases.
Note, species in phases with zero mole numbers are still considered active. Whether the phase pops back into existence or not is checked as part of the main iteration loop.
iph | Phase to be deleted |
Definition at line 2342 of file vcs_solve_TP.cpp.
References vcs_VolPhase::PhaseName, plogf, vcs_VolPhase::setTotalMoles(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, and VCS_SPECIES_ZEROEDMS.
|
private |
Change the concentration of a species by delta moles.
Make sure to conserve elements and keep track of the total kmoles in all phases.
kspec | The species index |
delta_ptr | pointer to the delta for the species. This may change during the calculation |
Definition at line 2110 of file vcs_solve_TP.cpp.
References plogendl, plogf, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, and VCS_STATECALC_OLD.
|
private |
Provide an estimate for the deleted species in phases that are not zeroed out.
Try to add back in all deleted species. An estimate of the kmol numbers are obtained and the species is added back into the equation system, into the old state vector.
This routine is called at the end of the calculation, just before returning to the user.
Definition at line 2635 of file vcs_solve_TP.cpp.
References vcs_VolPhase::m_singleSpecies, plogendl, plogf, vcs_VolPhase::sendToVCS_ActCoeff(), VCS_DATA_PTR, VCS_DELETE_MINORSPECIES_CUTOFF, VCS_STATECALC_NEW, and VCS_STATECALC_OLD.
|
private |
Recheck deleted species in multispecies phases.
We are checking the equation:
sum_u = sum_j_comp [ sigma_i_j * u_j ] = u_i_O + log((AC_i * W_i)/m_tPhaseMoles_old)
by first evaluating:
DG_i_O = u_i_O - sum_u.
Then, if TL is zero, the phase pops into existence if DG_i_O < 0. Also, if the phase exists, then we check to see if the species can have a mole number larger than VCS_DELETE_SPECIES_CUTOFF (default value = 1.0E-32).
Definition at line 2508 of file vcs_solve_TP.cpp.
References plogf, VCS_DATA_PTR, VCS_RELDELETE_SPECIES_CUTOFF, and VCS_STATECALC_NEW.
|
private |
Recheck deletion condition for multispecies phases.
We assume here that DG_i_0 has been calculated for deleted species correctly
m_feSpecies(I) = m_SSfeSpecies(I) + ln(ActCoeff[I]) - ln(Mnaught * m_units) + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase]; sum_u = sum_j_comp [ sigma_i_j * u_j ] = u_i_O + log((AC_i * W_i)/m_tPhaseMoles_old) DG_i_0 = m_feSpecies(I) - sum_m{ a_i_m DG_m }
by first evaluating:
DG_i_O = u_i_O - sum_u.
Then, the phase pops into existence iff
phaseDG = 1.0 - sum_i{exp(-DG_i_O)} < 0.0
This formula works for both single species phases and for multispecies phases. It's an overkill for single species phases.
iphase | Phase index number |
NOTE: this routine is currently not used in the code, and contains some basic changes that are incompatible.
assumptions:
Definition at line 2595 of file vcs_solve_TP.cpp.
References vcs_VolPhase::exists(), vcs_VolPhase::m_singleSpecies, vcs_VolPhase::nSpecies(), vcs_VolPhase::spGlobalIndexVCS(), VCS_PHASE_EXIST_NO, and VCS_PHASE_EXIST_ZEROEDPHASE.
|
private |
Minor species alternative calculation.
This is based upon the following approximation: The mole fraction changes due to these reactions don't affect the mole numbers of the component species. Therefore the following approximation is valid for a small component of an ideal phase:
0 = m_deltaGRxn_old(I) + log(molNum_new(I)/molNum_old(I))
m_deltaGRxn_old
contains the contribution from
m_feSpecies_old(I) = m_SSfeSpecies(I) + log(ActCoeff[i] * molNum_old(I) / m_tPhaseMoles_old(iph))
Thus,
molNum_new(I)= molNum_old(I) * EXP(-m_deltaGRxn_old(I))
Most of this section is mainly restricting the update to reasonable values. We restrict the update a factor of 1.0E10 up and 1.0E-10 down because we run into trouble with the addition operator due to roundoff if we go larger than ~1.0E15. Roundoff will then sometimes produce zero mole fractions.
Note: This routine was generalized to incorporate nonideal phases and phases on the molality basis
[in] | kspec | The current species and corresponding formation reaction number. |
[in] | irxn | The current species and corresponding formation reaction number. |
[out] | do_delete | BOOLEAN which if true on return, then we branch to the section that deletes a species from the current set of active species. |
Definition at line 1991 of file vcs_solve_TP.cpp.
References VCS_DELETE_MINORSPECIES_CUTOFF, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
|
private |
This routine optimizes the minimization of the total gibbs free energy by making sure the slope of the following functional stays negative:
The slope of the following functional is equivalent to the slope of the total Gibbs free energy of the system:
d_Gibbs/ds = sum_k( m_deltaGRxn * m_deltaMolNumSpecies[k] )
along the current direction m_deltaMolNumSpecies[], by choosing a value, al: (0<al<1) such that the a parabola approximation to Gibbs(al) fit to the end points al = 0 and al = 1 is minimized.
Only if there has been an inflection point (i.e., s1 < 0 and s2 > 0), does this code section kick in. It finds the point on the parabola where the slope is equal to zero.
Definition at line 2745 of file vcs_solve_TP.cpp.
References plogendl, plogf, VCS_DATA_PTR, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, and VCS_STATECALC_NEW.
|
private |
Switch rows and columns of a square matrix.
Switches the row and column of a matrix. So that after
J[k1][j] = J_old[k2][j] and J[j][k1] = J_old[j][k2] J[k2][j] = J_old[k1][j] and J[j][k2] = J_old[j][k1]
Jac | Double pointer to the Jacobian |
k1 | first row/column value to be switched |
k2 | second row/column value to be switched |
Definition at line 4357 of file vcs_solve_TP.cpp.
|
private |
Calculate the norm of a deltaGibbs free energy vector.
Positive DG for species which don't exist are ignored.
dg | Vector of local delta G's. |
Definition at line 4179 of file vcs_solve_TP.cpp.
References VCS_SPECIES_MAJOR, VCS_SPECIES_MINOR, and VCS_SPECIES_ZEROEDMS.
|
private |
Print out and check the elemental abundance vector.
Definition at line 4139 of file vcs_solve_TP.cpp.
References plogendl, plogf, and VCS_SPECIES_TYPE_INTERFACIALVOLTAGE.
|
private |
Estimate equilibrium compositions.
Algorithm covered in a section of Smith and Missen's Book.
Linear programming module is based on using dbolm.
aw | aw[i[ Mole fraction work space (ne in length) |
sa | sa[j] = Gram-Schmidt orthog work space (ne in length) |
sm | sm[i+j*ne] = QR matrix work space (ne*ne in length) |
ss | ss[j] = Gram-Schmidt orthog work space (ne in length) |
test | This is a small negative number. |
Definition at line 22 of file vcs_inest.cpp.
References VCS_SOLVE::m_debug_print_lvl, VCS_SOLVE::m_deltaGRxn_new, VCS_SOLVE::m_deltaMolNumPhase, VCS_SOLVE::m_deltaMolNumSpecies, VCS_SOLVE::m_deltaPhaseMoles, VCS_SOLVE::m_elemAbundancesGoal, VCS_SOLVE::m_elementActive, VCS_SOLVE::m_elementName, VCS_SOLVE::m_feSpecies_new, VCS_SOLVE::m_feSpecies_old, VCS_SOLVE::m_formulaMatrix, VCS_SOLVE::m_indexRxnToSpecies, VCS_SOLVE::m_molNumSpecies_new, VCS_SOLVE::m_molNumSpecies_old, VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numElemConstraints, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numRxnTot, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_phaseID, VCS_SOLVE::m_speciesName, VCS_SOLVE::m_speciesUnknownType, VCS_SOLVE::m_SSfeSpecies, VCS_SOLVE::m_SSPhase, VCS_SOLVE::m_stoichCoeffRxnMatrix, VCS_SOLVE::m_TmpPhase, VCS_SOLVE::m_TmpPhase2, VCS_SOLVE::m_tPhaseMoles_new, VCS_SOLVE::m_tPhaseMoles_old, VCS_SOLVE::m_VolPhaseList, plogendl, plogf, VCS_SOLVE::TPhInertMoles, VCS_SOLVE::vcs_basopt(), VCS_DATA_PTR, VCS_SOLVE::vcs_deltag(), VCS_SOLVE::vcs_dfe(), VCS_SOLVE::vcs_setMolesLinProg(), VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_SPECIES_TYPE_MOLNUM, VCS_STATECALC_NEW, VCS_STATECALC_OLD, and VCS_SOLVE::vcs_tmoles().
Referenced by VCS_SOLVE::vcs_inest_TP().
|
private |
Calculate the status of single species phases.
Definition at line 20 of file vcs_prep.cpp.
References VCS_SOLVE::m_numPhases, VCS_SOLVE::m_numSpeciesTot, VCS_SOLVE::m_phaseID, vcs_VolPhase::m_singleSpecies, VCS_SOLVE::m_SSPhase, VCS_SOLVE::m_VolPhaseList, vcs_VolPhase::setExistence(), and VCS_SOLVE::TPhInertMoles.
Referenced by VCS_SOLVE::vcs_prep_oneTime().
|
private |
This function recalculates the deltaG for reaction, irxn.
This function recalculates the deltaG for reaction irxn, given the mole numbers in molNum. It uses the temporary space mu_i, to hold the recalculated chemical potentials. It only recalculates the chemical potentials for species in phases which participate in the irxn reaction. This function is used by the vcs_line_search algorithm() and should not be used widely due to the unknown state it leaves the system.
[in] | irxn | Reaction number |
[in] | molNum | Current mole numbers of species to be used as input to the calculation (units = kmol), (length = totalNuMSpecies) |
[out] | ac | Activity coefficients (length = totalNumSpecies) Note this is only partially formed. Only species in phases that participate in the reaction will be updated |
[out] | mu_i | dimensionless chemical potentials (length totalNumSpecies) Note this is only partially formed. Only species in phases that participate in the reaction will be updated |
Definition at line 675 of file vcs_rxnadj.cpp.
References VCS_SOLVE::m_numComponents, VCS_SOLVE::m_numPhases, VCS_SOLVE::m_phaseParticipation, VCS_SOLVE::m_stoichCoeffRxnMatrix, and VCS_SOLVE::vcs_chemPotPhase().
|
private |
Delete memory that isn't just resizable STL containers.
This gets called by the destructor or by InitSizes().
Definition at line 223 of file vcs_solve.cpp.
|
private |
Initialize the internal counters.
Initialize the internal counters containing the subroutine call values and times spent in the subroutines.
ifunc = 0 Initialize only those counters appropriate for the top of vcs_solve_TP(). = 1 Initialize all counters.
Definition at line 972 of file vcs_solve.cpp.
|
private |
Create a report on the plog file containing timing and its information.
timing_print_lvl | If 0, just report the iteration count. If larger than zero, report the timing information |
Definition at line 396 of file vcs_report.cpp.
References VCS_SOLVE::m_VCount, plogf, VCS_COUNTERS::T_Basis_Opts, VCS_COUNTERS::T_Calls_Inest, VCS_COUNTERS::T_Calls_vcs_TP, VCS_COUNTERS::T_Its, VCS_COUNTERS::T_Time_basopt, VCS_COUNTERS::T_Time_inest, VCS_COUNTERS::T_Time_vcs, and VCS_COUNTERS::T_Time_vcs_TP.
|
private |
Update all underlying vcs_VolPhase objects.
Update the mole numbers and the phase voltages of all phases in the vcs problem
stateCalc | Location of the update (either VCS_STATECALC_NEW or VCS_STATECALC_OLD). |
Definition at line 5083 of file vcs_solve_TP.cpp.
References vcs_VolPhase::updateFromVCS_MoleNumbers().
int vcs_rank | ( | const double * | awtmp, |
size_t | numSpecies, | ||
const double * | matrix, | ||
size_t | numElemConstraints, | ||
std::vector< size_t > & | compRes, | ||
std::vector< size_t > & | elemComp, | ||
int *const | usedZeroedSpecies | ||
) | const |
Calculate the rank of a matrix and return the rows and columns that will generate an independent basis for that rank.
Definition at line 46 of file vcs_rank.cpp.
References plogendl, plogf, and VCS_DATA_PTR.
size_t NSPECIES0 |
value of the number of species used to malloc data structures
Definition at line 1491 of file vcs_solve.h.
size_t NPHASE0 |
value of the number of phases used to malloc data structures
Definition at line 1494 of file vcs_solve.h.
size_t m_numSpeciesTot |
Total number of species in the problems.
Definition at line 1497 of file vcs_solve.h.
Referenced by vcs_MultiPhaseEquil::getStoichVector(), VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elabPhase(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_fePrep_TP(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_rearrange(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_RxnStepSizes(), VCS_SOLVE::vcs_SSPhase(), and VCS_SOLVE::vcs_switch_elem_pos().
size_t m_numElemConstraints |
Number of element constraints in the problem.
This is typically equal to the number of elements in the problem
Definition at line 1503 of file vcs_solve.h.
Referenced by vcs_MultiPhaseEquil::numElemConstraints(), VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elabPhase(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_switch_elem_pos().
size_t m_numComponents |
Number of components calculated for the problem.
Definition at line 1506 of file vcs_solve.h.
Referenced by VCS_SOLVE::deltaG_Recalc_Rxn(), vcs_MultiPhaseEquil::numComponents(), VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_Hessian_actCoeff_diag(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
size_t m_numRxnTot |
Total number of non-component species in the problem.
Definition at line 1509 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_prep_oneTime(), and VCS_SOLVE::vcs_report().
size_t m_numSpeciesRdc |
Current number of species in the problems.
Species can be deleted if they aren't stable under the current conditions
Definition at line 1516 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_GibbsPhase(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_Total_Gibbs().
size_t m_numRxnRdc |
Current number of non-component species in the problem.
Species can be deleted if they aren't stable under the current conditions
Definition at line 1523 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
size_t m_numRxnMinorZeroed |
Number of active species which are currently either treated as minor species.
Definition at line 1527 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_rxn_adj_cg().
size_t m_numPhases |
Number of Phases in the problem.
Definition at line 1530 of file vcs_solve.h.
Referenced by VCS_SOLVE::deltaG_Recalc_Rxn(), VCS_SOLVE::vcs_CalcLnActCoeffJac(), VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), VCS_SOLVE::vcs_RxnStepSizes(), VCS_SOLVE::vcs_SSPhase(), VCS_SOLVE::vcs_switch_elem_pos(), and VCS_SOLVE::vcs_Total_Gibbs().
DoubleStarStar m_formulaMatrix |
Formula matrix for the problem.
FormulaMatrix[j][kspec] = Number of elements, j, in the kspec species
Both element and species indices are swapped.
Definition at line 1538 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elabPhase(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_inest(), and VCS_SOLVE::vcs_switch_elem_pos().
DoubleStarStar m_stoichCoeffRxnMatrix |
Stoichiometric coefficient matrix for the reaction mechanism expressed in Reduced Canonical Form.
This is the stoichiometric coefficient matrix for the reaction which forms species kspec from the component species. A stoichiometric coefficient of one is assumed for the species kspec in this mechanism.
NOTE: kspec = irxn + m_numComponents
m_stoichCoeffRxnMatrix[irxn][j] : j refers to the component number, and irxn refers to the irxn_th non-component species. The stoichiometric coefficients multilplied by the Formula coefficients of the component species add up to the negative value of the number of elements in the species kspec.
length = [nspecies0][nelements0]
Definition at line 1556 of file vcs_solve.h.
Referenced by VCS_SOLVE::deltaG_Recalc_Rxn(), vcs_MultiPhaseEquil::getStoichVector(), VCS_SOLVE::vcs_Hessian_actCoeff_diag(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<double> m_scSize |
Absolute size of the stoichiometric coefficients.
scSize[irxn] = abs(Size) of the stoichiometric coefficients. These are used to determine whether a given species should be handled by the alt_min treatment or should be handled as a major species.
Definition at line 1566 of file vcs_solve.h.
std::vector<double> m_spSize |
total size of the species
This is used as a multiplier to the mole number in figuring out which species should be components.
Definition at line 1573 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_prep_oneTime().
std::vector<double> m_SSfeSpecies |
Standard state chemical potentials for species K at the current temperature and pressure.
The first NC entries are for components. The following NR entries are for the current non-component species in the mechanism.
Definition at line 1581 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_fePrep_TP(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().
std::vector<double> m_feSpecies_old |
Free energy vector from the start of the current iteration.
The free energies are saved at the start of the current iteration. Length = number of species
Definition at line 1588 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_fePrep_TP(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().
std::vector<double> m_feSpecies_new |
Dimensionless new free energy for all the species in the mechanism at the new tentatite T, P, and mole numbers.
The first NC entries are for components. The following NR entries are for the current non-component species in the mechanism. Length = number of species
Definition at line 1597 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_fePrep_TP(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_inest_TP(), and VCS_SOLVE::vcs_prep().
int m_doEstimateEquil |
Setting for whether to do an initial estimate.
Initial estimate:
Definition at line 1609 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_prep_oneTime(), and VCS_SOLVE::vcs_TP().
std::vector<double> m_molNumSpecies_old |
Total moles of the species.
Total number of moles of the kth species. Length = Total number of species = m
Definition at line 1616 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elabPhase(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<int> m_speciesUnknownType |
Specifies the species unknown type.
There are two types. One is the straightforward species, with the mole number w[k], as the unknown. The second is the an interfacial voltage where w[k] refers to the interfacial voltage in volts.
These species types correspond to metallic electrons corresponding to electrodes. The voltage and other interfacial conditions sets up an interfacial current, which is set to zero in this initial treatment. Later we may have non-zero interfacial currents.
Definition at line 1629 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabPhase(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_GibbsPhase(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_RxnStepSizes(), and VCS_SOLVE::vcs_Total_Gibbs().
DoubleStarStar m_deltaMolNumPhase |
Change in the number of moles of phase, iphase, due to the noncomponent formation reaction, irxn, for species, k:
m_deltaMolNumPhase[irxn][iphase] = k = nc + irxn
Definition at line 1636 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_prep(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
IntStarStar m_phaseParticipation |
This is 1 if the phase, iphase, participates in the formation reaction irxn, and zero otherwise.
PhaseParticipation[irxn][iphase]
Definition at line 1640 of file vcs_solve.h.
Referenced by VCS_SOLVE::deltaG_Recalc_Rxn(), and VCS_SOLVE::vcs_prep().
std::vector<double> m_phasePhi |
electric potential of the iph phase
Definition at line 1643 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_report().
std::vector<double> m_molNumSpecies_new |
Tentative value of the mole number vector.
It's also used to store the mole fraction vector.
Definition at line 1648 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_prep(), and VCS_SOLVE::vcs_report().
std::vector<double> m_deltaGRxn_new |
Delta G(irxn) for the noncomponent species in the mechanism.
Computed by the subroutine deltaG. m_deltaGRxn is the free energy change for the reaction which forms species K from the component species. This vector has length equal to the number of noncomponent species in the mechanism. It starts with the first current noncomponent species in the mechanism.
Definition at line 1658 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<double> m_deltaGRxn_old |
Last deltag[irxn] from the previous step.
Definition at line 1661 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_nondim_TP(), and VCS_SOLVE::vcs_redim_TP().
std::vector<double> m_deltaGRxn_Deficient |
Last deltag[irxn] from the previous step with additions for possible births of zeroed phases.
Definition at line 1665 of file vcs_solve.h.
std::vector<double> m_deltaGRxn_tmp |
Temporary vector of Rxn DeltaG's.
This is used from time to time, for printing purposes
Definition at line 1671 of file vcs_solve.h.
std::vector<double> m_deltaMolNumSpecies |
Reaction Adjustments for each species during the current step.
delta Moles for each species during the current step. Length = number of species
Definition at line 1678 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<double> m_elemAbundances |
Element abundances vector.
Vector of moles of each element actually in the solution vector. Except for certain parts of the algorithm, this is a constant. Note other constraint conditions are added to this vector. This is input from the input file and is considered a constant from thereon. units = kmoles
Definition at line 1690 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elab(), VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_switch_elem_pos().
std::vector<double> m_elemAbundancesGoal |
Element abundances vector Goals.
Vector of moles of each element that are the goals of the simulation. This is a constant in the problem. Note other constraint conditions are added to this vector. This is input from the input file and is considered a constant from thereon. units = kmoles
Definition at line 1699 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_switch_elem_pos().
double m_totalMolNum |
Total number of kmoles in all phases.
This number includes the inerts. -> Don't use this except for scaling purposes
Definition at line 1706 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_RxnStepSizes().
std::vector<double> m_tPhaseMoles_old |
Total kmols of species in each phase.
This contains the total number of moles of species in each phase
Length = number of phases
Definition at line 1714 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_Hessian_actCoeff_diag(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<double> m_tPhaseMoles_new |
total kmols of species in each phase in the tentative soln vector
This contains the total number of moles of species in each phase in the tentative solution vector
Length = number of phases
Definition at line 1723 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest(), and VCS_SOLVE::vcs_prep().
|
mutable |
Temporary vector of length NPhase.
Definition at line 1726 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest().
|
mutable |
Temporary vector of length NPhase.
Definition at line 1729 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest().
std::vector<double> m_deltaPhaseMoles |
Change in the total moles in each phase.
Length number of phases.
Definition at line 1735 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest(), and VCS_SOLVE::vcs_prep().
double m_temperature |
Temperature (Kelvin)
Definition at line 1738 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_TP().
double m_pressurePA |
Pressure (units are determined by m_VCS_UnitsFormat.
Values | units |
---|---|
-1: | atm |
0: | atm |
1: | atm |
2: | atm |
3: | Pa |
Units being changed to Pa
Definition at line 1752 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_GibbsPhase(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_Total_Gibbs(), and VCS_SOLVE::vcs_TP().
std::vector<double> TPhInertMoles |
Total kmoles of inert to add to each phase.
TPhInertMoles[iph] = Total kmoles of inert to add to each phase length = number of phases
Definition at line 1759 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_GibbsPhase(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_SSPhase(), and VCS_SOLVE::vcs_Total_Gibbs().
double m_tolmaj |
Tolerance requirement for major species.
Definition at line 1762 of file vcs_solve.h.
double m_tolmin |
Tolerance requirements for minor species.
Definition at line 1765 of file vcs_solve.h.
double m_tolmaj2 |
Below this, major species aren't refined any more.
Definition at line 1768 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
double m_tolmin2 |
Below this, minor species aren't refined any more.
Definition at line 1771 of file vcs_solve.h.
std::vector<size_t> m_speciesMapIndex |
Index vector that keeps track of the species vector rearrangement.
At the end of each run, the species vector and associated data gets put back in the original order.
Example
k = m_speciesMapIndex[kspec] kspec = current order in the vcs_solve object k = original order in the vcs_prob object and in the MultiPhase object
Definition at line 1785 of file vcs_solve.h.
Referenced by vcs_MultiPhaseEquil::component(), vcs_MultiPhaseEquil::getStoichVector(), and VCS_SOLVE::vcs_rearrange().
std::vector<size_t> m_speciesLocalPhaseIndex |
Index that keeps track of the index of the species within the local phase.
This returns the local index of the species within the phase. Its argument is the global species index within the VCS problem.
k = m_speciesLocalPhaseIndex[kspec]
k varies between 0 and the nSpecies in the phase
Length = number of species
Definition at line 1799 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_prep_oneTime().
std::vector<size_t> m_elementMapIndex |
Index vector that keeps track of the rearrangement of the elements.
At the end of each run, the element vector and associated data gets put back in the original order.
Example
e = m_elementMapIndex[eNum] eNum = current order in the vcs_solve object e = original order in the vcs_prob object and in the MultiPhase object
Definition at line 1812 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_switch_elem_pos().
std::vector<size_t> m_indexRxnToSpecies |
Mapping between the species index for noncomponent species and the full species index.
ir[irxn] = Mapping between the reaction index for noncomponent formation reaction of a species and the full species index.
Initially set to a value of K = NC + I This vector has length equal to number of noncomponent species in the mechanism. It starts with the first current noncomponent species in the mechanism. kspec = ir[irxn]
Definition at line 1825 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_Hessian_actCoeff_diag(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<int> m_speciesStatus |
Major -Minor status vector for the species in the problem.
The index for this is species. The reaction that this is referring to is kspec = irxn + m_numComponents
. For possible values and their meanings, see vcs_evaluate_speciesType().
Definition at line 1833 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<size_t> m_phaseID |
Mapping from the species number to the phase number.
Definition at line 1836 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elabPhase(), VCS_SOLVE::vcs_GibbsPhase(), VCS_SOLVE::vcs_Hessian_actCoeff_diag(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), VCS_SOLVE::vcs_RxnStepSizes(), and VCS_SOLVE::vcs_SSPhase().
std::vector<char> m_SSPhase |
Boolean indicating whether a species belongs to a single-species phase.
Definition at line 1840 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_fePrep_TP(), VCS_SOLVE::vcs_Hessian_actCoeff_diag(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_rxn_adj_cg(), VCS_SOLVE::vcs_RxnStepSizes(), and VCS_SOLVE::vcs_SSPhase().
std::vector<std::string> m_speciesName |
Species string name for the kth species.
Definition at line 1843 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), and VCS_SOLVE::vcs_RxnStepSizes().
std::vector<std::string> m_elementName |
Vector of strings containing the element names.
ElName[j] = String containing element names
Definition at line 1849 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_switch_elem_pos().
std::vector<int> m_elType |
Type of the element constraint.
Definition at line 1864 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_switch_elem_pos().
std::vector<int> m_elementActive |
Specifies whether an element constraint is active.
The default is true. Length = nelements
Definition at line 1870 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elabcheck(), VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_switch_elem_pos().
std::vector<vcs_VolPhase*> m_VolPhaseList |
Array of Phase Structures. Length = number of phases.
Definition at line 1873 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_CalcLnActCoeffJac(), VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_GibbsPhase(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_redim_TP(), VCS_SOLVE::vcs_report(), VCS_SOLVE::vcs_rxn_adj_cg(), VCS_SOLVE::vcs_RxnStepSizes(), VCS_SOLVE::vcs_SSPhase(), VCS_SOLVE::vcs_switch_elem_pos(), and VCS_SOLVE::vcs_Total_Gibbs().
std::string m_title |
String containing the title of the run.
Definition at line 1876 of file vcs_solve.h.
char m_unitsState |
This specifies the current state of units for the Gibbs free energy properties in the program.
The default is to have this unitless
Definition at line 1883 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().
double m_totalMoleScale |
Multiplier for the mole numbers within the nondimensionless formulation.
All numbers within the main routine are on an absolute basis. This presents some problems wrt very large and very small mole numbers. We get around this by using a multiplier coming into and coming out of the equilibrium routines
Definition at line 1892 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().
std::vector<int> m_actConventionSpecies |
specifies the activity convention of the phase containing the species
length = number of species
Definition at line 1901 of file vcs_solve.h.
std::vector<int> m_phaseActConvention |
specifies the activity convention of the phase.
length = number of phases
Definition at line 1910 of file vcs_solve.h.
std::vector<double> m_lnMnaughtSpecies |
specifies the ln(Mnaught) used to calculate the chemical potentials
For molar based activity conventions this will be equal to 0.0. length = number of species.
Definition at line 1917 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_report().
std::vector<double> m_actCoeffSpecies_new |
Molar-based Activity Coefficients for Species.
Length = number of species
Definition at line 1921 of file vcs_solve.h.
std::vector<double> m_actCoeffSpecies_old |
Molar-based Activity Coefficients for Species based on old mole numbers.
These activity coefficients are based on the m_molNumSpecies_old values Molar based activity coeffients. Length = number of species
Definition at line 1928 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_report().
DoubleStarStar m_np_dLnActCoeffdMolNum |
Change in the log of the activity coefficient with respect to the mole number multiplied by the phase mole number.
length = [nspecies][nspecies]
This is a temporary array that gets regenerated every time it's needed. It is not swapped wrt species.
Definition at line 1938 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_CalcLnActCoeffJac(), and VCS_SOLVE::vcs_Hessian_actCoeff_diag().
std::vector<double> m_wtSpecies |
Molecular weight of each species.
units = kg/kmol. length = number of species.
note: this is a candidate for removal. I don't think we use it.
Definition at line 1946 of file vcs_solve.h.
std::vector<double> m_chargeSpecies |
Charge of each species. Length = number of species.
Definition at line 1949 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_report().
std::vector<VCS_SPECIES_THERMO*> m_speciesThermoList |
Vector of pointers to thermostructures which identify the model and parameters for evaluating the thermodynamic functions for that particular species.
SpeciesThermo[k] pointer to the thermo information for the kth species
Definition at line 1959 of file vcs_solve.h.
int m_useActCoeffJac |
Choice of Hessians.
If this is true, then we will use a better approximation to the Hessian based on Jacobian of the ln(ActCoeff) with respect to mole numbers
Definition at line 1967 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_RxnStepSizes().
double m_totalVol |
Total volume of all phases. Units are m^3.
Definition at line 1970 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_report().
std::vector<double> m_PMVolumeSpecies |
Partial molar volumes of the species.
units = mks (m^3/kmol) -determined by m_VCS_UnitsFormat Length = number of species
Definition at line 1977 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_report().
double m_Faraday_dim |
dimensionless value of Faraday's constant, F / RT (1/volt)
Definition at line 1980 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().
VCS_COUNTERS* m_VCount |
Timing and iteration counters for the vcs object.
Definition at line 1983 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_report(), and VCS_SOLVE::vcs_TCounters_report().
int m_debug_print_lvl |
Debug printing lvl.
Levels correspond to the following guidlines
Levels of printing above 4 are only accessible when DEBUG_MODE is turned on
Definition at line 1999 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_elcorr(), VCS_SOLVE::vcs_elem_rearrange(), VCS_SOLVE::vcs_inest(), VCS_SOLVE::vcs_inest_TP(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_prep_oneTime(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_RxnStepSizes().
int m_timing_print_lvl |
printing level of timing information
Definition at line 2006 of file vcs_solve.h.
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability(), vcs_MultiPhaseEquil::equilibrate_TP(), and VCS_SOLVE::vcs_report().
int m_VCS_UnitsFormat |
Units for the chemical potential data.
Value | chemical potential units | pressure units |
---|---|---|
-1 | kcal/mol | Pa |
0 | MU/RT | Pa |
1 | kJ/mol | Pa |
2 | Kelvin | Pa |
3 | J / kmol | Pa |
Definition at line 2018 of file vcs_solve.h.
Referenced by VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_nondim_TP(), VCS_SOLVE::vcs_redim_TP(), and VCS_SOLVE::vcs_report().