Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
vcs_VolPhase.cpp
Go to the documentation of this file.
1 /**
2  * @file vcs_VolPhase.cpp
3  */
4 /*
5  * Copyright (2005) Sandia Corporation. Under the terms of
6  * Contract DE-AC04-94AL85000 with Sandia Corporation, the
7  * U.S. Government retains certain rights in this software.
8  */
12 
16 
17 #include <sstream>
18 #include <cstdio>
19 
20 namespace Cantera
21 {
22 
23 vcs_VolPhase::vcs_VolPhase(VCS_SOLVE* owningSolverObject) :
24  m_owningSolverObject(0),
25  VP_ID_(npos),
26  m_singleSpecies(true),
27  m_gasPhase(false),
28  m_eqnState(VCS_EOS_CONSTANT),
29  ChargeNeutralityElement(npos),
30  p_VCS_UnitsFormat(VCS_UNITS_MKS),
31  p_activityConvention(0),
32  m_numElemConstraints(0),
33  m_elemGlobalIndex(0),
34  m_numSpecies(0),
35  m_totalMolesInert(0.0),
36  m_isIdealSoln(false),
37  m_existence(VCS_PHASE_EXIST_NO),
38  m_MFStartIndex(0),
39  IndSpecies(0),
40  m_useCanteraCalls(false),
41  TP_ptr(0),
42  v_totalMoles(0.0),
43  m_phiVarIndex(npos),
44  m_totalVol(0.0),
45  m_vcsStateStatus(VCS_STATECALC_OLD),
46  m_phi(0.0),
47  m_UpToDate(false),
48  m_UpToDate_AC(false),
49  m_UpToDate_VolStar(false),
50  m_UpToDate_VolPM(false),
51  m_UpToDate_GStar(false),
52  m_UpToDate_G0(false),
53  Temp_(273.15),
54  Pres_(1.01325E5)
55 {
56  m_owningSolverObject = owningSolverObject;
57 }
58 
59 vcs_VolPhase::~vcs_VolPhase()
60 {
61  for (size_t k = 0; k < m_numSpecies; k++) {
62  vcs_SpeciesProperties* sp = ListSpeciesPtr[k];
63  delete sp;
64  sp = 0;
65  }
66 }
67 
68 vcs_VolPhase::vcs_VolPhase(const vcs_VolPhase& b) :
69  m_owningSolverObject(b.m_owningSolverObject),
70  VP_ID_(b.VP_ID_),
71  m_singleSpecies(b.m_singleSpecies),
72  m_gasPhase(b.m_gasPhase),
73  m_eqnState(b.m_eqnState),
74  ChargeNeutralityElement(b.ChargeNeutralityElement),
75  p_VCS_UnitsFormat(b.p_VCS_UnitsFormat),
76  p_activityConvention(b.p_activityConvention),
77  m_numElemConstraints(b.m_numElemConstraints),
78  m_numSpecies(b.m_numSpecies),
79  m_totalMolesInert(b.m_totalMolesInert),
80  m_isIdealSoln(b.m_isIdealSoln),
81  m_existence(b.m_existence),
82  m_MFStartIndex(b.m_MFStartIndex),
83  m_useCanteraCalls(b.m_useCanteraCalls),
84  TP_ptr(b.TP_ptr),
85  v_totalMoles(b.v_totalMoles),
86  creationMoleNumbers_(0),
87  creationGlobalRxnNumbers_(0),
88  m_phiVarIndex(npos),
89  m_totalVol(b.m_totalVol),
90  m_vcsStateStatus(VCS_STATECALC_OLD),
91  m_phi(b.m_phi),
92  m_UpToDate(false),
93  m_UpToDate_AC(false),
94  m_UpToDate_VolStar(false),
95  m_UpToDate_VolPM(false),
96  m_UpToDate_GStar(false),
97  m_UpToDate_G0(false),
98  Temp_(b.Temp_),
99  Pres_(b.Pres_)
100 {
101  //! Objects that are owned by this object are deep copied here, except for
102  //! the ThermoPhase object. The assignment operator does most of the work.
103  *this = b;
104 }
105 
106 vcs_VolPhase& vcs_VolPhase::operator=(const vcs_VolPhase& b)
107 {
108  if (&b != this) {
109  size_t old_num = m_numSpecies;
110 
111  // Note: we comment this out for the assignment operator
112  // specifically, because it isn't true for the assignment
113  // operator but is true for a copy constructor
114  // m_owningSolverObject = b.m_owningSolverObject;
115 
116  VP_ID_ = b.VP_ID_;
117  m_singleSpecies = b.m_singleSpecies;
118  m_gasPhase = b.m_gasPhase;
119  m_eqnState = b.m_eqnState;
120  ChargeNeutralityElement = b.ChargeNeutralityElement;
121  p_VCS_UnitsFormat = b.p_VCS_UnitsFormat;
122  p_activityConvention= b.p_activityConvention;
123  m_numSpecies = b.m_numSpecies;
124  m_numElemConstraints = b.m_numElemConstraints;
125  m_elementNames.resize(b.m_numElemConstraints);
126  for (size_t e = 0; e < b.m_numElemConstraints; e++) {
127  m_elementNames[e] = b.m_elementNames[e];
128  }
129  m_elementActive = b.m_elementActive;
130  m_elementType = b.m_elementType;
131  m_formulaMatrix = b.m_formulaMatrix;
132  m_speciesUnknownType = b.m_speciesUnknownType;
133  m_elemGlobalIndex = b.m_elemGlobalIndex;
134  PhaseName = b.PhaseName;
135  m_totalMolesInert = b.m_totalMolesInert;
136  m_isIdealSoln = b.m_isIdealSoln;
137  m_existence = b.m_existence;
138  m_MFStartIndex = b.m_MFStartIndex;
139  /*
140  * Do a shallow copy because we haven' figured this out.
141  */
142  IndSpecies = b.IndSpecies;
143 
144  for (size_t k = 0; k < old_num; k++) {
145  if (ListSpeciesPtr[k]) {
146  delete ListSpeciesPtr[k];
147  ListSpeciesPtr[k] = 0;
148  }
149  }
150  ListSpeciesPtr.resize(m_numSpecies, 0);
151  for (size_t k = 0; k < m_numSpecies; k++) {
152  ListSpeciesPtr[k] =
153  new vcs_SpeciesProperties(*(b.ListSpeciesPtr[k]));
154  }
155  m_useCanteraCalls = b.m_useCanteraCalls;
156  /*
157  * Do a shallow copy of the ThermoPhase object pointer.
158  * We don't duplicate the object.
159  * Um, there is no reason we couldn't do a
160  * duplicateMyselfAsThermoPhase() call here. This will
161  * have to be looked into.
162  */
163  TP_ptr = b.TP_ptr;
164  v_totalMoles = b.v_totalMoles;
165  Xmol_ = b.Xmol_;
166  creationMoleNumbers_ = b.creationMoleNumbers_;
167  creationGlobalRxnNumbers_ = b.creationGlobalRxnNumbers_;
168  m_phiVarIndex = b.m_phiVarIndex;
169  m_totalVol = b.m_totalVol;
170  SS0ChemicalPotential = b.SS0ChemicalPotential;
171  StarChemicalPotential = b.StarChemicalPotential;
172  StarMolarVol = b.StarMolarVol;
173  PartialMolarVol = b.PartialMolarVol;
174  ActCoeff = b.ActCoeff;
175  np_dLnActCoeffdMolNumber = b.np_dLnActCoeffdMolNumber;
176  m_vcsStateStatus = b.m_vcsStateStatus;
177  m_phi = b.m_phi;
178  m_UpToDate = false;
179  m_UpToDate_AC = false;
180  m_UpToDate_VolStar = false;
181  m_UpToDate_VolPM = false;
182  m_UpToDate_GStar = false;
183  m_UpToDate_G0 = false;
184  Temp_ = b.Temp_;
185  Pres_ = b.Pres_;
186 
187  setState_TP(Temp_, Pres_);
188  _updateMoleFractionDependencies();
189  }
190  return *this;
191 }
192 
193 void vcs_VolPhase::resize(const size_t phaseNum, const size_t nspecies,
194  const size_t numElem, const char* const phaseName,
195  const double molesInert)
196 {
197  AssertThrowMsg(nspecies > 0, "vcs_VolPhase::resize", "nspecies Error");
198  setTotalMolesInert(molesInert);
199  m_phi = 0.0;
200  m_phiVarIndex = npos;
201 
202  if (phaseNum == VP_ID_) {
203  if (strcmp(PhaseName.c_str(), phaseName)) {
204  throw CanteraError("vcs_VolPhase::resize",
205  "Strings are different: " + PhaseName + " " +
206  phaseName + " :unknown situation");
207  }
208  } else {
209  VP_ID_ = phaseNum;
210  if (!phaseName) {
211  std::stringstream sstmp;
212  sstmp << "Phase_" << VP_ID_;
213  PhaseName = sstmp.str();
214  } else {
215  PhaseName = phaseName;
216  }
217  }
218  if (nspecies > 1) {
219  m_singleSpecies = false;
220  } else {
221  m_singleSpecies = true;
222  }
223 
224  if (m_numSpecies == nspecies && numElem == m_numElemConstraints) {
225  return;
226  }
227 
228  m_numSpecies = nspecies;
229  if (nspecies > 1) {
230  m_singleSpecies = false;
231  }
232 
233 
234  IndSpecies.resize(nspecies, npos);
235 
236  if (ListSpeciesPtr.size() >= m_numSpecies) {
237  for (size_t i = 0; i < m_numSpecies; i++) {
238  if (ListSpeciesPtr[i]) {
239  delete ListSpeciesPtr[i];
240  ListSpeciesPtr[i] = 0;
241  }
242  }
243  }
244  ListSpeciesPtr.resize(nspecies, 0);
245  for (size_t i = 0; i < nspecies; i++) {
246  ListSpeciesPtr[i] = new vcs_SpeciesProperties(phaseNum, i, this);
247  }
248 
249  Xmol_.resize(nspecies, 0.0);
250  creationMoleNumbers_.resize(nspecies, 0.0);
251  creationGlobalRxnNumbers_.resize(nspecies, npos);
252  for (size_t i = 0; i < nspecies; i++) {
253  Xmol_[i] = 1.0/nspecies;
254  creationMoleNumbers_[i] = 1.0/nspecies;
255  if (IndSpecies[i] >= m_numElemConstraints) {
256  creationGlobalRxnNumbers_[i] = IndSpecies[i] - m_numElemConstraints;
257  } else {
258  creationGlobalRxnNumbers_[i] = npos;
259  }
260  }
261 
262  SS0ChemicalPotential.resize(nspecies, -1.0);
263  StarChemicalPotential.resize(nspecies, -1.0);
264  StarMolarVol.resize(nspecies, -1.0);
265  PartialMolarVol.resize(nspecies, -1.0);
266  ActCoeff.resize(nspecies, 1.0);
267  np_dLnActCoeffdMolNumber.resize(nspecies, nspecies, 0.0);
268 
269 
270  m_speciesUnknownType.resize(nspecies, VCS_SPECIES_TYPE_MOLNUM);
271  m_UpToDate = false;
272  m_vcsStateStatus = VCS_STATECALC_OLD;
273  m_UpToDate_AC = false;
274  m_UpToDate_VolStar = false;
275  m_UpToDate_VolPM = false;
276  m_UpToDate_GStar = false;
277  m_UpToDate_G0 = false;
278 
279 
280  elemResize(numElem);
281 
282 }
283 
284 void vcs_VolPhase::elemResize(const size_t numElemConstraints)
285 {
286 
287  m_elementNames.resize(numElemConstraints);
288 
289  m_elementActive.resize(numElemConstraints+1, 1);
290  m_elementType.resize(numElemConstraints, VCS_ELEM_TYPE_ABSPOS);
291  m_formulaMatrix.resize(m_numSpecies, numElemConstraints, 0.0);
292 
293  m_elementNames.resize(numElemConstraints, "");
294  m_elemGlobalIndex.resize(numElemConstraints, npos);
295 
296  m_numElemConstraints = numElemConstraints;
297 }
298 
299 void vcs_VolPhase::_updateActCoeff() const
300 {
301  if (m_isIdealSoln) {
302  m_UpToDate_AC = true;
303  return;
304  }
305  if (m_useCanteraCalls) {
306  TP_ptr->getActivityCoefficients(VCS_DATA_PTR(ActCoeff));
307  } else {
308  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
309  "deprecated and will not work after Cantera 2.2.");
310  }
311  m_UpToDate_AC = true;
312 }
313 
314 double vcs_VolPhase::AC_calc_one(size_t kspec) const
315 {
316  if (! m_UpToDate_AC) {
317  _updateActCoeff();
318  }
319  return ActCoeff[kspec];
320 }
321 
322 void vcs_VolPhase::_updateG0() const
323 {
324  if (m_useCanteraCalls) {
325  TP_ptr->getGibbs_ref(VCS_DATA_PTR(SS0ChemicalPotential));
326  } else {
327  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
328  "deprecated and will not work after Cantera 2.2.");
329  double R = vcsUtil_gasConstant(p_VCS_UnitsFormat);
330  for (size_t k = 0; k < m_numSpecies; k++) {
331  size_t kglob = IndSpecies[k];
332  vcs_SpeciesProperties* sProp = ListSpeciesPtr[k];
333  VCS_SPECIES_THERMO* sTherm = sProp->SpeciesThermo;
334  SS0ChemicalPotential[k] =
335  R * (sTherm->G0_R_calc(kglob, Temp_));
336  }
337  }
338  m_UpToDate_G0 = true;
339 }
340 
341 double vcs_VolPhase::G0_calc_one(size_t kspec) const
342 {
343  if (!m_UpToDate_G0) {
344  _updateG0();
345  }
346  return SS0ChemicalPotential[kspec];
347 }
348 
349 void vcs_VolPhase::_updateGStar() const
350 {
351  if (m_useCanteraCalls) {
352  TP_ptr->getStandardChemPotentials(VCS_DATA_PTR(StarChemicalPotential));
353  } else {
354  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
355  "deprecated and will not work after Cantera 2.2.");
356  double R = vcsUtil_gasConstant(p_VCS_UnitsFormat);
357  for (size_t k = 0; k < m_numSpecies; k++) {
358  size_t kglob = IndSpecies[k];
359  vcs_SpeciesProperties* sProp = ListSpeciesPtr[k];
360  VCS_SPECIES_THERMO* sTherm = sProp->SpeciesThermo;
361  StarChemicalPotential[k] =
362  R * (sTherm->GStar_R_calc(kglob, Temp_, Pres_));
363  }
364  }
365  m_UpToDate_GStar = true;
366 }
367 
368 double vcs_VolPhase::GStar_calc_one(size_t kspec) const
369 {
370  if (!m_UpToDate_GStar) {
371  _updateGStar();
372  }
373  return StarChemicalPotential[kspec];
374 }
375 
376 void vcs_VolPhase::setMoleFractions(const double* const xmol)
377 {
378  double sum = -1.0;
379  for (size_t k = 0; k < m_numSpecies; k++) {
380  Xmol_[k] = xmol[k];
381  sum+= xmol[k];
382  }
383  if (std::fabs(sum) > 1.0E-13) {
384  for (size_t k = 0; k < m_numSpecies; k++) {
385  Xmol_[k] /= sum;
386  }
387  }
388  _updateMoleFractionDependencies();
389  m_UpToDate = false;
390  m_vcsStateStatus = VCS_STATECALC_TMP;
391 }
392 
393 void vcs_VolPhase::_updateMoleFractionDependencies()
394 {
395  if (m_useCanteraCalls) {
396  if (TP_ptr) {
397  TP_ptr->setState_PX(Pres_, &(Xmol_[m_MFStartIndex]));
398  }
399  } else {
400  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
401  "deprecated and will not work after Cantera 2.2.");
402  }
403  if (!m_isIdealSoln) {
404  m_UpToDate_AC = false;
405  m_UpToDate_VolPM = false;
406  }
407 }
408 
409 const std::vector<double> & vcs_VolPhase::moleFractions() const
410 {
411  return Xmol_;
412 }
413 
414 double vcs_VolPhase::moleFraction(size_t k) const
415 {
416  return Xmol_[k];
417 }
418 
419 void vcs_VolPhase::setMoleFractionsState(const double totalMoles,
420  const double* const moleFractions,
421  const int vcsStateStatus)
422 {
423 
424  if (totalMoles != 0.0) {
425  // There are other ways to set the mole fractions when VCS_STATECALC
426  // is set to a normal settting.
427  if (vcsStateStatus != VCS_STATECALC_TMP) {
428  throw CanteraError("vcs_VolPhase::setMolesFractionsState",
429  "inappropriate usage");
430  }
431  m_UpToDate = false;
432  m_vcsStateStatus = VCS_STATECALC_TMP;
433  if (m_existence == VCS_PHASE_EXIST_ZEROEDPHASE) {
434  throw CanteraError("vcs_VolPhase::setMolesFractionsState",
435  "inappropriate usage");
436  }
437  m_existence = VCS_PHASE_EXIST_YES;
438  } else {
439  m_UpToDate = true;
440  m_vcsStateStatus = vcsStateStatus;
441  m_existence = std::min(m_existence, VCS_PHASE_EXIST_NO);
442  }
443  double fractotal = 1.0;
444  v_totalMoles = totalMoles;
445  if (m_totalMolesInert > 0.0) {
446  if (m_totalMolesInert > v_totalMoles) {
447  throw CanteraError("vcs_VolPhase::setMolesFractionsState",
448  "inerts greater than total: " + fp2str(v_totalMoles) + " " +
449  fp2str(m_totalMolesInert));
450  }
451  fractotal = 1.0 - m_totalMolesInert/v_totalMoles;
452  }
453  double sum = 0.0;
454  for (size_t k = 0; k < m_numSpecies; k++) {
455  Xmol_[k] = moleFractions[k];
456  sum += moleFractions[k];
457  }
458  if (sum == 0.0) {
459  throw CanteraError("vcs_VolPhase::setMolesFractionsState",
460  "inappropriate usage");
461  }
462  if (sum != fractotal) {
463  for (size_t k = 0; k < m_numSpecies; k++) {
464  Xmol_[k] *= (fractotal /sum);
465  }
466  }
467  _updateMoleFractionDependencies();
468 
469 }
470 
471 void vcs_VolPhase::setMolesFromVCS(const int stateCalc,
472  const double* molesSpeciesVCS)
473 {
474  v_totalMoles = m_totalMolesInert;
475 
476  if (molesSpeciesVCS == 0) {
477  AssertThrowMsg(m_owningSolverObject, "vcs_VolPhase::setMolesFromVCS",
478  "shouldn't be here");
479  if (stateCalc == VCS_STATECALC_OLD) {
480  molesSpeciesVCS = VCS_DATA_PTR(m_owningSolverObject->m_molNumSpecies_old);
481  } else if (stateCalc == VCS_STATECALC_NEW) {
482  molesSpeciesVCS = VCS_DATA_PTR(m_owningSolverObject->m_molNumSpecies_new);
483  } else if (DEBUG_MODE_ENABLED) {
484  throw CanteraError("vcs_VolPhase::setMolesFromVCS", "shouldn't be here"); }
485  } else if (DEBUG_MODE_ENABLED && m_owningSolverObject) {
486  if (stateCalc == VCS_STATECALC_OLD) {
487  if (molesSpeciesVCS != VCS_DATA_PTR(m_owningSolverObject->m_molNumSpecies_old)) {
488  throw CanteraError("vcs_VolPhase::setMolesFromVCS", "shouldn't be here");
489  }
490  } else if (stateCalc == VCS_STATECALC_NEW) {
491  if (molesSpeciesVCS != VCS_DATA_PTR(m_owningSolverObject->m_molNumSpecies_new)) {
492  throw CanteraError("vcs_VolPhase::setMolesFromVCS", "shouldn't be here");
493  }
494  }
495  }
496 
497  for (size_t k = 0; k < m_numSpecies; k++) {
498  if (m_speciesUnknownType[k] != VCS_SPECIES_TYPE_INTERFACIALVOLTAGE) {
499  size_t kglob = IndSpecies[k];
500  v_totalMoles += std::max(0.0, molesSpeciesVCS[kglob]);
501  }
502  }
503  if (v_totalMoles > 0.0) {
504  for (size_t k = 0; k < m_numSpecies; k++) {
505  if (m_speciesUnknownType[k] != VCS_SPECIES_TYPE_INTERFACIALVOLTAGE) {
506  size_t kglob = IndSpecies[k];
507  double tmp = std::max(0.0, molesSpeciesVCS[kglob]);
508  Xmol_[k] = tmp / v_totalMoles;
509  }
510  }
511  m_existence = VCS_PHASE_EXIST_YES;
512  } else {
513  // This is where we will start to store a better approximation
514  // for the mole fractions, when the phase doesn't exist.
515  // This is currently unimplemented.
516  m_existence = VCS_PHASE_EXIST_NO;
517  }
518  /*
519  * Update the electric potential if it is a solution variable
520  * in the equation system
521  */
522  if (m_phiVarIndex != npos) {
523  size_t kglob = IndSpecies[m_phiVarIndex];
524  if (m_numSpecies == 1) {
525  Xmol_[m_phiVarIndex] = 1.0;
526  } else {
527  Xmol_[m_phiVarIndex] = 0.0;
528  }
529  double phi = molesSpeciesVCS[kglob];
530  setElectricPotential(phi);
531  if (m_numSpecies == 1) {
532  m_existence = VCS_PHASE_EXIST_YES;
533  }
534  }
535  _updateMoleFractionDependencies();
536  if (m_totalMolesInert > 0.0) {
537  m_existence = VCS_PHASE_EXIST_ALWAYS;
538  }
539 
540  /*
541  * If stateCalc is old and the total moles is positive,
542  * then we have a valid state. If the phase went away, it would
543  * be a valid starting point for F_k's. So, save the state.
544  */
545  if (stateCalc == VCS_STATECALC_OLD) {
546  if (v_totalMoles > 0.0) {
547  creationMoleNumbers_ = Xmol_;
548  }
549  }
550 
551  /*
552  * Set flags indicating we are up to date with the VCS state vector.
553  */
554  m_UpToDate = true;
555  m_vcsStateStatus = stateCalc;
556 }
557 
558 void vcs_VolPhase::setMolesFromVCSCheck(const int vcsStateStatus,
559  const double* molesSpeciesVCS,
560  const double* const TPhMoles)
561 {
562  setMolesFromVCS(vcsStateStatus, molesSpeciesVCS);
563  /*
564  * Check for consistency with TPhMoles[]
565  */
566  double Tcheck = TPhMoles[VP_ID_];
567  if (Tcheck != v_totalMoles) {
568  if (vcs_doubleEqual(Tcheck, v_totalMoles)) {
569  Tcheck = v_totalMoles;
570  } else {
571  throw CanteraError("vcs_VolPhase::setMolesFromVCSCheck",
572  "We have a consistency problem: " + fp2str(Tcheck) + " " +
573  fp2str(v_totalMoles));
574  }
575  }
576 }
577 
578 void vcs_VolPhase::updateFromVCS_MoleNumbers(const int vcsStateStatus)
579 {
580  if (!m_UpToDate || (vcsStateStatus != m_vcsStateStatus)) {
581  if (vcsStateStatus == VCS_STATECALC_OLD || vcsStateStatus == VCS_STATECALC_NEW) {
582  if (m_owningSolverObject) {
583  setMolesFromVCS(vcsStateStatus);
584  }
585  }
586  }
587 }
588 
589 void vcs_VolPhase::sendToVCS_ActCoeff(const int vcsStateStatus,
590  double* const AC)
591 {
592  updateFromVCS_MoleNumbers(vcsStateStatus);
593  if (!m_UpToDate_AC) {
594  _updateActCoeff();
595  }
596  for (size_t k = 0; k < m_numSpecies; k++) {
597  size_t kglob = IndSpecies[k];
598  AC[kglob] = ActCoeff[k];
599  }
600 }
601 
602 double vcs_VolPhase::sendToVCS_VolPM(double* const VolPM) const
603 {
604  if (!m_UpToDate_VolPM) {
605  (void) _updateVolPM();
606  }
607  for (size_t k = 0; k < m_numSpecies; k++) {
608  size_t kglob = IndSpecies[k];
609  VolPM[kglob] = PartialMolarVol[k];
610  }
611  return m_totalVol;
612 }
613 
614 void vcs_VolPhase::sendToVCS_GStar(double* const gstar) const
615 {
616  if (!m_UpToDate_GStar) {
617  _updateGStar();
618  }
619  for (size_t k = 0; k < m_numSpecies; k++) {
620  size_t kglob = IndSpecies[k];
621  gstar[kglob] = StarChemicalPotential[k];
622  }
623 }
624 
625 void vcs_VolPhase::setElectricPotential(const double phi)
626 {
627  m_phi = phi;
628  if (m_useCanteraCalls) {
629  TP_ptr->setElectricPotential(m_phi);
630  } else {
631  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
632  "deprecated and will not work after Cantera 2.2.");
633 }
634  // We have changed the state variable. Set uptodate flags to false
635  m_UpToDate_AC = false;
636  m_UpToDate_VolStar = false;
637  m_UpToDate_VolPM = false;
638  m_UpToDate_GStar = false;
639 }
640 
641 double vcs_VolPhase::electricPotential() const
642 {
643  return m_phi;
644 }
645 
646 void vcs_VolPhase::setState_TP(const double temp, const double pres)
647 {
648  if (Temp_ == temp) {
649  if (Pres_ == pres) {
650  return;
651  }
652  }
653  if (m_useCanteraCalls) {
654  TP_ptr->setElectricPotential(m_phi);
655  TP_ptr->setState_TP(temp, pres);
656  } else {
657  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
658  "deprecated and will not work after Cantera 2.2.");
659  }
660  Temp_ = temp;
661  Pres_ = pres;
662  m_UpToDate_AC = false;
663  m_UpToDate_VolStar = false;
664  m_UpToDate_VolPM = false;
665  m_UpToDate_GStar = false;
666  m_UpToDate_G0 = false;
667 }
668 
669 void vcs_VolPhase::setState_T(const double temp)
670 {
671  setState_TP(temp, Pres_);
672 }
673 
674 void vcs_VolPhase::_updateVolStar() const
675 {
676  if (m_useCanteraCalls) {
677  TP_ptr->getStandardVolumes(VCS_DATA_PTR(StarMolarVol));
678  } else {
679  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
680  "deprecated and will not work after Cantera 2.2.");
681  for (size_t k = 0; k < m_numSpecies; k++) {
682  size_t kglob = IndSpecies[k];
683  vcs_SpeciesProperties* sProp = ListSpeciesPtr[k];
684  VCS_SPECIES_THERMO* sTherm = sProp->SpeciesThermo;
685  StarMolarVol[k] = (sTherm->VolStar_calc(kglob, Temp_, Pres_));
686  }
687  }
688  m_UpToDate_VolStar = true;
689 }
690 
691 double vcs_VolPhase::VolStar_calc_one(size_t kspec) const
692 {
693  if (!m_UpToDate_VolStar) {
694  _updateVolStar();
695  }
696  return StarMolarVol[kspec];
697 }
698 
699 double vcs_VolPhase::_updateVolPM() const
700 {
701  if (m_useCanteraCalls) {
702  TP_ptr->getPartialMolarVolumes(VCS_DATA_PTR(PartialMolarVol));
703  } else {
704  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
705  "deprecated and will not work after Cantera 2.2.");
706  for (size_t k = 0; k < m_numSpecies; k++) {
707  size_t kglob = IndSpecies[k];
708  vcs_SpeciesProperties* sProp = ListSpeciesPtr[k];
709  VCS_SPECIES_THERMO* sTherm = sProp->SpeciesThermo;
710  StarMolarVol[k] = (sTherm->VolStar_calc(kglob, Temp_, Pres_));
711  }
712  for (size_t k = 0; k < m_numSpecies; k++) {
713  PartialMolarVol[k] = StarMolarVol[k];
714  }
715  }
716 
717  m_totalVol = 0.0;
718  for (size_t k = 0; k < m_numSpecies; k++) {
719  m_totalVol += PartialMolarVol[k] * Xmol_[k];
720  }
721  m_totalVol *= v_totalMoles;
722 
723  if (m_totalMolesInert > 0.0) {
724  if (m_gasPhase) {
725  double volI = m_totalMolesInert * GasConstant * Temp_ / Pres_;
726  m_totalVol += volI;
727  } else {
728  throw CanteraError("vcs_VolPhase::_updateVolPM", "unknown situation");
729  }
730  }
731  m_UpToDate_VolPM = true;
732  return m_totalVol;
733 }
734 
735 void vcs_VolPhase::_updateLnActCoeffJac()
736 {
737  double phaseTotalMoles = v_totalMoles;
738  if (phaseTotalMoles < 1.0E-14) {
739  phaseTotalMoles = 1.0;
740  }
741 
742  /*
743  * Evaluate the current base activity coefficients if necessary
744  */
745  if (!m_UpToDate_AC) {
746  _updateActCoeff();
747  }
748  if (!TP_ptr) {
749  return;
750  }
751  TP_ptr->getdlnActCoeffdlnN(m_numSpecies, &np_dLnActCoeffdMolNumber(0,0));
752  for (size_t j = 0; j < m_numSpecies; j++) {
753  double moles_j_base = phaseTotalMoles * Xmol_[j];
754  double* const np_lnActCoeffCol = np_dLnActCoeffdMolNumber.ptrColumn(j);
755  if (moles_j_base < 1.0E-200) {
756  moles_j_base = 1.0E-7 * moles_j_base + 1.0E-13 * phaseTotalMoles + 1.0E-150;
757  }
758  for (size_t k = 0; k < m_numSpecies; k++) {
759  np_lnActCoeffCol[k] = np_lnActCoeffCol[k] * phaseTotalMoles / moles_j_base;
760  }
761  }
762 
763  double deltaMoles_j = 0.0;
764  // Make copies of ActCoeff and Xmol_ for use in taking differences
765  std::vector<double> ActCoeff_Base(ActCoeff);
766  std::vector<double> Xmol_Base(Xmol_);
767  double TMoles_base = phaseTotalMoles;
768 
769  /*
770  * Loop over the columns species to be deltad
771  */
772  for (size_t j = 0; j < m_numSpecies; j++) {
773  /*
774  * Calculate a value for the delta moles of species j
775  * -> Note Xmol_[] and Tmoles are always positive or zero
776  * quantities.
777  */
778  double moles_j_base = phaseTotalMoles * Xmol_Base[j];
779  deltaMoles_j = 1.0E-7 * moles_j_base + 1.0E-13 * phaseTotalMoles + 1.0E-150;
780  /*
781  * Now, update the total moles in the phase and all of the
782  * mole fractions based on this.
783  */
784  phaseTotalMoles = TMoles_base + deltaMoles_j;
785  for (size_t k = 0; k < m_numSpecies; k++) {
786  Xmol_[k] = Xmol_Base[k] * TMoles_base / phaseTotalMoles;
787  }
788  Xmol_[j] = (moles_j_base + deltaMoles_j) / phaseTotalMoles;
789 
790  /*
791  * Go get new values for the activity coefficients.
792  * -> Note this calls setState_PX();
793  */
794  _updateMoleFractionDependencies();
795  _updateActCoeff();
796  /*
797  * Revert to the base case Xmol_, v_totalMoles
798  */
799  v_totalMoles = TMoles_base;
800  Xmol_ = Xmol_Base;
801  }
802  /*
803  * Go get base values for the activity coefficients.
804  * -> Note this calls setState_TPX() again;
805  * -> Just wanted to make sure that cantera is in sync
806  * with VolPhase after this call.
807  */
808  setMoleFractions(VCS_DATA_PTR(Xmol_Base));
809  _updateMoleFractionDependencies();
810  _updateActCoeff();
811 }
812 
813 void vcs_VolPhase::sendToVCS_LnActCoeffJac(Array2D& np_LnACJac_VCS)
814 {
815  /*
816  * update the Ln Act Coeff Jacobian entries with respect to the
817  * mole number of species in the phase -> we always assume that
818  * they are out of date.
819  */
820  _updateLnActCoeffJac();
821 
822  /*
823  * Now copy over the values
824  */
825  for (size_t j = 0; j < m_numSpecies; j++) {
826  size_t jglob = IndSpecies[j];
827  for (size_t k = 0; k < m_numSpecies; k++) {
828  size_t kglob = IndSpecies[k];
829  np_LnACJac_VCS(kglob,jglob) = np_dLnActCoeffdMolNumber(k,j);
830  }
831  }
832 }
833 
834 void vcs_VolPhase::setPtrThermoPhase(ThermoPhase* tp_ptr)
835 {
836  TP_ptr = tp_ptr;
837  if (TP_ptr) {
838  m_useCanteraCalls = true;
839  Temp_ = TP_ptr->temperature();
840  Pres_ = TP_ptr->pressure();
841  setState_TP(Temp_, Pres_);
842  p_VCS_UnitsFormat = VCS_UNITS_MKS;
843  m_phi = TP_ptr->electricPotential();
844  size_t nsp = TP_ptr->nSpecies();
845  size_t nelem = TP_ptr->nElements();
846  if (nsp != m_numSpecies) {
847  if (m_numSpecies != 0) {
848  plogf("Warning Nsp != NVolSpeces: %d %d \n", nsp, m_numSpecies);
849  }
850  resize(VP_ID_, nsp, nelem, PhaseName.c_str());
851  }
852  TP_ptr->getMoleFractions(VCS_DATA_PTR(Xmol_));
853  creationMoleNumbers_ = Xmol_;
854  _updateMoleFractionDependencies();
855 
856  /*
857  * figure out ideal solution tag
858  */
859  if (nsp == 1) {
860  m_isIdealSoln = true;
861  } else {
862  int eos = TP_ptr->eosType();
863  switch (eos) {
864  case cIdealGas:
865  case cIncompressible:
866  case cSurf:
867  case cMetal:
868  case cStoichSubstance:
869  case cSemiconductor:
870  case cLatticeSolid:
871  case cLattice:
872  case cEdge:
874  m_isIdealSoln = true;
875  break;
876  default:
877  m_isIdealSoln = false;
878  };
879  }
880  } else {
881  m_useCanteraCalls = false;
882  warn_deprecated("m_useCanteraCalls", "Setting this flag to 'false' is "
883  "deprecated and will not work after Cantera 2.2.");
884 
885  }
886 }
887 
888 const ThermoPhase* vcs_VolPhase::ptrThermoPhase() const
889 {
890  return TP_ptr;
891 }
892 
893 double vcs_VolPhase::totalMoles() const
894 {
895  return v_totalMoles;
896 }
897 
898 double vcs_VolPhase::molefraction(size_t k) const
899 {
900  return Xmol_[k];
901 }
902 
903 void vcs_VolPhase::setCreationMoleNumbers(const double* const n_k,
904  const std::vector<size_t> &creationGlobalRxnNumbers)
905 {
906  creationMoleNumbers_.assign(n_k, n_k+m_numSpecies);
907  for (size_t k = 0; k < m_numSpecies; k++) {
908  creationGlobalRxnNumbers_[k] = creationGlobalRxnNumbers[k];
909  }
910 }
911 
912 const std::vector<double> & vcs_VolPhase::creationMoleNumbers(std::vector<size_t> &creationGlobalRxnNumbers) const
913 {
914  creationGlobalRxnNumbers = creationGlobalRxnNumbers_;
915  return creationMoleNumbers_;
916 }
917 
918 void vcs_VolPhase::setTotalMoles(const double totalMols)
919 {
920  v_totalMoles = totalMols;
921  if (m_totalMolesInert > 0.0) {
922  m_existence = VCS_PHASE_EXIST_ALWAYS;
923  AssertThrowMsg(totalMols >= m_totalMolesInert,
924  "vcs_VolPhase::setTotalMoles",
925  "totalMoles less than inert moles: " +
926  fp2str(totalMols) + " " + fp2str(m_totalMolesInert));
927  } else {
928  if (m_singleSpecies && (m_phiVarIndex == 0)) {
929  m_existence = VCS_PHASE_EXIST_ALWAYS;
930  } else {
931  if (totalMols > 0.0) {
932  m_existence = VCS_PHASE_EXIST_YES;
933  } else {
934  m_existence = VCS_PHASE_EXIST_NO;
935  }
936  }
937  }
938 }
939 
940 void vcs_VolPhase::setMolesOutOfDate(int stateCalc)
941 {
942  m_UpToDate = false;
943  if (stateCalc != -1) {
944  m_vcsStateStatus = stateCalc;
945  }
946 }
947 
948 void vcs_VolPhase::setMolesCurrent(int stateCalc)
949 {
950  m_UpToDate = true;
951  m_vcsStateStatus = stateCalc;
952 }
953 
954 std::string string16_EOSType(int EOSType)
955 {
956  char st[32];
957  st[16] = '\0';
958  switch (EOSType) {
959  case VCS_EOS_CONSTANT:
960  sprintf(st,"Constant ");
961  break;
962  case VCS_EOS_IDEAL_GAS:
963  sprintf(st,"Ideal Gas ");
964  break;
965  case VCS_EOS_STOICH_SUB:
966  sprintf(st,"Stoich Sub ");
967  break;
968  case VCS_EOS_IDEAL_SOLN:
969  sprintf(st,"Ideal Soln ");
970  break;
971  case VCS_EOS_DEBEYE_HUCKEL:
972  sprintf(st,"Debeye Huckel ");
973  break;
974  case VCS_EOS_REDLICK_KWONG:
975  sprintf(st,"Redlick_Kwong ");
976  break;
977  case VCS_EOS_REGULAR_SOLN:
978  sprintf(st,"Regular Soln ");
979  break;
980  default:
981  sprintf(st,"UnkType: %-7d", EOSType);
982  break;
983  }
984  st[16] = '\0';
985  return st;
986 }
987 
988 bool vcs_VolPhase::isIdealSoln() const
989 {
990  return m_isIdealSoln;
991 }
992 
993 bool vcs_VolPhase::usingCanteraCalls() const
994 {
995  return m_useCanteraCalls;
996 }
997 
998 size_t vcs_VolPhase::phiVarIndex() const
999 {
1000  return m_phiVarIndex;
1001 }
1002 
1003 void vcs_VolPhase::setPhiVarIndex(size_t phiVarIndex)
1004 {
1005  m_phiVarIndex = phiVarIndex;
1006  m_speciesUnknownType[m_phiVarIndex] = VCS_SPECIES_TYPE_INTERFACIALVOLTAGE;
1007  if (m_singleSpecies) {
1008  if (m_phiVarIndex == 0) {
1009  m_existence = VCS_PHASE_EXIST_ALWAYS;
1010  }
1011  }
1012 }
1013 
1014 vcs_SpeciesProperties* vcs_VolPhase::speciesProperty(const size_t kindex)
1015 {
1016  return ListSpeciesPtr[kindex];
1017 }
1018 
1019 int vcs_VolPhase::exists() const
1020 {
1021  return m_existence;
1022 }
1023 
1024 void vcs_VolPhase::setExistence(const int existence)
1025 {
1026  if (existence == VCS_PHASE_EXIST_NO || existence == VCS_PHASE_EXIST_ZEROEDPHASE) {
1027  if (v_totalMoles != 0.0) {
1028  if (DEBUG_MODE_ENABLED) {
1029  throw CanteraError("vcs_VolPhase::setExistence",
1030  "setting false existence for phase with moles");
1031  } else {
1032  v_totalMoles = 0.0;
1033  }
1034  }
1035  } else if (DEBUG_MODE_ENABLED && m_totalMolesInert == 0.0) {
1036  if (v_totalMoles == 0.0) {
1037  if (!m_singleSpecies || m_phiVarIndex != 0) {
1038  throw CanteraError("vcs_VolPhase::setExistence",
1039  "setting true existence for phase with no moles");
1040  }
1041  }
1042  }
1043  if (DEBUG_MODE_ENABLED && m_singleSpecies) {
1044  if (m_phiVarIndex == 0) {
1045  if (existence == VCS_PHASE_EXIST_NO || existence == VCS_PHASE_EXIST_ZEROEDPHASE) {
1046  throw CanteraError("vcs_VolPhase::setExistence",
1047  "Trying to set existence of an electron phase to false");
1048  }
1049  }
1050  }
1051  m_existence = existence;
1052 }
1053 
1054 size_t vcs_VolPhase::spGlobalIndexVCS(const size_t spIndex) const
1055 {
1056  return IndSpecies[spIndex];
1057 }
1058 
1059 void vcs_VolPhase::setSpGlobalIndexVCS(const size_t spIndex,
1060  const size_t spGlobalIndex)
1061 {
1062  IndSpecies[spIndex] = spGlobalIndex;
1063  if (spGlobalIndex >= m_numElemConstraints) {
1064  creationGlobalRxnNumbers_[spIndex] = spGlobalIndex - m_numElemConstraints;
1065  }
1066 }
1067 
1068 void vcs_VolPhase::setTotalMolesInert(const double tMolesInert)
1069 {
1070  if (m_totalMolesInert != tMolesInert) {
1071  m_UpToDate = false;
1072  m_UpToDate_AC = false;
1073  m_UpToDate_VolStar = false;
1074  m_UpToDate_VolPM = false;
1075  m_UpToDate_GStar = false;
1076  m_UpToDate_G0 = false;
1077  v_totalMoles += (tMolesInert - m_totalMolesInert);
1078  m_totalMolesInert = tMolesInert;
1079  }
1080  if (m_totalMolesInert > 0.0) {
1081  m_existence = VCS_PHASE_EXIST_ALWAYS;
1082  } else if (m_singleSpecies && (m_phiVarIndex == 0)) {
1083  m_existence = VCS_PHASE_EXIST_ALWAYS;
1084  } else {
1085  if (v_totalMoles > 0.0) {
1086  m_existence = VCS_PHASE_EXIST_YES;
1087  } else {
1088  m_existence = VCS_PHASE_EXIST_NO;
1089  }
1090  }
1091 }
1092 
1093 double vcs_VolPhase::totalMolesInert() const
1094 {
1095  return m_totalMolesInert;
1096 }
1097 
1098 size_t vcs_VolPhase::elemGlobalIndex(const size_t e) const
1099 {
1100  AssertThrow(e < m_numElemConstraints, " vcs_VolPhase::elemGlobalIndex");
1101  return m_elemGlobalIndex[e];
1102 }
1103 
1104 void vcs_VolPhase::setElemGlobalIndex(const size_t eLocal, const size_t eGlobal)
1105 {
1106  AssertThrow(eLocal < m_numElemConstraints,
1107  "vcs_VolPhase::setElemGlobalIndex");
1108  m_elemGlobalIndex[eLocal] = eGlobal;
1109 }
1110 
1111 size_t vcs_VolPhase::nElemConstraints() const
1112 {
1113  return m_numElemConstraints;
1114 }
1115 
1116 std::string vcs_VolPhase::elementName(const size_t e) const
1117 {
1118  return m_elementNames[e];
1119 }
1120 
1121 //! This function decides whether a phase has charged species or not.
1122 static bool hasChargedSpecies(const ThermoPhase* const tPhase)
1123 {
1124  for (size_t k = 0; k < tPhase->nSpecies(); k++) {
1125  if (tPhase->charge(k) != 0.0) {
1126  return true;
1127  }
1128  }
1129  return false;
1130 }
1131 
1132 /*!
1133  * This utility routine decides whether a Cantera ThermoPhase needs
1134  * a constraint equation representing the charge neutrality of the
1135  * phase. It does this by searching for charged species. If it
1136  * finds one, and if the phase needs one, then it returns true.
1137  */
1138 static bool chargeNeutralityElement(const ThermoPhase* const tPhase)
1139 {
1140  int hasCharge = hasChargedSpecies(tPhase);
1141  if (tPhase->chargeNeutralityNecessary()) {
1142  if (hasCharge) {
1143  return true;
1144  }
1145  }
1146  return false;
1147 }
1148 
1149 size_t vcs_VolPhase::transferElementsFM(const ThermoPhase* const tPhase)
1150 {
1151  size_t nebase = tPhase->nElements();
1152  size_t ne = nebase;
1153  size_t ns = tPhase->nSpecies();
1154 
1155  /*
1156  * Decide whether we need an extra element constraint for charge
1157  * neutrality of the phase
1158  */
1159  bool cne = chargeNeutralityElement(tPhase);
1160  if (cne) {
1161  ChargeNeutralityElement = ne;
1162  ne++;
1163  }
1164 
1165  /*
1166  * Assign and malloc structures
1167  */
1168  elemResize(ne);
1169 
1170 
1171  if (ChargeNeutralityElement != npos) {
1172  m_elementType[ChargeNeutralityElement] = VCS_ELEM_TYPE_CHARGENEUTRALITY;
1173  }
1174 
1175  size_t eFound = npos;
1176  if (hasChargedSpecies(tPhase)) {
1177  if (cne) {
1178  /*
1179  * We need a charge neutrality constraint.
1180  * We also have an Electron Element. These are
1181  * duplicates of each other. To avoid trouble with
1182  * possible range error conflicts, sometimes we eliminate
1183  * the Electron condition. Flag that condition for elimination
1184  * by toggling the ElActive variable. If we find we need it
1185  * later, we will retoggle ElActive to true.
1186  */
1187  for (size_t eT = 0; eT < nebase; eT++) {
1188  if (tPhase->elementName(eT) == "E") {
1189  eFound = eT;
1190  m_elementActive[eT] = 0;
1191  m_elementType[eT] = VCS_ELEM_TYPE_ELECTRONCHARGE;
1192  }
1193  }
1194  } else {
1195  for (size_t eT = 0; eT < nebase; eT++) {
1196  if (tPhase->elementName(eT) == "E") {
1197  eFound = eT;
1198  m_elementType[eT] = VCS_ELEM_TYPE_ELECTRONCHARGE;
1199  }
1200  }
1201  }
1202  if (eFound == npos) {
1203  eFound = ne;
1204  m_elementType[ne] = VCS_ELEM_TYPE_ELECTRONCHARGE;
1205  m_elementActive[ne] = 0;
1206  std::string ename = "E";
1207  m_elementNames[ne] = ename;
1208  ne++;
1209  elemResize(ne);
1210  }
1211 
1212  }
1213 
1214  m_formulaMatrix.resize(ns, ne, 0.0);
1215 
1216  m_speciesUnknownType.resize(ns, VCS_SPECIES_TYPE_MOLNUM);
1217 
1218  elemResize(ne);
1219 
1220  size_t e = 0;
1221  for (size_t eT = 0; eT < nebase; eT++) {
1222  m_elementNames[e] = tPhase->elementName(eT);
1223  m_elementType[e] = tPhase->elementType(eT);
1224  e++;
1225  }
1226 
1227  if (cne) {
1228  std::string pname = tPhase->id();
1229  if (pname == "") {
1230  std::stringstream sss;
1231  sss << "phase" << VP_ID_;
1232  pname = sss.str();
1233  }
1234  e = ChargeNeutralityElement;
1235  m_elementNames[e] = "cn_" + pname;
1236  }
1237 
1238  for (size_t k = 0; k < ns; k++) {
1239  e = 0;
1240  for (size_t eT = 0; eT < nebase; eT++) {
1241  m_formulaMatrix(k,e) = tPhase->nAtoms(k, eT);
1242  e++;
1243  }
1244  if (eFound != npos) {
1245  m_formulaMatrix(k,eFound) = - tPhase->charge(k);
1246  }
1247  }
1248 
1249  if (cne) {
1250  for (size_t k = 0; k < ns; k++) {
1251  m_formulaMatrix(k,ChargeNeutralityElement) = tPhase->charge(k);
1252  }
1253  }
1254 
1255 
1256  /*
1257  * Here, we figure out what is the species types are
1258  * The logic isn't set in stone, and is just for a particular type
1259  * of problem that I'm solving first.
1260  */
1261  if (ns == 1) {
1262  if (tPhase->charge(0) != 0.0) {
1263  m_speciesUnknownType[0] = VCS_SPECIES_TYPE_INTERFACIALVOLTAGE;
1264  setPhiVarIndex(0);
1265  }
1266  }
1267 
1268  return ne;
1269 }
1270 
1271 int vcs_VolPhase::elementType(const size_t e) const
1272 {
1273  return m_elementType[e];
1274 }
1275 
1276 void vcs_VolPhase::setElementType(const size_t e, const int eType)
1277 {
1278  m_elementType[e] = eType;
1279 }
1280 
1281 const Array2D& vcs_VolPhase::getFormulaMatrix() const
1282 {
1283  return m_formulaMatrix;
1284 }
1285 
1286 int vcs_VolPhase::speciesUnknownType(const size_t k) const
1287 {
1288  return m_speciesUnknownType[k];
1289 }
1290 
1291 int vcs_VolPhase::elementActive(const size_t e) const
1292 {
1293  return m_elementActive[e];
1294 }
1295 
1296 size_t vcs_VolPhase::nSpecies() const
1297 {
1298  return m_numSpecies;
1299 }
1300 
1301 }
#define VCS_ELEM_TYPE_ELECTRONCHARGE
This refers to conservation of electrons.
Definition: vcs_defs.h:296
doublereal nAtoms(size_t k, size_t m) const
Number of atoms of element m in species k.
Definition: Phase.cpp:243
#define VCS_PHASE_EXIST_NO
Phase doesn't currently exist in the mixture.
Definition: vcs_defs.h:225
bool vcs_doubleEqual(double d1, double d2)
Simple routine to check whether two doubles are equal up to roundoff error.
Definition: vcs_util.cpp:185
size_t nElements() const
Number of elements.
Definition: Phase.cpp:167
virtual double GStar_R_calc(size_t kspec, double TKelvin, double pres)
This function calculates the standard state Gibbs free energy for species, kspec, at the temperature ...
#define VCS_PHASE_EXIST_ALWAYS
Always exists because it contains inerts which can't exist in any other phase.
Definition: vcs_defs.h:212
#define VCS_PHASE_EXIST_ZEROEDPHASE
Phase currently is zeroed due to a programmatic issue.
Definition: vcs_defs.h:231
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:165
#define VCS_DATA_PTR(vvv)
Points to the data in a std::vector<> object.
Definition: vcs_internal.h:18
virtual double VolStar_calc(size_t kglob, double TKelvin, double Pres)
This function calculates the standard state molar volume for species, kspec, at the temperature TKelv...
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
Definition: global.cpp:78
double vcsUtil_gasConstant(int mu_units)
Returns the value of the gas constant in the units specified by parameter.
Definition: vcs_util.cpp:76
static bool hasChargedSpecies(const ThermoPhase *const tPhase)
This function decides whether a phase has charged species or not.
static bool chargeNeutralityElement(const ThermoPhase *const tPhase)
#define VCS_PHASE_EXIST_YES
Phase is a normal phase that currently exists.
Definition: vcs_defs.h:215
#define AssertThrowMsg(expr, procedure, message)
Assertion must be true or an error is thrown.
Definition: ctexceptions.h:301
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:29
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:97
const int cSurf
A surface phase. Used by class SurfPhase.
Definition: mix_defs.h:40
Header file for the internal object that holds the vcs equilibrium problem (see Class VCS_SOLVE and E...
Header for the object representing each phase within vcs.
#define VCS_SPECIES_TYPE_MOLNUM
Unknown refers to mole number of a single species.
Definition: vcs_defs.h:343
#define VCS_ELEM_TYPE_ABSPOS
Normal element constraint consisting of positive coefficients for the formula matrix.
Definition: vcs_defs.h:290
std::string fp2str(const double x, const std::string &fmt)
Convert a double into a c++ string.
Definition: stringUtils.cpp:28
Properties of a single species.
bool chargeNeutralityNecessary() const
Returns the chargeNeutralityNecessity boolean.
Definition: ThermoPhase.h:218
const int cMetal
A metal phase.
Definition: mix_defs.h:43
#define VCS_ELEM_TYPE_CHARGENEUTRALITY
This refers to a charge neutrality of a single phase.
Definition: vcs_defs.h:302
#define AssertThrow(expr, procedure)
Assertion must be true or an error is thrown.
Definition: ctexceptions.h:283
std::string id() const
Return the string id for the phase.
Definition: Phase.cpp:147
#define VCS_STATECALC_NEW
State Calculation based on the new or tentative mole numbers.
Definition: vcs_defs.h:368
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:99
const int cEdge
An edge between two 2D surfaces.
Definition: mix_defs.h:58
std::string string16_EOSType(int EOSType)
Return a string representing the equation of state.
virtual double G0_R_calc(size_t kglob, double TKelvin)
This function calculates the standard state Gibbs free energy for species, kspec, at the temperature ...
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:265
doublereal temperature() const
Temperature (K).
Definition: Phase.h:602
#define VCS_STATECALC_OLD
State Calculation based on the old or base mole numbers.
Definition: vcs_defs.h:365
Contains declarations for string manipulation functions within Cantera.
std::string elementName(size_t m) const
Name of the element with index m.
Definition: Phase.cpp:186
#define plogf
define this Cantera function to replace printf
Definition: vcs_internal.h:24
int elementType(size_t m) const
Return the element constraint type Possible types include:
Definition: Phase.cpp:231
#define VCS_SPECIES_TYPE_INTERFACIALVOLTAGE
Unknown refers to the voltage level of a phase.
Definition: vcs_defs.h:352
#define VCS_STATECALC_TMP
State Calculation based on a temporary set of mole numbers.
Definition: vcs_defs.h:376
Header file for class ThermoPhase, the base class for phases with thermodynamic properties, and the text for the Module thermoprops (see Thermodynamic Properties and class ThermoPhase).
const int cIdealGas
Equation of state types:
Definition: mix_defs.h:37
const int cIdealSolidSolnPhase
Constant partial molar volume solution IdealSolidSolnPhase.h.
Definition: mix_defs.h:64
doublereal charge(size_t k) const
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the...
Definition: Phase.h:578