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