Cantera  3.3.0a1
Loading...
Searching...
No Matches
vcs_VolPhase.cpp
Go to the documentation of this file.
1/**
2 * @file vcs_VolPhase.cpp
3 */
4
5// This file is part of Cantera. See License.txt in the top-level directory or
6// at https://cantera.org/license.txt for license and copyright information.
7
10
13
14#include <cstdio>
15
16namespace Cantera
17{
18
19vcs_VolPhase::vcs_VolPhase(VCS_SOLVE* owningSolverObject,
20 ThermoPhase* thermoPhase, size_t phaseNum)
21 : m_owningSolverObject(owningSolverObject)
22{
23 if (!m_owningSolverObject) {
24 throw CanteraError("vcs_VolPhase::vcs_VolPhase",
25 "owningSolverObject must not be null");
26 }
27 if (!thermoPhase) {
28 throw CanteraError("vcs_VolPhase::vcs_VolPhase",
29 "thermoPhase must not be null");
30 }
31
32 VP_ID_ = phaseNum;
33 PhaseName = thermoPhase->name().empty() ?
34 fmt::format("Phase_{}", VP_ID_) : thermoPhase->name();
35
36 m_numSpecies = thermoPhase->nSpecies();
37 AssertThrowMsg(m_numSpecies > 0, "vcs_VolPhase::vcs_VolPhase",
38 "nspecies Error");
39 m_singleSpecies = (m_numSpecies == 1);
40 m_phiVarIndex = npos;
41
42 p_activityConvention = thermoPhase->activityConvention();
43
44 IndSpecies.assign(m_numSpecies, npos);
45 Xmol_.assign(m_numSpecies, 0.0);
46 creationMoleNumbers_.assign(m_numSpecies, 0.0);
47 creationGlobalRxnNumbers_.assign(m_numSpecies, npos);
48 for (size_t k = 0; k < m_numSpecies; k++) {
49 double defaultFrac = 1.0 / m_numSpecies;
50 Xmol_[k] = defaultFrac;
51 creationMoleNumbers_[k] = defaultFrac;
52 }
53
54 SS0ChemicalPotential.assign(m_numSpecies, -1.0);
55 StarChemicalPotential.assign(m_numSpecies, -1.0);
56 StarMolarVol.assign(m_numSpecies, -1.0);
57 PartialMolarVol.assign(m_numSpecies, -1.0);
58 ActCoeff.assign(m_numSpecies, 1.0);
59 np_dLnActCoeffdMolNumber.resize(m_numSpecies, m_numSpecies, 0.0);
60
61 m_speciesUnknownType.assign(m_numSpecies, VCS_SPECIES_TYPE_MOLNUM);
62 m_UpToDate = false;
63 m_vcsStateStatus = VCS_STATECALC_OLD;
64 m_UpToDate_AC = false;
65 m_UpToDate_VolStar = false;
66 m_UpToDate_VolPM = false;
67 m_UpToDate_GStar = false;
68 m_UpToDate_G0 = false;
69
70 TP_ptr = thermoPhase;
71 Temp_ = TP_ptr->temperature();
72 Pres_ = TP_ptr->pressure();
73 m_phi = TP_ptr->electricPotential();
74
75 transferElementsFM(TP_ptr);
76 setTotalMoles(0.0);
77
78 // Determine the global index for each of the elements in the phase
79 size_t neVP = nElemConstraints();
80 for (size_t eVP = 0; eVP < neVP; eVP++) {
81 string enVP = elementName(eVP);
82 size_t foundPos = m_owningSolverObject->elementIndex(enVP);
83 if (foundPos == npos) {
84 int elType = elementType(eVP);
85 int elactive = elementActive(eVP);
86 foundPos = m_owningSolverObject->addElement(enVP.c_str(), elType, elactive);
87 }
88 setElemGlobalIndex(eVP, foundPos);
89 }
90
91 setState_TP(Temp_, Pres_);
92 TP_ptr->getMoleFractions(&Xmol_[0]);
93 creationMoleNumbers_ = Xmol_;
94 _updateMoleFractionDependencies();
95
96 m_isIdealSoln = m_singleSpecies ? true : TP_ptr->isIdeal();
97}
98
99void vcs_VolPhase::elemResize(const size_t numElemConstraints)
100{
101 m_elementNames.resize(numElemConstraints);
102 m_elementActive.resize(numElemConstraints+1, 1);
103 m_elementType.resize(numElemConstraints, VCS_ELEM_TYPE_ABSPOS);
104 m_formulaMatrix.resize(m_numSpecies, numElemConstraints, 0.0);
105 m_elementNames.resize(numElemConstraints, "");
106 m_elemGlobalIndex.resize(numElemConstraints, npos);
107 m_numElemConstraints = numElemConstraints;
108}
109
110void vcs_VolPhase::_updateActCoeff() const
111{
112 if (m_isIdealSoln) {
113 m_UpToDate_AC = true;
114 return;
115 }
116 TP_ptr->getActivityCoefficients(&ActCoeff[0]);
117 m_UpToDate_AC = true;
118}
119
120void vcs_VolPhase::_updateG0() const
121{
122 TP_ptr->getGibbs_ref(&SS0ChemicalPotential[0]);
123 m_UpToDate_G0 = true;
124}
125
126double vcs_VolPhase::G0_calc_one(size_t kspec) const
127{
128 if (!m_UpToDate_G0) {
129 _updateG0();
130 }
131 return SS0ChemicalPotential[kspec];
132}
133
134void vcs_VolPhase::_updateGStar() const
135{
136 TP_ptr->getStandardChemPotentials(&StarChemicalPotential[0]);
137 m_UpToDate_GStar = true;
138}
139
140double vcs_VolPhase::GStar_calc_one(size_t kspec) const
141{
142 if (!m_UpToDate_GStar) {
143 _updateGStar();
144 }
145 return StarChemicalPotential[kspec];
146}
147
148void vcs_VolPhase::setMoleFractions(const double* const xmol)
149{
150 double sum = -1.0;
151 for (size_t k = 0; k < m_numSpecies; k++) {
152 Xmol_[k] = xmol[k];
153 sum+= xmol[k];
154 }
155 if (std::fabs(sum) > 1.0E-13) {
156 for (size_t k = 0; k < m_numSpecies; k++) {
157 Xmol_[k] /= sum;
158 }
159 }
160 _updateMoleFractionDependencies();
161 m_UpToDate = false;
162 m_vcsStateStatus = VCS_STATECALC_TMP;
163}
164
165void vcs_VolPhase::_updateMoleFractionDependencies()
166{
167 if (TP_ptr) {
168 TP_ptr->setMoleFractions(&Xmol_[m_MFStartIndex]);
169 TP_ptr->setPressure(Pres_);
170 }
171 if (!m_isIdealSoln) {
172 m_UpToDate_AC = false;
173 m_UpToDate_VolPM = false;
174 }
175}
176
177const vector<double> & vcs_VolPhase::moleFractions() const
178{
179 return Xmol_;
180}
181
182double vcs_VolPhase::moleFraction(size_t k) const
183{
184 return Xmol_[k];
185}
186
187void vcs_VolPhase::setMoleFractionsState(const double totalMoles,
188 const double* const moleFractions,
189 const int vcsStateStatus)
190{
191 if (totalMoles != 0.0) {
192 // There are other ways to set the mole fractions when VCS_STATECALC
193 // is set to a normal settting.
194 if (vcsStateStatus != VCS_STATECALC_TMP) {
195 throw CanteraError("vcs_VolPhase::setMolesFractionsState",
196 "inappropriate usage");
197 }
198 m_UpToDate = false;
199 m_vcsStateStatus = VCS_STATECALC_TMP;
200 if (m_existence == VCS_PHASE_EXIST_ZEROEDPHASE) {
201 throw CanteraError("vcs_VolPhase::setMolesFractionsState",
202 "inappropriate usage");
203 }
204 m_existence = VCS_PHASE_EXIST_YES;
205 } else {
206 m_UpToDate = true;
207 m_vcsStateStatus = vcsStateStatus;
208 m_existence = std::min(m_existence, VCS_PHASE_EXIST_NO);
209 }
210 double fractotal = 1.0;
211 v_totalMoles = totalMoles;
212 double sum = 0.0;
213 for (size_t k = 0; k < m_numSpecies; k++) {
214 Xmol_[k] = moleFractions[k];
215 sum += moleFractions[k];
216 }
217 if (sum == 0.0) {
218 throw CanteraError("vcs_VolPhase::setMolesFractionsState",
219 "inappropriate usage");
220 }
221 if (sum != fractotal) {
222 for (size_t k = 0; k < m_numSpecies; k++) {
223 Xmol_[k] *= (fractotal /sum);
224 }
225 }
226 _updateMoleFractionDependencies();
227}
228
229void vcs_VolPhase::setMolesFromVCS(const int stateCalc,
230 const double* molesSpeciesVCS)
231{
232 v_totalMoles = 0.0;
233
234 if (molesSpeciesVCS == 0) {
235 if (stateCalc == VCS_STATECALC_OLD) {
236 molesSpeciesVCS = &m_owningSolverObject->m_molNumSpecies_old[0];
237 } else if (stateCalc == VCS_STATECALC_NEW) {
238 molesSpeciesVCS = &m_owningSolverObject->m_molNumSpecies_new[0];
239 } else {
240 throw CanteraError("vcs_VolPhase::setMolesFromVCS", "shouldn't be here");
241 }
242 }
243
244 for (size_t k = 0; k < m_numSpecies; k++) {
245 if (m_speciesUnknownType[k] != VCS_SPECIES_TYPE_INTERFACIALVOLTAGE) {
246 size_t kglob = IndSpecies[k];
247 v_totalMoles += std::max(0.0, molesSpeciesVCS[kglob]);
248 }
249 }
250 if (v_totalMoles > 0.0) {
251 for (size_t k = 0; k < m_numSpecies; k++) {
252 if (m_speciesUnknownType[k] != VCS_SPECIES_TYPE_INTERFACIALVOLTAGE) {
253 size_t kglob = IndSpecies[k];
254 double tmp = std::max(0.0, molesSpeciesVCS[kglob]);
255 Xmol_[k] = tmp / v_totalMoles;
256 }
257 }
258 m_existence = VCS_PHASE_EXIST_YES;
259 } else {
260 // This is where we will start to store a better approximation
261 // for the mole fractions, when the phase doesn't exist.
262 // This is currently unimplemented.
263 m_existence = VCS_PHASE_EXIST_NO;
264 }
265
266 // Update the electric potential if it is a solution variable in the
267 // equation system
268 if (m_phiVarIndex != npos) {
269 size_t kglob = IndSpecies[m_phiVarIndex];
270 if (m_numSpecies == 1) {
271 Xmol_[m_phiVarIndex] = 1.0;
272 } else {
273 Xmol_[m_phiVarIndex] = 0.0;
274 }
275 double phi = molesSpeciesVCS[kglob];
276 setElectricPotential(phi);
277 if (m_numSpecies == 1) {
278 m_existence = VCS_PHASE_EXIST_YES;
279 }
280 }
281 _updateMoleFractionDependencies();
282
283 // If stateCalc is old and the total moles is positive, then we have a valid
284 // state. If the phase went away, it would be a valid starting point for
285 // F_k's. So, save the state.
286 if (stateCalc == VCS_STATECALC_OLD && v_totalMoles > 0.0) {
287 creationMoleNumbers_ = Xmol_;
288 }
289
290 // Set flags indicating we are up to date with the VCS state vector.
291 m_UpToDate = true;
292 m_vcsStateStatus = stateCalc;
293}
294
295void vcs_VolPhase::setMolesFromVCSCheck(const int vcsStateStatus,
296 const double* molesSpeciesVCS,
297 const double* const TPhMoles)
298{
299 setMolesFromVCS(vcsStateStatus, molesSpeciesVCS);
300
301 // Check for consistency with TPhMoles[]
302 double Tcheck = TPhMoles[VP_ID_];
303 if (Tcheck != v_totalMoles) {
304 if (vcs_doubleEqual(Tcheck, v_totalMoles)) {
305 Tcheck = v_totalMoles;
306 } else {
307 throw CanteraError("vcs_VolPhase::setMolesFromVCSCheck",
308 "We have a consistency problem: {} {}", Tcheck, v_totalMoles);
309 }
310 }
311}
312
313void vcs_VolPhase::updateFromVCS_MoleNumbers(const int vcsStateStatus)
314{
315 if ((!m_UpToDate || vcsStateStatus != m_vcsStateStatus) &&
316 (vcsStateStatus == VCS_STATECALC_OLD || vcsStateStatus == VCS_STATECALC_NEW)) {
317 setMolesFromVCS(vcsStateStatus);
318 }
319}
320
321void vcs_VolPhase::sendToVCS_ActCoeff(const int vcsStateStatus,
322 double* const AC)
323{
324 updateFromVCS_MoleNumbers(vcsStateStatus);
325 if (!m_UpToDate_AC) {
326 _updateActCoeff();
327 }
328 for (size_t k = 0; k < m_numSpecies; k++) {
329 size_t kglob = IndSpecies[k];
330 AC[kglob] = ActCoeff[k];
331 }
332}
333
334double vcs_VolPhase::sendToVCS_VolPM(double* const VolPM) const
335{
336 if (!m_UpToDate_VolPM) {
337 _updateVolPM();
338 }
339 for (size_t k = 0; k < m_numSpecies; k++) {
340 size_t kglob = IndSpecies[k];
341 VolPM[kglob] = PartialMolarVol[k];
342 }
343 return m_totalVol;
344}
345
346void vcs_VolPhase::sendToVCS_GStar(double* const gstar) const
347{
348 if (!m_UpToDate_GStar) {
349 _updateGStar();
350 }
351 for (size_t k = 0; k < m_numSpecies; k++) {
352 size_t kglob = IndSpecies[k];
353 gstar[kglob] = StarChemicalPotential[k];
354 }
355}
356
357void vcs_VolPhase::setElectricPotential(const double phi)
358{
359 m_phi = phi;
360 TP_ptr->setElectricPotential(m_phi);
361 // We have changed the state variable. Set uptodate flags to false
362 m_UpToDate_AC = false;
363 m_UpToDate_VolStar = false;
364 m_UpToDate_VolPM = false;
365 m_UpToDate_GStar = false;
366}
367
368double vcs_VolPhase::electricPotential() const
369{
370 return m_phi;
371}
372
373void vcs_VolPhase::setState_TP(const double temp, const double pres)
374{
375 if (Temp_ == temp && Pres_ == pres) {
376 return;
377 }
378 TP_ptr->setElectricPotential(m_phi);
379 TP_ptr->setState_TP(temp, pres);
380 Temp_ = temp;
381 Pres_ = pres;
382 m_UpToDate_AC = false;
383 m_UpToDate_VolStar = false;
384 m_UpToDate_VolPM = false;
385 m_UpToDate_GStar = false;
386 m_UpToDate_G0 = false;
387}
388
389void vcs_VolPhase::_updateVolStar() const
390{
391 TP_ptr->getStandardVolumes(&StarMolarVol[0]);
392 m_UpToDate_VolStar = true;
393}
394
395double vcs_VolPhase::_updateVolPM() const
396{
397 TP_ptr->getPartialMolarVolumes(&PartialMolarVol[0]);
398 m_totalVol = 0.0;
399 for (size_t k = 0; k < m_numSpecies; k++) {
400 m_totalVol += PartialMolarVol[k] * Xmol_[k];
401 }
402 m_totalVol *= v_totalMoles;
403 m_UpToDate_VolPM = true;
404 return m_totalVol;
405}
406
407void vcs_VolPhase::_updateLnActCoeffJac()
408{
409 double phaseTotalMoles = v_totalMoles;
410 if (phaseTotalMoles < 1.0E-14) {
411 phaseTotalMoles = 1.0;
412 }
413
414 // Evaluate the current base activity coefficients if necessary
415 if (!m_UpToDate_AC) {
416 _updateActCoeff();
417 }
418 if (!TP_ptr) {
419 return;
420 }
421 TP_ptr->getdlnActCoeffdlnN(m_numSpecies, &np_dLnActCoeffdMolNumber(0,0));
422 for (size_t j = 0; j < m_numSpecies; j++) {
423 double moles_j_base = phaseTotalMoles * Xmol_[j];
424 double* const np_lnActCoeffCol = np_dLnActCoeffdMolNumber.ptrColumn(j);
425 if (moles_j_base < 1.0E-200) {
426 moles_j_base = 1.0E-7 * moles_j_base + 1.0E-13 * phaseTotalMoles + 1.0E-150;
427 }
428 for (size_t k = 0; k < m_numSpecies; k++) {
429 np_lnActCoeffCol[k] = np_lnActCoeffCol[k] * phaseTotalMoles / moles_j_base;
430 }
431 }
432
433 double deltaMoles_j = 0.0;
434 // Make copies of ActCoeff and Xmol_ for use in taking differences
435 vector<double> ActCoeff_Base(ActCoeff);
436 vector<double> Xmol_Base(Xmol_);
437 double TMoles_base = phaseTotalMoles;
438
439 // Loop over the columns species to be deltad
440 for (size_t j = 0; j < m_numSpecies; j++) {
441 // Calculate a value for the delta moles of species j. Note Xmol_[] and
442 // Tmoles are always positive or zero quantities.
443 double moles_j_base = phaseTotalMoles * Xmol_Base[j];
444 deltaMoles_j = 1.0E-7 * moles_j_base + 1.0E-13 * phaseTotalMoles + 1.0E-150;
445
446 // Now, update the total moles in the phase and all of the mole
447 // fractions based on this.
448 phaseTotalMoles = TMoles_base + deltaMoles_j;
449 for (size_t k = 0; k < m_numSpecies; k++) {
450 Xmol_[k] = Xmol_Base[k] * TMoles_base / phaseTotalMoles;
451 }
452 Xmol_[j] = (moles_j_base + deltaMoles_j) / phaseTotalMoles;
453
454 // Go get new values for the activity coefficients. Note this calls
455 // setState_PX();
456 _updateMoleFractionDependencies();
457 _updateActCoeff();
458
459 // Revert to the base case Xmol_, v_totalMoles
460 v_totalMoles = TMoles_base;
461 Xmol_ = Xmol_Base;
462 }
463
464 // Go get base values for the activity coefficients. Note this calls
465 // setState_TPX() again; Just wanted to make sure that cantera is in sync
466 // with VolPhase after this call.
467 setMoleFractions(&Xmol_Base[0]);
468 _updateMoleFractionDependencies();
469 _updateActCoeff();
470}
471
472void vcs_VolPhase::sendToVCS_LnActCoeffJac(Array2D& np_LnACJac_VCS)
473{
474 // update the Ln Act Coeff Jacobian entries with respect to the mole number
475 // of species in the phase -> we always assume that they are out of date.
476 _updateLnActCoeffJac();
477
478 // Now copy over the values
479 for (size_t j = 0; j < m_numSpecies; j++) {
480 size_t jglob = IndSpecies[j];
481 for (size_t k = 0; k < m_numSpecies; k++) {
482 size_t kglob = IndSpecies[k];
483 np_LnACJac_VCS(kglob,jglob) = np_dLnActCoeffdMolNumber(k,j);
484 }
485 }
486}
487
488double vcs_VolPhase::totalMoles() const
489{
490 return v_totalMoles;
491}
492
493double vcs_VolPhase::molefraction(size_t k) const
494{
495 return Xmol_[k];
496}
497
498void vcs_VolPhase::setCreationMoleNumbers(const double* const n_k,
499 const vector<size_t> &creationGlobalRxnNumbers)
500{
501 creationMoleNumbers_.assign(n_k, n_k+m_numSpecies);
502 for (size_t k = 0; k < m_numSpecies; k++) {
503 creationGlobalRxnNumbers_[k] = creationGlobalRxnNumbers[k];
504 }
505}
506
507const vector<double>& vcs_VolPhase::creationMoleNumbers(
508 vector<size_t> &creationGlobalRxnNumbers) const
509{
510 creationGlobalRxnNumbers = creationGlobalRxnNumbers_;
511 return creationMoleNumbers_;
512}
513
514void vcs_VolPhase::setTotalMoles(const double totalMols)
515{
516 v_totalMoles = totalMols;
517 if (m_singleSpecies && (m_phiVarIndex == 0)) {
518 m_existence = VCS_PHASE_EXIST_ALWAYS;
519 } else {
520 if (totalMols > 0.0) {
521 m_existence = VCS_PHASE_EXIST_YES;
522 } else {
523 m_existence = VCS_PHASE_EXIST_NO;
524 }
525 }
526}
527
528void vcs_VolPhase::setMolesOutOfDate(int stateCalc)
529{
530 m_UpToDate = false;
531 if (stateCalc != -1) {
532 m_vcsStateStatus = stateCalc;
533 }
534}
535
536void vcs_VolPhase::setMolesCurrent(int stateCalc)
537{
538 m_UpToDate = true;
539 m_vcsStateStatus = stateCalc;
540}
541
542bool vcs_VolPhase::isIdealSoln() const
543{
544 return m_isIdealSoln;
545}
546
547size_t vcs_VolPhase::phiVarIndex() const
548{
549 return m_phiVarIndex;
550}
551
552void vcs_VolPhase::setPhiVarIndex(size_t phiVarIndex)
553{
554 m_phiVarIndex = phiVarIndex;
555 m_speciesUnknownType[m_phiVarIndex] = VCS_SPECIES_TYPE_INTERFACIALVOLTAGE;
556 if (m_singleSpecies && m_phiVarIndex == 0) {
557 m_existence = VCS_PHASE_EXIST_ALWAYS;
558 }
559}
560
561int vcs_VolPhase::exists() const
562{
563 return m_existence;
564}
565
566void vcs_VolPhase::setExistence(const int existence)
567{
568 if (existence == VCS_PHASE_EXIST_NO || existence == VCS_PHASE_EXIST_ZEROEDPHASE) {
569 if (v_totalMoles != 0.0) {
570 throw CanteraError("vcs_VolPhase::setExistence",
571 "setting false existence for phase with moles");
572 }
573 } else if (v_totalMoles == 0.0 && (!m_singleSpecies || m_phiVarIndex != 0)) {
574 throw CanteraError("vcs_VolPhase::setExistence",
575 "setting true existence for phase with no moles");
576 }
577 if (m_singleSpecies && m_phiVarIndex == 0 && (existence == VCS_PHASE_EXIST_NO || existence == VCS_PHASE_EXIST_ZEROEDPHASE)) {
578 throw CanteraError("vcs_VolPhase::setExistence",
579 "Trying to set existence of an electron phase to false");
580 }
581 m_existence = existence;
582}
583
584size_t vcs_VolPhase::spGlobalIndexVCS(const size_t spIndex) const
585{
586 return IndSpecies[spIndex];
587}
588
589void vcs_VolPhase::setSpGlobalIndexVCS(const size_t spIndex,
590 const size_t spGlobalIndex)
591{
592 IndSpecies[spIndex] = spGlobalIndex;
593 if (spGlobalIndex >= m_numElemConstraints) {
594 creationGlobalRxnNumbers_[spIndex] = spGlobalIndex - m_numElemConstraints;
595 }
596}
597
598size_t vcs_VolPhase::elemGlobalIndex(const size_t e) const
599{
600 AssertThrow(e < m_numElemConstraints, " vcs_VolPhase::elemGlobalIndex");
601 return m_elemGlobalIndex[e];
602}
603
604void vcs_VolPhase::setElemGlobalIndex(const size_t eLocal, const size_t eGlobal)
605{
606 AssertThrow(eLocal < m_numElemConstraints,
607 "vcs_VolPhase::setElemGlobalIndex");
608 m_elemGlobalIndex[eLocal] = eGlobal;
609}
610
611size_t vcs_VolPhase::nElemConstraints() const
612{
613 return m_numElemConstraints;
614}
615
616string vcs_VolPhase::elementName(size_t e) const
617{
618
619 if (e < m_elementNames.size()) {
620 return m_elementNames[e];
621 }
622 throw IndexError("vcs_VolPhase::elementName", "element", e, m_elementNames.size());
623}
624
625//! This function decides whether a phase has charged species or not.
626static bool hasChargedSpecies(const ThermoPhase* const tPhase)
627{
628 for (size_t k = 0; k < tPhase->nSpecies(); k++) {
629 if (tPhase->charge(k) != 0.0) {
630 return true;
631 }
632 }
633 return false;
634}
635
636//! This utility routine decides whether a Cantera ThermoPhase needs
637//! a constraint equation representing the charge neutrality of the
638//! phase. It does this by searching for charged species. If it
639//! finds one, and if the phase needs one, then it returns true.
640static bool chargeNeutralityElement(const ThermoPhase* const tPhase)
641{
642 int hasCharge = hasChargedSpecies(tPhase);
643 if (tPhase->chargeNeutralityNecessary() && hasCharge) {
644 return true;
645 }
646 return false;
647}
648
649size_t vcs_VolPhase::transferElementsFM(const ThermoPhase* const tPhase)
650{
651 size_t nebase = tPhase->nElements();
652 size_t ne = nebase;
653 size_t ns = tPhase->nSpecies();
654
655 // Decide whether we need an extra element constraint for charge
656 // neutrality of the phase
657 bool cne = chargeNeutralityElement(tPhase);
658 if (cne) {
659 ChargeNeutralityElement = ne;
660 ne++;
661 }
662
663 // Assign and resize structures
664 elemResize(ne);
665
666 if (ChargeNeutralityElement != npos) {
667 m_elementType[ChargeNeutralityElement] = VCS_ELEM_TYPE_CHARGENEUTRALITY;
668 }
669
670 size_t eFound = npos;
671 if (hasChargedSpecies(tPhase)) {
672 if (cne) {
673 // We need a charge neutrality constraint. We also have an Electron
674 // Element. These are duplicates of each other. To avoid trouble
675 // with possible range error conflicts, sometimes we eliminate the
676 // Electron condition. Flag that condition for elimination by
677 // toggling the ElActive variable. If we find we need it later, we
678 // will retoggle ElActive to true.
679 for (size_t eT = 0; eT < nebase; eT++) {
680 if (tPhase->elementName(eT) == "E") {
681 eFound = eT;
682 m_elementActive[eT] = 0;
683 m_elementType[eT] = VCS_ELEM_TYPE_ELECTRONCHARGE;
684 }
685 }
686 } else {
687 for (size_t eT = 0; eT < nebase; eT++) {
688 if (tPhase->elementName(eT) == "E") {
689 eFound = eT;
690 m_elementType[eT] = VCS_ELEM_TYPE_ELECTRONCHARGE;
691 }
692 }
693 }
694 if (eFound == npos) {
695 eFound = ne;
696 m_elementType[ne] = VCS_ELEM_TYPE_ELECTRONCHARGE;
697 m_elementActive[ne] = 0;
698 string ename = "E";
699 m_elementNames[ne] = ename;
700 ne++;
701 elemResize(ne);
702 }
703 }
704
705 m_formulaMatrix.resize(ns, ne, 0.0);
706 m_speciesUnknownType.resize(ns, VCS_SPECIES_TYPE_MOLNUM);
707 elemResize(ne);
708
709 size_t e = 0;
710 for (size_t eT = 0; eT < nebase; eT++) {
711 m_elementNames[e] = tPhase->elementName(eT);
712 m_elementType[e] = tPhase->elementType(eT);
713 e++;
714 }
715
716 if (cne) {
717 string pname = tPhase->name();
718 if (pname == "") {
719 pname = fmt::format("phase{}", VP_ID_);
720 }
721 e = ChargeNeutralityElement;
722 m_elementNames[e] = "cn_" + pname;
723 }
724
725 for (size_t k = 0; k < ns; k++) {
726 e = 0;
727 for (size_t eT = 0; eT < nebase; eT++) {
728 m_formulaMatrix(k,e) = tPhase->nAtoms(k, eT);
729 e++;
730 }
731 if (eFound != npos) {
732 m_formulaMatrix(k,eFound) = - tPhase->charge(k);
733 }
734 }
735
736 if (cne) {
737 for (size_t k = 0; k < ns; k++) {
738 m_formulaMatrix(k,ChargeNeutralityElement) = tPhase->charge(k);
739 }
740 }
741
742 // Here, we figure out what is the species types are The logic isn't set in
743 // stone, and is just for a particular type of problem that I'm solving
744 // first.
745 if (ns == 1 && tPhase->charge(0) != 0.0) {
746 m_speciesUnknownType[0] = VCS_SPECIES_TYPE_INTERFACIALVOLTAGE;
747 setPhiVarIndex(0);
748 }
749
750 return ne;
751}
752
753int vcs_VolPhase::elementType(const size_t e) const
754{
755 return m_elementType[e];
756}
757
758const Array2D& vcs_VolPhase::getFormulaMatrix() const
759{
760 return m_formulaMatrix;
761}
762
763int vcs_VolPhase::speciesUnknownType(const size_t k) const
764{
765 return m_speciesUnknownType[k];
766}
767
768int vcs_VolPhase::elementActive(const size_t e) const
769{
770 return m_elementActive[e];
771}
772
773size_t vcs_VolPhase::nSpecies() const
774{
775 return m_numSpecies;
776}
777
778string vcs_VolPhase::eos_name() const
779{
780 return TP_ptr->type();
781}
782}
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition Array.h:32
Base class for exceptions thrown by Cantera classes.
An array index is out of range.
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:245
int elementType(size_t m) const
Return the element constraint type Possible types include:
Definition Phase.cpp:89
double nAtoms(size_t k, size_t m) const
Number of atoms of element m in species k.
Definition Phase.cpp:101
size_t nElements() const
Number of elements.
Definition Phase.cpp:30
string elementName(size_t m) const
Name of the element with index m.
Definition Phase.cpp:43
double 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:574
string name() const
Return the name of the phase.
Definition Phase.cpp:20
Base class for a phase with thermodynamic properties.
bool chargeNeutralityNecessary() const
Returns the chargeNeutralityNecessity boolean.
#define AssertThrow(expr, procedure)
Assertion must be true or an error is thrown.
#define AssertThrowMsg(expr, procedure,...)
Assertion must be true or an error is thrown.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
static bool chargeNeutralityElement(const ThermoPhase *const tPhase)
This utility routine decides whether a Cantera ThermoPhase needs a constraint equation representing t...
static bool hasChargedSpecies(const ThermoPhase *const tPhase)
This function decides whether a phase has charged species or not.
bool vcs_doubleEqual(double d1, double d2)
Simple routine to check whether two doubles are equal up to roundoff error.
Definition vcs_util.cpp:89
Contains declarations for string manipulation functions within Cantera.
Header for the object representing each phase within vcs.
#define VCS_SPECIES_TYPE_INTERFACIALVOLTAGE
Unknown refers to the voltage level of a phase.
Definition vcs_defs.h:197
#define VCS_STATECALC_OLD
State Calculation based on the old or base mole numbers.
Definition vcs_defs.h:206
#define VCS_STATECALC_TMP
State Calculation based on a temporary set of mole numbers.
Definition vcs_defs.h:217
#define VCS_PHASE_EXIST_NO
Phase doesn't currently exist in the mixture.
Definition vcs_defs.h:147
#define VCS_PHASE_EXIST_ALWAYS
These defines are valid values for the phase existence flag.
Definition vcs_defs.h:141
#define VCS_ELEM_TYPE_ABSPOS
Normal element constraint consisting of positive coefficients for the formula matrix.
Definition vcs_defs.h:170
#define VCS_STATECALC_NEW
State Calculation based on the new or tentative mole numbers.
Definition vcs_defs.h:209
#define VCS_SPECIES_TYPE_MOLNUM
Unknown refers to mole number of a single species.
Definition vcs_defs.h:189
#define VCS_ELEM_TYPE_ELECTRONCHARGE
This refers to conservation of electrons.
Definition vcs_defs.h:176
#define VCS_ELEM_TYPE_CHARGENEUTRALITY
This refers to a charge neutrality of a single phase.
Definition vcs_defs.h:182
#define VCS_PHASE_EXIST_YES
Phase is a normal phase that currently exists.
Definition vcs_defs.h:144
#define VCS_PHASE_EXIST_ZEROEDPHASE
Phase currently is zeroed due to a programmatic issue.
Definition vcs_defs.h:153
Header file for the internal object that holds the vcs equilibrium problem (see Class VCS_SOLVE and C...