Cantera  2.4.0
vcs_VolPhase.h
Go to the documentation of this file.
1 /**
2  * @file vcs_VolPhase.h
3  * Header for the object representing each phase within vcs
4  */
5 
6 // This file is part of Cantera. See License.txt in the top-level directory or
7 // at http://www.cantera.org/license.txt for license and copyright information.
8 
9 #ifndef VCS_VOLPHASE_H
10 #define VCS_VOLPHASE_H
11 
13 #include "cantera/base/Array.h"
14 
15 namespace Cantera
16 {
17 
18 class ThermoPhase;
19 
20 //! Models for the standard state volume of each species
21 #define VCS_SSVOL_IDEALGAS 0
22 #define VCS_SSVOL_CONSTANT 1
23 
24 /*
25  * DEFINITIONS FOR THE vcs_VolPhase structure
26  *
27  * Equation of State Types
28  * - Permissible values for the EqnState variable in CPC_PHASE structure
29  */
30 #define VCS_EOS_CONSTANT 0
31 #define VCS_EOS_IDEAL_GAS 1
32 #define VCS_EOS_STOICH_SUB 5
33 #define VCS_EOS_IDEAL_SOLN 22
34 #define VCS_EOS_DEBEYE_HUCKEL 23
35 #define VCS_EOS_REDLICH_KWONG 24
36 #define VCS_EOS_REGULAR_SOLN 25
37 #define VCS_EOS_UNK_CANTERA -1
38 
39 struct VCS_SPECIES;
40 class vcs_SpeciesProperties;
41 class VCS_SOLVE;
42 
43 //! Phase information and Phase calculations for vcs.
44 /*!
45  * Each phase in a vcs calculation has a vcs_VolPhase object associated with it.
46  * This object helps to coordinate property evaluations for species within the
47  * phase. Usually these evaluations must be carried out on a per phase basis.
48  * However, vcs frequently needs per species quantities. Therefore, we need an
49  * interface layer between vcs and Cantera's ThermoPhase.
50  *
51  * The species stay in the same ordering within this structure. The vcs
52  * algorithm will change the ordering of species in the global species list.
53  * However, the indexing of species in this list stays the same. This structure
54  * contains structures that point to the species belonging to this phase in the
55  * global vcs species list.
56  *
57  * This object is considered not to own the underlying Cantera ThermoPhase
58  * object for the phase.
59  *
60  * This object contains an idea of the temperature and pressure. It checks to
61  * see if if the temperature and pressure has changed before calling underlying
62  * property evaluation routines.
63  *
64  * The object contains values for the electric potential of a phase. It
65  * coordinates the evaluation of properties wrt when the electric potential of a
66  * phase has changed.
67  *
68  * The object knows about the mole fractions of the phase. It controls the
69  * values of mole fractions, and coordinates the property evaluation wrt to
70  * changes in the mole fractions. It also will keep track of the likely values
71  * of mole fractions in multicomponent phases even when the phase doesn't
72  * actually exist within the thermo program.
73  *
74  * The object knows about the total moles of a phase. It checks to see if the
75  * phase currently exists or not, and modifies its behavior accordingly.
76  *
77  * Activity coefficients and volume calculations are lagged. They are only
78  * called when they are needed (and when the state has changed so that they need
79  * to be recalculated).
80  */
82 {
83 public:
84  vcs_VolPhase(VCS_SOLVE* owningSolverObject = 0);
85 
86  vcs_VolPhase(const vcs_VolPhase& b) = delete;
87  vcs_VolPhase& operator=(const vcs_VolPhase& b) = delete;
88  ~vcs_VolPhase();
89 
90  //! The resize() function fills in all of the initial information if it
91  //! is not given in the constructor.
92  /*!
93  * @param phaseNum index of the phase in the vcs problem
94  * @param numSpecies Number of species in the phase
95  * @param phaseName String name for the phase
96  * @param molesInert kmoles of inert in the phase (defaults to zero)
97  */
98  void resize(const size_t phaseNum, const size_t numSpecies,
99  const size_t numElem, const char* const phaseName,
100  const double molesInert = 0.0);
101 
102  void elemResize(const size_t numElemConstraints);
103 
104  //! Set the moles and/or mole fractions within the phase
105  /*!
106  * @param molNum total moles in the phase
107  * @param moleFracVec Vector of input mole fractions
108  * @param vcsStateStatus Status flag for this update
109  */
110  void setMoleFractionsState(const double molNum, const double* const moleFracVec,
111  const int vcsStateStatus);
112 
113  //! Set the moles within the phase
114  /*!
115  * This function takes as input the mole numbers in vcs format, and then
116  * updates this object with their values. This is essentially a gather
117  * routine.
118  *
119  * @param molesSpeciesVCS Array of mole numbers. Note, the indices for
120  * species in this array may not be contiguous. IndSpecies[] is needed
121  * to gather the species into the local contiguous vector format.
122  */
123  void setMolesFromVCS(const int stateCalc,
124  const double* molesSpeciesVCS = 0);
125 
126  //! Set the moles within the phase
127  /*!
128  * This function takes as input the mole numbers in vcs format, and then
129  * updates this object with their values. This is essentially a gather
130  * routine.
131  *
132  * Additionally it checks to see that the total moles value in
133  * TPhMoles[iplace] is equal to the internally computed value. If this isn't
134  * the case, an error exit is carried out.
135  *
136  * @param vcsStateStatus State calc value either `VCS_STATECALC_OLD` or
137  * `VCS_STATECALC_NEW`. With any other value nothing is done.
138  * @param molesSpeciesVCS array of mole numbers. Note, the indices for
139  * species in this array may not be contiguous. IndSpecies[] is needed
140  * to gather the species into the local contiguous vector format.
141  * @param TPhMoles VCS's array containing the number of moles in each phase.
142  */
143  void setMolesFromVCSCheck(const int vcsStateStatus,
144  const double* molesSpeciesVCS,
145  const double* const TPhMoles);
146 
147  //! Update the moles within the phase, if necessary
148  /*!
149  * This function takes as input the stateCalc value, which determines where
150  * within VCS_SOLVE to fetch the mole numbers. It then updates this object
151  * with their values. This is essentially a gather routine.
152  *
153  * @param stateCalc State calc value either VCS_STATECALC_OLD or
154  * VCS_STATECALC_NEW. With any other value nothing is done.
155  */
156  void updateFromVCS_MoleNumbers(const int stateCalc);
157 
158  //! Fill in an activity coefficients vector within a VCS_SOLVE object
159  /*!
160  * This routine will calculate the activity coefficients for the current
161  * phase, and fill in the corresponding entries in the VCS activity
162  * coefficients vector.
163  *
164  * @param AC vector of activity coefficients for all of the species in all
165  * of the phases in a VCS problem. Only the entries for the current
166  * phase are filled in.
167  */
168  void sendToVCS_ActCoeff(const int stateCalc, double* const AC);
169 
170  //! set the electric potential of the phase
171  /*!
172  * @param phi electric potential (volts)
173  */
174  void setElectricPotential(const double phi);
175 
176  //! Returns the electric field of the phase
177  /*!
178  * Units are potential
179  */
180  double electricPotential() const;
181 
182  //! Gibbs free energy calculation for standard state of one species
183  /*!
184  * Calculate the Gibbs free energies for the standard state of the kth
185  * species. The results are held internally within the object.
186  *
187  * @param kspec Species number (within the phase)
188  * @returns the Gibbs free energy for the standard state of the kth species.
189  */
190  double GStar_calc_one(size_t kspec) const;
191 
192  //! Gibbs free energy calculation at a temperature for the reference state
193  //! of a species, return a value for one species
194  /*!
195  * @param kspec species index
196  * @returns value of the Gibbs free energy
197  */
198  double G0_calc_one(size_t kspec) const;
199 
200  //! Molar volume calculation for standard state of one species
201  /*!
202  * Calculate the molar volume for the standard states. The results are held
203  * internally within the object.
204  *
205  * @param kspec Species number (within the phase)
206  * @return molar volume of the kspec species's standard state (m**3/kmol)
207  */
208  double VolStar_calc_one(size_t kspec) const;
209 
210  //! Fill in the partial molar volume vector for VCS
211  /*!
212  * This routine will calculate the partial molar volumes for the current
213  * phase (if needed), and fill in the corresponding entries in the VCS
214  * partial molar volumes vector.
215  *
216  * @param VolPM vector of partial molar volumes for all of the species in
217  * all of the phases in a VCS problem. Only the entries for the current
218  * phase are filled in.
219  */
220  double sendToVCS_VolPM(double* const VolPM) const;
221 
222  //! Fill in the partial molar volume vector for VCS
223  /*!
224  * This routine will calculate the partial molar volumes for the
225  * current phase (if needed), and fill in the corresponding entries in the
226  * VCS partial molar volumes vector.
227  *
228  * @param VolPM vector of partial molar volumes for all of the species in
229  * all of the phases in a VCS problem. Only the entries for the current
230  * phase are filled in.
231  *
232  * @todo This function's documentation is incorrect.
233  */
234  void sendToVCS_GStar(double* const gstar) const;
235 
236  //! Sets the temperature and pressure in this object and underlying
237  //! ThermoPhase objects
238  /*!
239  * @param temperature_Kelvin (Kelvin)
240  * @param pressure_PA Pressure (MKS units - Pascal)
241  */
242  void setState_TP(const double temperature_Kelvin, const double pressure_PA);
243 
244  //! Sets the temperature in this object and underlying ThermoPhase objects
245  /*!
246  * @param temperature_Kelvin (Kelvin)
247  */
248  void setState_T(const double temperature_Kelvin);
249 
250  // Downloads the ln ActCoeff Jacobian into the VCS version of the
251  // ln ActCoeff Jacobian.
252  /*
253  * This is essentially a scatter operation.
254  *
255  * @param LnAcJac_VCS Jacobian parameter
256  * The Jacobians are actually d( lnActCoeff) / d (MolNumber);
257  * dLnActCoeffdMolNumber(k,j)
258  *
259  * j = id of the species mole number
260  * k = id of the species activity coefficient
261  */
262  void sendToVCS_LnActCoeffJac(Array2D& LnACJac_VCS);
263 
264  //! Set the pointer for Cantera's ThermoPhase parameter
265  /*!
266  * When we first initialize the ThermoPhase object, we read the state of the
267  * ThermoPhase into vcs_VolPhase object.
268  *
269  * @param tp_ptr Pointer to the ThermoPhase object corresponding
270  * to this phase.
271  */
272  void setPtrThermoPhase(ThermoPhase* tp_ptr);
273 
274  //! Return a const ThermoPhase pointer corresponding to this phase
275  /*!
276  * @return pointer to the ThermoPhase.
277  */
278  const ThermoPhase* ptrThermoPhase() const;
279 
280  //! Return the total moles in the phase
281  double totalMoles() const;
282 
283  //! Returns the mole fraction of the kspec species
284  /*!
285  * @param kspec Index of the species in the phase
286  *
287  * @return Value of the mole fraction
288  */
289  double molefraction(size_t kspec) const;
290 
291  //! Sets the total moles in the phase
292  /*!
293  * We don't have to flag the internal state as changing here because we have
294  * just changed the total moles.
295  *
296  * @param totalMols Total moles in the phase (kmol)
297  */
298  void setTotalMoles(const double totalMols);
299 
300  //! Sets the mole flag within the object to out of date
301  /*!
302  * This will trigger the object to go get the current mole numbers when it
303  * needs it.
304  */
305  void setMolesOutOfDate(int stateCalc = -1);
306 
307  //! Sets the mole flag within the object to be current
308  void setMolesCurrent(int vcsStateStatus);
309 
310 private:
311  //! Set the mole fractions from a conventional mole fraction vector
312  /*!
313  * @param xmol Value of the mole fractions for the species in the phase.
314  * These are contiguous.
315  */
316  void setMoleFractions(const double* const xmol);
317 
318 public:
319  //! Return a const reference to the mole fractions stored in the object.
320  const vector_fp & moleFractions() const;
321 
322  double moleFraction(size_t klocal) const;
323 
324  //! Sets the creationMoleNum's within the phase object
325  /*!
326  * @param F_k Pointer to a vector of n_k's
327  */
328  void setCreationMoleNumbers(const double* const n_k, const std::vector<size_t> &creationGlobalRxnNumbers);
329 
330  //! Return a const reference to the creationMoleNumbers stored in the object.
331  /*!
332  * @returns a const reference to the vector of creationMoleNumbers
333  */
334  const vector_fp & creationMoleNumbers(std::vector<size_t> &creationGlobalRxnNumbers) const;
335 
336  //! Returns whether the phase is an ideal solution phase
337  bool isIdealSoln() const;
338 
339  //! Return the index of the species that represents the the voltage of the
340  //! phase
341  size_t phiVarIndex() const;
342 
343  void setPhiVarIndex(size_t phiVarIndex);
344 
345  //! Retrieve the kth Species structure for the species belonging to this
346  //! phase
347  /*!
348  * The index into this vector is the species index within the phase.
349  *
350  * @param kindex kth species index.
351  */
352  vcs_SpeciesProperties* speciesProperty(const size_t kindex);
353 
354  //! int indicating whether the phase exists or not
355  /*!
356  * returns the m_existence int for the phase
357  *
358  * - VCS_PHASE_EXIST_ZEROEDPHASE = -6: Set to not exist by fiat from a
359  * higher level. This is used in phase stability boundary calculations
360  * - VCS_PHASE_EXIST_NO = 0: Doesn't exist currently
361  * - VCS_PHASE_EXIST_MINORCONC = 1: Exists, but the concentration is so low
362  * that an alternate method is used to calculate the total phase
363  * concentrations.
364  * - VCS_PHASE_EXIST_YES = 2 : Does exist currently
365  * - VCS_PHASE_EXIST_ALWAYS = 3: Always exists because it contains inerts
366  * which can't exist in any other phase. Or, the phase exists always
367  * because it consists of a single species, which is identified with the
368  * voltage, i.e., it's an electron metal phase.
369  */
370  int exists() const;
371 
372  //! Set the existence flag in the object
373  /*!
374  * Note the total moles of the phase must have been set appropriately before
375  * calling this routine.
376  *
377  * @param existence Phase existence flag
378  *
379  * @note try to eliminate this routine
380  */
381  void setExistence(const int existence);
382 
383  //! Return the Global VCS index of the kth species in the phase
384  /*!
385  * @param spIndex local species index (0 to the number of species in the
386  * phase)
387  *
388  * @returns the VCS_SOLVE species index of the species. This changes as
389  * rearrangements are carried out.
390  */
391  size_t spGlobalIndexVCS(const size_t spIndex) const;
392 
393 
394  //! set the Global VCS index of the kth species in the phase
395  /*!
396  * @param spIndex local species index (0 to the number of species
397  * in the phase)
398  *
399  * @returns the VCS_SOLVE species index of the that species This changes as
400  * rearrangements are carried out.
401  */
402  void setSpGlobalIndexVCS(const size_t spIndex, const size_t spGlobalIndex);
403 
404  //! Sets the total moles of inert in the phase
405  /*!
406  * @param tMolesInert Value of the total kmols of inert species in the
407  * phase.
408  */
409  void setTotalMolesInert(const double tMolesInert);
410 
411  //! Returns the value of the total kmol of inert in the phase
412  double totalMolesInert() const;
413 
414  //! Returns the global index of the local element index for the phase
415  size_t elemGlobalIndex(const size_t e) const;
416 
417  //! sets a local phase element to a global index value
418  /*!
419  * @param eLocal Local phase element index
420  * @param eGlobal Global phase element index
421  */
422  void setElemGlobalIndex(const size_t eLocal, const size_t eGlobal);
423 
424  //! Returns the number of element constraints
425  size_t nElemConstraints() const;
426 
427  //! Name of the element constraint with index \c e.
428  /*!
429  * @param e Element index.
430  */
431  std::string elementName(const size_t e) const;
432 
433  //! Type of the element constraint with index \c e.
434  /*!
435  * @param e Element index.
436  */
437  int elementType(const size_t e) const;
438 
439  //! Transfer all of the element information from the ThermoPhase object to
440  //! the vcs_VolPhase object.
441  /*!
442  * Also decide whether we need a new charge neutrality element in the phase
443  * to enforce a charge neutrality constraint.
444  *
445  * @param tPhase Pointer to the ThermoPhase object
446  */
447  size_t transferElementsFM(const ThermoPhase* const tPhase);
448 
449  //! Get a constant form of the Species Formula Matrix
450  /*!
451  * Returns a `double**` pointer such that `fm[e][f]` is the formula
452  * matrix entry for element `e` for species `k`
453  */
454  const Array2D& getFormulaMatrix() const;
455 
456  //! Returns the type of the species unknown
457  /*!
458  * @param k species index
459  * @return the SpeciesUnknownType[k] = type of species
460  * - Normal -> VCS_SPECIES_TYPE_MOLUNK (unknown is the mole number in
461  * the phase)
462  * - metal electron -> VCS_SPECIES_INTERFACIALVOLTAGE (unknown is the
463  * interfacial voltage (volts))
464  */
465  int speciesUnknownType(const size_t k) const;
466 
467  int elementActive(const size_t e) const;
468 
469  //! Return the number of species in the phase
470  size_t nSpecies() const;
471 
472  //! Return the name corresponding to the equation of state
473  std::string eos_name() const;
474 
475 private:
476  //! Evaluate the activity coefficients at the current conditions
477  /*!
478  * We carry out a calculation whenever #m_UpToDate_AC is false. Specifically
479  * whenever a phase goes zero, we do not carry out calculations on it.
480  */
481  void _updateActCoeff() const;
482 
483  //! Gibbs free energy calculation for standard states
484  /*!
485  * Calculate the Gibbs free energies for the standard states The results are
486  * held internally within the object.
487  */
488  void _updateGStar() const;
489 
490  //! Gibbs free energy calculation at a temperature for the reference state
491  //! of each species
492  void _updateG0() const;
493 
494  //! Molar volume calculation for standard states
495  /*!
496  * Calculate the molar volume for the standard states. The results are held
497  * internally within the object. Units are in m**3/kmol.
498  */
499  void _updateVolStar() const;
500 
501  //! Calculate the partial molar volumes of all species and return the
502  //! total volume
503  /*!
504  * Calculates these quantities internally and then stores them
505  *
506  * @return total volume [m^3]
507  */
508  double _updateVolPM() const;
509 
510  //! Evaluation of Activity Coefficient Jacobians
511  /*!
512  * This is the derivative of the ln of the activity coefficient with respect
513  * to mole number of jth species. (temp, pressure, and other mole numbers
514  * held constant)
515  *
516  * We employ a finite difference derivative approach here. Because we have
517  * to change the mole numbers, this is not a const function, even though the
518  * paradigm would say that it should be.
519  */
520  void _updateLnActCoeffJac();
521 
522  //! Updates the mole fraction dependencies
523  /*!
524  * Whenever the mole fractions change, this routine should be called.
525  */
527 
528 private:
529  //! Backtrack value of VCS_SOLVE *
531 
532 public:
533  //! Original ID of the phase in the problem.
534  /*!
535  * If a non-ideal phase splits into two due to a miscibility gap, these
536  * numbers will stay the same after the split.
537  */
538  size_t VP_ID_;
539 
540  //! If true, this phase consists of a single species
542 
543  //! If true, this phase is a gas-phase like phase
544  /*!
545  * A RTlog(p/1atm) term is added onto the chemical potential for inert
546  * species if this is true.
547  */
549 
550  //! Type of the equation of state
551  /*!
552  * The known types are listed at the top of this file.
553  */
555 
556  //! This is the element number for the charge neutrality condition of the
557  //! phase
558  /*!
559  * If it has one. If it does not have a charge neutrality
560  * constraint, then this value is equal to -1
561  */
563 
564  //! Convention for the activity formulation
565  /*!
566  * * 0 = molar based activities (default)
567  * * 1 = Molality based activities, mu = mu_0 + ln a_molality. Standard
568  * state is based on unity molality
569  */
571 
572 private:
573  //! Number of element constraints within the problem
574  /*!
575  * This is usually equal to the number of elements.
576  */
578 
579  //! vector of strings containing the element constraint names
580  /*!
581  * Length = nElemConstraints
582  */
583  std::vector<std::string> m_elementNames;
584 
585  //! boolean indicating whether an element constraint is active
586  //! for the current problem
588 
589  //! Type of the element constraint
590  /*!
591  * m_elType[j] = type of the element:
592  * * 0 VCS_ELEM_TYPE_ABSPOS Normal element that is positive or zero in
593  * all species.
594  * * 1 VCS_ELEM_TYPE_ELECTRONCHARGE element dof that corresponds to the
595  * charge DOF.
596  * * 2 VCS_ELEM_TYPE_OTHERCONSTRAINT Other constraint which may mean that
597  * a species has neg 0 or pos value of that constraint (other than
598  * charge)
599  */
601 
602  //! Formula Matrix for the phase
603  /*!
604  * FormulaMatrix(kspec,j) = Formula Matrix for the species
605  * Number of elements, j, in the kspec species
606  */
608 
609  //! Type of the species unknown
610  /*!
611  * SpeciesUnknownType[k] = type of species
612  * - Normal -> VCS_SPECIES_TYPE_MOLUNK.
613  * (unknown is the mole number in the phase)
614  * - metal electron -> VCS_SPECIES_INTERFACIALVOLTAGE.
615  * (unknown is the interfacial voltage (volts))
616  */
618 
619  //! Index of the element number in the global list of elements stored in VCS_SOLVE
620  std::vector<size_t> m_elemGlobalIndex;
621 
622  //! Number of species in the phase
623  size_t m_numSpecies;
624 
625 public:
626  //! String name for the phase
627  std::string PhaseName;
628 
629 private:
630  //! Total moles of inert in the phase
632 
633  //! Boolean indicating whether the phase is an ideal solution
634  //! and therefore its molar-based activity coefficients are
635  //! uniformly equal to one.
637 
638  //! Current state of existence:
639  /*!
640  * - VCS_PHASE_EXIST_ZEROEDPHASE = -6: Set to not exist by fiat from a
641  * higher level. This is used in phase stability boundary calculations
642  * - VCS_PHASE_EXIST_NO = 0: Doesn't exist currently
643  * - VCS_PHASE_EXIST_MINORCONC = 1: Exists, but the concentration is so
644  * low that an alternate method is used to calculate the total phase
645  * concentrations.
646  * - VCS_PHASE_EXIST_YES = 2 : Does exist currently
647  * - VCS_PHASE_EXIST_ALWAYS = 3: Always exists because it contains inerts
648  * which can't exist in any other phase. Or, the phase exists always
649  * because it consists of a single species, which is identified with the
650  * voltage, i.e., its an electron metal phase.
651  */
653 
654  // Index of the first MF species in the list of unknowns for this phase
655  /*!
656  * This is always equal to zero.
657  * Am anticipating the case where the phase potential is species # 0,
658  * for multiphase phases. Right now we have the phase potential equal
659  * to 0 for single species phases, where we set by hand the mole fraction
660  * of species 0 to one.
661  */
663 
664  //! Index into the species vectors
665  /*!
666  * Maps the phase species number into the global species number.
667  * Note, as part of the vcs algorithm, the order of the species
668  * vector is changed during the algorithm
669  */
670  std::vector<size_t> IndSpecies;
671 
672  //! Vector of Species structures for the species belonging to this phase
673  /*!
674  * The index into this vector is the species index within the phase.
675  */
676  std::vector<vcs_SpeciesProperties*> ListSpeciesPtr;
677 
678  /**
679  * If we are using Cantera, this is the pointer to the ThermoPhase
680  * object. If not, this is null.
681  */
683 
684  //! Total mols in the phase. units are kmol
685  double v_totalMoles;
686 
687  //! Vector of the current mole fractions for species in the phase
689 
690  //! Vector of current creationMoleNumbers_
691  /*!
692  * These are the actual unknowns in the phase stability problem
693  */
695 
696  //! Vector of creation global reaction numbers for the phase stability problem
697  /*!
698  * The phase stability problem requires a global reaction number for each
699  * species in the phase. Usually this is the krxn = kglob - M for species in
700  * the phase that are not components. For component species, the choice of
701  * the reaction is one which maximizes the chance that the phase pops into
702  * (or remains in) existence.
703  *
704  * The index here is the local phase species index. the value of the
705  * variable is the global vcs reaction number. Note, that the global
706  * reaction number will go out of order when the species positions are
707  * swapped. So, this number has to be recalculated.
708  *
709  * Length = number of species in phase
710  */
711  std::vector<size_t> creationGlobalRxnNumbers_;
712 
713  //! If the potential is a solution variable in VCS, it acts as a species.
714  //! This is the species index in the phase for the potential
716 
717  //! Total Volume of the phase. Units are m**3.
718  mutable double m_totalVol;
719 
720  //! Vector of calculated SS0 chemical potentials for the
721  //! current Temperature.
722  /*!
723  * Note, This is the chemical potential derived strictly from the polynomial
724  * in temperature. Pressure effects have to be added in to get to the
725  * standard state. Units are J/kmol.
726  */
728 
729  //! Vector of calculated Star chemical potentials for the
730  //! current Temperature and pressure.
731  /*!
732  * Note, This is the chemical potential at unit activity. Thus, we can call
733  * it the standard state chemical potential as well. Units are J/kmol.
734  */
736 
737  //! Vector of the Star molar Volumes of the species. units m3 / kmol
739 
740  //! Vector of the Partial molar Volumes of the species. units m3 / kmol
742 
743  //! Vector of calculated activity coefficients for the current state
744  /*!
745  * Whether or not this vector is current is determined by the bool
746  * #m_UpToDate_AC.
747  */
749 
750  //! Vector of the derivatives of the ln activity coefficient wrt to the
751  //! current mole number multiplied by the current phase moles
752  /*!
753  * np_dLnActCoeffdMolNumber(k,j);
754  * - j = id of the species mole number
755  * - k = id of the species activity coefficient
756  */
758 
759  //! Status
760  /*!
761  * valid values are
762  * - VCS_STATECALC_OLD
763  * - VCS_STATECALC_NEW
764  * - VCS_STATECALC_TMP
765  */
767 
768  //! Value of the potential for the phase (Volts)
769  double m_phi;
770 
771  //! Boolean indicating whether the object has an up-to-date mole number vector
772  //! and potential with respect to the current vcs state calc status
774 
775  //! Boolean indicating whether activity coefficients are up to date.
776  /*!
777  * Activity coefficients and volume calculations are lagged. They are only
778  * called when they are needed (and when the state has changed so that they
779  * need to be recalculated).
780  */
781  mutable bool m_UpToDate_AC;
782 
783  //! Boolean indicating whether Star volumes are up to date.
784  /*!
785  * Activity coefficients and volume calculations are lagged. They are only
786  * called when they are needed (and when the state has changed so that they
787  * need to be recalculated). Star volumes are sensitive to temperature and
788  * pressure
789  */
790  mutable bool m_UpToDate_VolStar;
791 
792  //! Boolean indicating whether partial molar volumes are up to date.
793  /*!
794  * Activity coefficients and volume calculations are lagged. They are only
795  * called when they are needed (and when the state has changed so that they
796  * need to be recalculated). partial molar volumes are sensitive to
797  * everything
798  */
799  mutable bool m_UpToDate_VolPM;
800 
801  //! Boolean indicating whether GStar is up to date.
802  /*!
803  * GStar is sensitive to the temperature and the pressure, only
804  */
805  mutable bool m_UpToDate_GStar;
806 
807  //! Boolean indicating whether G0 is up to date.
808  /*!
809  * G0 is sensitive to the temperature and the pressure, only
810  */
811  mutable bool m_UpToDate_G0;
812 
813  //! Current value of the temperature for this object, and underlying objects
814  double Temp_;
815 
816  //! Current value of the pressure for this object, and underlying objects
817  double Pres_;
818 };
819 
820 }
821 
822 #endif
void setCreationMoleNumbers(const double *const n_k, const std::vector< size_t > &creationGlobalRxnNumbers)
Sets the creationMoleNum&#39;s within the phase object.
bool isIdealSoln() const
Returns whether the phase is an ideal solution phase.
void sendToVCS_GStar(double *const gstar) const
Fill in the partial molar volume vector for VCS.
void _updateMoleFractionDependencies()
Updates the mole fraction dependencies.
bool m_UpToDate_VolStar
Boolean indicating whether Star volumes are up to date.
Definition: vcs_VolPhase.h:790
void _updateVolStar() const
Molar volume calculation for standard states.
double m_phi
Value of the potential for the phase (Volts)
Definition: vcs_VolPhase.h:769
bool m_singleSpecies
If true, this phase consists of a single species.
Definition: vcs_VolPhase.h:541
void updateFromVCS_MoleNumbers(const int stateCalc)
Update the moles within the phase, if necessary.
double totalMolesInert() const
Returns the value of the total kmol of inert in the phase.
void setMoleFractionsState(const double molNum, const double *const moleFracVec, const int vcsStateStatus)
Set the moles and/or mole fractions within the phase.
vector_fp StarChemicalPotential
Vector of calculated Star chemical potentials for the current Temperature and pressure.
Definition: vcs_VolPhase.h:735
int m_vcsStateStatus
Status.
Definition: vcs_VolPhase.h:766
size_t nElemConstraints() const
Returns the number of element constraints.
void setMolesFromVCSCheck(const int vcsStateStatus, const double *molesSpeciesVCS, const double *const TPhMoles)
Set the moles within the phase.
void _updateActCoeff() const
Evaluate the activity coefficients at the current conditions.
int m_existence
Current state of existence:
Definition: vcs_VolPhase.h:652
int p_activityConvention
Convention for the activity formulation.
Definition: vcs_VolPhase.h:570
Array2D m_formulaMatrix
Formula Matrix for the phase.
Definition: vcs_VolPhase.h:607
double VolStar_calc_one(size_t kspec) const
Molar volume calculation for standard state of one species.
vector_fp StarMolarVol
Vector of the Star molar Volumes of the species. units m3 / kmol.
Definition: vcs_VolPhase.h:738
double electricPotential() const
Returns the electric field of the phase.
const Array2D & getFormulaMatrix() const
Get a constant form of the Species Formula Matrix.
size_t spGlobalIndexVCS(const size_t spIndex) const
Return the Global VCS index of the kth species in the phase.
size_t transferElementsFM(const ThermoPhase *const tPhase)
Transfer all of the element information from the ThermoPhase object to the vcs_VolPhase object...
const vector_fp & creationMoleNumbers(std::vector< size_t > &creationGlobalRxnNumbers) const
Return a const reference to the creationMoleNumbers stored in the object.
std::string elementName(const size_t e) const
Name of the element constraint with index e.
std::string PhaseName
String name for the phase.
Definition: vcs_VolPhase.h:627
size_t ChargeNeutralityElement
This is the element number for the charge neutrality condition of the phase.
Definition: vcs_VolPhase.h:562
vector_int m_elementActive
boolean indicating whether an element constraint is active for the current problem ...
Definition: vcs_VolPhase.h:587
VCS_SOLVE * m_owningSolverObject
Backtrack value of VCS_SOLVE *.
Definition: vcs_VolPhase.h:530
size_t elemGlobalIndex(const size_t e) const
Returns the global index of the local element index for the phase.
vector_fp SS0ChemicalPotential
Vector of calculated SS0 chemical potentials for the current Temperature.
Definition: vcs_VolPhase.h:727
void setExistence(const int existence)
Set the existence flag in the object.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:31
void setPtrThermoPhase(ThermoPhase *tp_ptr)
Set the pointer for Cantera&#39;s ThermoPhase parameter.
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:93
bool m_UpToDate_VolPM
Boolean indicating whether partial molar volumes are up to date.
Definition: vcs_VolPhase.h:799
std::vector< vcs_SpeciesProperties * > ListSpeciesPtr
Vector of Species structures for the species belonging to this phase.
Definition: vcs_VolPhase.h:676
Header file for class Cantera::Array2D.
size_t VP_ID_
Original ID of the phase in the problem.
Definition: vcs_VolPhase.h:538
double m_totalVol
Total Volume of the phase. Units are m**3.
Definition: vcs_VolPhase.h:718
void setMolesFromVCS(const int stateCalc, const double *molesSpeciesVCS=0)
Set the moles within the phase.
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:159
double GStar_calc_one(size_t kspec) const
Gibbs free energy calculation for standard state of one species.
void _updateGStar() const
Gibbs free energy calculation for standard states.
double v_totalMoles
Total mols in the phase. units are kmol.
Definition: vcs_VolPhase.h:685
std::string eos_name() const
Return the name corresponding to the equation of state.
Properties of a single species.
size_t m_numElemConstraints
Number of element constraints within the problem.
Definition: vcs_VolPhase.h:577
int m_eqnState
Type of the equation of state.
Definition: vcs_VolPhase.h:554
void setElectricPotential(const double phi)
set the electric potential of the phase
size_t nSpecies() const
Return the number of species in the phase.
int elementType(const size_t e) const
Type of the element constraint with index e.
void sendToVCS_ActCoeff(const int stateCalc, double *const AC)
Fill in an activity coefficients vector within a VCS_SOLVE object.
void setMoleFractions(const double *const xmol)
Set the mole fractions from a conventional mole fraction vector.
void _updateLnActCoeffJac()
Evaluation of Activity Coefficient Jacobians.
std::vector< size_t > creationGlobalRxnNumbers_
Vector of creation global reaction numbers for the phase stability problem.
Definition: vcs_VolPhase.h:711
void resize(const size_t phaseNum, const size_t numSpecies, const size_t numElem, const char *const phaseName, const double molesInert=0.0)
The resize() function fills in all of the initial information if it is not given in the constructor...
bool m_UpToDate_G0
Boolean indicating whether G0 is up to date.
Definition: vcs_VolPhase.h:811
bool m_gasPhase
If true, this phase is a gas-phase like phase.
Definition: vcs_VolPhase.h:548
double totalMoles() const
Return the total moles in the phase.
double Temp_
Current value of the temperature for this object, and underlying objects.
Definition: vcs_VolPhase.h:814
double Pres_
Current value of the pressure for this object, and underlying objects.
Definition: vcs_VolPhase.h:817
size_t m_numSpecies
Number of species in the phase.
Definition: vcs_VolPhase.h:623
size_t phiVarIndex() const
Return the index of the species that represents the the voltage of the phase.
void setMolesOutOfDate(int stateCalc=-1)
Sets the mole flag within the object to out of date.
Phase information and Phase calculations for vcs.
Definition: vcs_VolPhase.h:81
bool m_UpToDate_AC
Boolean indicating whether activity coefficients are up to date.
Definition: vcs_VolPhase.h:781
ThermoPhase * TP_ptr
If we are using Cantera, this is the pointer to the ThermoPhase object.
Definition: vcs_VolPhase.h:682
vector_fp PartialMolarVol
Vector of the Partial molar Volumes of the species. units m3 / kmol.
Definition: vcs_VolPhase.h:741
vector_fp ActCoeff
Vector of calculated activity coefficients for the current state.
Definition: vcs_VolPhase.h:748
int exists() const
int indicating whether the phase exists or not
std::vector< size_t > IndSpecies
Index into the species vectors.
Definition: vcs_VolPhase.h:670
bool m_UpToDate_GStar
Boolean indicating whether GStar is up to date.
Definition: vcs_VolPhase.h:805
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:157
size_t m_phiVarIndex
If the potential is a solution variable in VCS, it acts as a species.
Definition: vcs_VolPhase.h:715
const ThermoPhase * ptrThermoPhase() const
Return a const ThermoPhase pointer corresponding to this phase.
double m_totalMolesInert
Total moles of inert in the phase.
Definition: vcs_VolPhase.h:631
vcs_SpeciesProperties * speciesProperty(const size_t kindex)
Retrieve the kth Species structure for the species belonging to this phase.
void setSpGlobalIndexVCS(const size_t spIndex, const size_t spGlobalIndex)
set the Global VCS index of the kth species in the phase
double sendToVCS_VolPM(double *const VolPM) const
Fill in the partial molar volume vector for VCS.
vector_fp Xmol_
Vector of the current mole fractions for species in the phase.
Definition: vcs_VolPhase.h:688
void setMolesCurrent(int vcsStateStatus)
Sets the mole flag within the object to be current.
bool m_UpToDate
Boolean indicating whether the object has an up-to-date mole number vector and potential with respect...
Definition: vcs_VolPhase.h:773
vector_int m_elementType
Type of the element constraint.
Definition: vcs_VolPhase.h:600
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:8
void setState_T(const double temperature_Kelvin)
Sets the temperature in this object and underlying ThermoPhase objects.
void _updateG0() const
Gibbs free energy calculation at a temperature for the reference state of each species.
const vector_fp & moleFractions() const
Return a const reference to the mole fractions stored in the object.
void setState_TP(const double temperature_Kelvin, const double pressure_PA)
Sets the temperature and pressure in this object and underlying ThermoPhase objects.
bool m_isIdealSoln
Boolean indicating whether the phase is an ideal solution and therefore its molar-based activity coef...
Definition: vcs_VolPhase.h:636
std::vector< std::string > m_elementNames
vector of strings containing the element constraint names
Definition: vcs_VolPhase.h:583
Array2D np_dLnActCoeffdMolNumber
Vector of the derivatives of the ln activity coefficient wrt to the current mole number multiplied by...
Definition: vcs_VolPhase.h:757
vector_int m_speciesUnknownType
Type of the species unknown.
Definition: vcs_VolPhase.h:617
int speciesUnknownType(const size_t k) const
Returns the type of the species unknown.
This is the main structure used to hold the internal data used in vcs_solve_TP(), and to solve TP sys...
Definition: vcs_solve.h:44
double molefraction(size_t kspec) const
Returns the mole fraction of the kspec species.
double _updateVolPM() const
Calculate the partial molar volumes of all species and return the total volume.
void setTotalMolesInert(const double tMolesInert)
Sets the total moles of inert in the phase.
std::vector< size_t > m_elemGlobalIndex
Index of the element number in the global list of elements stored in VCS_SOLVE.
Definition: vcs_VolPhase.h:620
double G0_calc_one(size_t kspec) const
Gibbs free energy calculation at a temperature for the reference state of a species, return a value for one species.
void setTotalMoles(const double totalMols)
Sets the total moles in the phase.
void setElemGlobalIndex(const size_t eLocal, const size_t eGlobal)
sets a local phase element to a global index value
vector_fp creationMoleNumbers_
Vector of current creationMoleNumbers_.
Definition: vcs_VolPhase.h:694