Cantera  3.2.0
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
11
14
15#include <cstdio>
16
17namespace Cantera
18{
19
20vcs_VolPhase::vcs_VolPhase(VCS_SOLVE* owningSolverObject) :
21 m_owningSolverObject(owningSolverObject)
22{
23}
24
25vcs_VolPhase::~vcs_VolPhase()
26{
27 for (size_t k = 0; k < m_numSpecies; k++) {
28 delete ListSpeciesPtr[k];
29 }
30}
31
32void vcs_VolPhase::resize(const size_t phaseNum, const size_t nspecies,
33 const size_t numElem, const char* const phaseName,
34 const double molesInert)
35{
36 AssertThrowMsg(nspecies > 0, "vcs_VolPhase::resize", "nspecies Error");
37 setTotalMolesInert(molesInert);
38 m_phi = 0.0;
40
41 if (phaseNum == VP_ID_) {
42 if (strcmp(PhaseName.c_str(), phaseName)) {
43 throw CanteraError("vcs_VolPhase::resize",
44 "Strings are different: " + PhaseName + " " +
45 phaseName + " :unknown situation");
46 }
47 } else {
48 VP_ID_ = phaseNum;
49 if (!phaseName) {
50 PhaseName = fmt::format("Phase_{}", VP_ID_);
51 } else {
52 PhaseName = phaseName;
53 }
54 }
55 if (nspecies > 1) {
56 m_singleSpecies = false;
57 } else {
58 m_singleSpecies = true;
59 }
60
61 if (m_numSpecies == nspecies && numElem == m_numElemConstraints) {
62 return;
63 }
64
65 m_numSpecies = nspecies;
66 if (nspecies > 1) {
67 m_singleSpecies = false;
68 }
69
70 IndSpecies.resize(nspecies, npos);
71
72 if (ListSpeciesPtr.size() >= m_numSpecies) {
73 for (size_t i = 0; i < m_numSpecies; i++) {
74 if (ListSpeciesPtr[i]) {
75 delete ListSpeciesPtr[i];
76 ListSpeciesPtr[i] = 0;
77 }
78 }
79 }
80 ListSpeciesPtr.resize(nspecies, 0);
81 for (size_t i = 0; i < nspecies; i++) {
82 ListSpeciesPtr[i] = new vcs_SpeciesProperties(phaseNum, i, this);
83 }
84
85 Xmol_.resize(nspecies, 0.0);
86 creationMoleNumbers_.resize(nspecies, 0.0);
87 creationGlobalRxnNumbers_.resize(nspecies, npos);
88 for (size_t i = 0; i < nspecies; i++) {
89 Xmol_[i] = 1.0/nspecies;
90 creationMoleNumbers_[i] = 1.0/nspecies;
93 } else {
95 }
96 }
97
98 SS0ChemicalPotential.resize(nspecies, -1.0);
99 StarChemicalPotential.resize(nspecies, -1.0);
100 StarMolarVol.resize(nspecies, -1.0);
101 PartialMolarVol.resize(nspecies, -1.0);
102 ActCoeff.resize(nspecies, 1.0);
103 np_dLnActCoeffdMolNumber.resize(nspecies, nspecies, 0.0);
104
106 m_UpToDate = false;
108 m_UpToDate_AC = false;
109 m_UpToDate_VolStar = false;
110 m_UpToDate_VolPM = false;
111 m_UpToDate_GStar = false;
112 m_UpToDate_G0 = false;
113
114 elemResize(numElem);
115}
116
117void vcs_VolPhase::elemResize(const size_t numElemConstraints)
118{
119 m_elementNames.resize(numElemConstraints);
120 m_elementActive.resize(numElemConstraints+1, 1);
121 m_elementType.resize(numElemConstraints, VCS_ELEM_TYPE_ABSPOS);
122 m_formulaMatrix.resize(m_numSpecies, numElemConstraints, 0.0);
123 m_elementNames.resize(numElemConstraints, "");
124 m_elemGlobalIndex.resize(numElemConstraints, npos);
125 m_numElemConstraints = numElemConstraints;
126}
127
129{
130 if (m_isIdealSoln) {
131 m_UpToDate_AC = true;
132 return;
133 }
134 TP_ptr->getActivityCoefficients(&ActCoeff[0]);
135 m_UpToDate_AC = true;
136}
137
139{
140 TP_ptr->getGibbs_ref(&SS0ChemicalPotential[0]);
141 m_UpToDate_G0 = true;
142}
143
144double vcs_VolPhase::G0_calc_one(size_t kspec) const
145{
146 if (!m_UpToDate_G0) {
147 _updateG0();
148 }
149 return SS0ChemicalPotential[kspec];
150}
151
153{
154 TP_ptr->getStandardChemPotentials(&StarChemicalPotential[0]);
155 m_UpToDate_GStar = true;
156}
157
158double vcs_VolPhase::GStar_calc_one(size_t kspec) const
159{
160 if (!m_UpToDate_GStar) {
161 _updateGStar();
162 }
163 return StarChemicalPotential[kspec];
164}
165
166void vcs_VolPhase::setMoleFractions(const double* const xmol)
167{
168 double sum = -1.0;
169 for (size_t k = 0; k < m_numSpecies; k++) {
170 Xmol_[k] = xmol[k];
171 sum+= xmol[k];
172 }
173 if (std::fabs(sum) > 1.0E-13) {
174 for (size_t k = 0; k < m_numSpecies; k++) {
175 Xmol_[k] /= sum;
176 }
177 }
179 m_UpToDate = false;
181}
182
184{
185 if (TP_ptr) {
186 TP_ptr->setMoleFractions(&Xmol_[m_MFStartIndex]);
187 TP_ptr->setPressure(Pres_);
188 }
189 if (!m_isIdealSoln) {
190 m_UpToDate_AC = false;
191 m_UpToDate_VolPM = false;
192 }
193}
194
195const vector<double> & vcs_VolPhase::moleFractions() const
196{
197 return Xmol_;
198}
199
200double vcs_VolPhase::moleFraction(size_t k) const
201{
202 return Xmol_[k];
203}
204
206 const double* const moleFractions,
207 const int vcsStateStatus)
208{
209 if (totalMoles != 0.0) {
210 // There are other ways to set the mole fractions when VCS_STATECALC
211 // is set to a normal settting.
212 if (vcsStateStatus != VCS_STATECALC_TMP) {
213 throw CanteraError("vcs_VolPhase::setMolesFractionsState",
214 "inappropriate usage");
215 }
216 m_UpToDate = false;
219 throw CanteraError("vcs_VolPhase::setMolesFractionsState",
220 "inappropriate usage");
221 }
223 } else {
224 m_UpToDate = true;
225 m_vcsStateStatus = vcsStateStatus;
227 }
228 double fractotal = 1.0;
230 if (m_totalMolesInert > 0.0) {
232 throw CanteraError("vcs_VolPhase::setMolesFractionsState",
233 "inerts greater than total: {} {}",
235 }
236 fractotal = 1.0 - m_totalMolesInert/v_totalMoles;
237 }
238 double sum = 0.0;
239 for (size_t k = 0; k < m_numSpecies; k++) {
240 Xmol_[k] = moleFractions[k];
241 sum += moleFractions[k];
242 }
243 if (sum == 0.0) {
244 throw CanteraError("vcs_VolPhase::setMolesFractionsState",
245 "inappropriate usage");
246 }
247 if (sum != fractotal) {
248 for (size_t k = 0; k < m_numSpecies; k++) {
249 Xmol_[k] *= (fractotal /sum);
250 }
251 }
253}
254
255void vcs_VolPhase::setMolesFromVCS(const int stateCalc,
256 const double* molesSpeciesVCS)
257{
259
260 if (molesSpeciesVCS == 0) {
261 AssertThrowMsg(m_owningSolverObject, "vcs_VolPhase::setMolesFromVCS",
262 "shouldn't be here");
263 if (stateCalc == VCS_STATECALC_OLD) {
264 molesSpeciesVCS = &m_owningSolverObject->m_molNumSpecies_old[0];
265 } else if (stateCalc == VCS_STATECALC_NEW) {
266 molesSpeciesVCS = &m_owningSolverObject->m_molNumSpecies_new[0];
267 } else {
268 throw CanteraError("vcs_VolPhase::setMolesFromVCS", "shouldn't be here");
269 }
270 } else if (m_owningSolverObject) {
271 // if (stateCalc == VCS_STATECALC_OLD) {
272 // if (molesSpeciesVCS != &m_owningSolverObject->m_molNumSpecies_old[0]) {
273 // throw CanteraError("vcs_VolPhase::setMolesFromVCS", "shouldn't be here");
274 // }
275 // } else if (stateCalc == VCS_STATECALC_NEW) {
276 // if (molesSpeciesVCS != &m_owningSolverObject->m_molNumSpecies_new[0]) {
277 // throw CanteraError("vcs_VolPhase::setMolesFromVCS", "shouldn't be here");
278 // }
279 // }
280 }
281
282 for (size_t k = 0; k < m_numSpecies; k++) {
284 size_t kglob = IndSpecies[k];
285 v_totalMoles += std::max(0.0, molesSpeciesVCS[kglob]);
286 }
287 }
288 if (v_totalMoles > 0.0) {
289 for (size_t k = 0; k < m_numSpecies; k++) {
291 size_t kglob = IndSpecies[k];
292 double tmp = std::max(0.0, molesSpeciesVCS[kglob]);
293 Xmol_[k] = tmp / v_totalMoles;
294 }
295 }
297 } else {
298 // This is where we will start to store a better approximation
299 // for the mole fractions, when the phase doesn't exist.
300 // This is currently unimplemented.
302 }
303
304 // Update the electric potential if it is a solution variable in the
305 // equation system
306 if (m_phiVarIndex != npos) {
307 size_t kglob = IndSpecies[m_phiVarIndex];
308 if (m_numSpecies == 1) {
309 Xmol_[m_phiVarIndex] = 1.0;
310 } else {
311 Xmol_[m_phiVarIndex] = 0.0;
312 }
313 double phi = molesSpeciesVCS[kglob];
315 if (m_numSpecies == 1) {
317 }
318 }
320 if (m_totalMolesInert > 0.0) {
322 }
323
324 // If stateCalc is old and the total moles is positive, then we have a valid
325 // state. If the phase went away, it would be a valid starting point for
326 // F_k's. So, save the state.
327 if (stateCalc == VCS_STATECALC_OLD && v_totalMoles > 0.0) {
329 }
330
331 // Set flags indicating we are up to date with the VCS state vector.
332 m_UpToDate = true;
333 m_vcsStateStatus = stateCalc;
334}
335
336void vcs_VolPhase::setMolesFromVCSCheck(const int vcsStateStatus,
337 const double* molesSpeciesVCS,
338 const double* const TPhMoles)
339{
340 setMolesFromVCS(vcsStateStatus, molesSpeciesVCS);
341
342 // Check for consistency with TPhMoles[]
343 double Tcheck = TPhMoles[VP_ID_];
344 if (Tcheck != v_totalMoles) {
345 if (vcs_doubleEqual(Tcheck, v_totalMoles)) {
346 Tcheck = v_totalMoles;
347 } else {
348 throw CanteraError("vcs_VolPhase::setMolesFromVCSCheck",
349 "We have a consistency problem: {} {}", Tcheck, v_totalMoles);
350 }
351 }
352}
353
354void vcs_VolPhase::updateFromVCS_MoleNumbers(const int vcsStateStatus)
355{
356 if ((!m_UpToDate || vcsStateStatus != m_vcsStateStatus) && m_owningSolverObject &&
357 (vcsStateStatus == VCS_STATECALC_OLD || vcsStateStatus == VCS_STATECALC_NEW)) {
358 setMolesFromVCS(vcsStateStatus);
359 }
360}
361
362void vcs_VolPhase::sendToVCS_ActCoeff(const int vcsStateStatus,
363 double* const AC)
364{
365 updateFromVCS_MoleNumbers(vcsStateStatus);
366 if (!m_UpToDate_AC) {
368 }
369 for (size_t k = 0; k < m_numSpecies; k++) {
370 size_t kglob = IndSpecies[k];
371 AC[kglob] = ActCoeff[k];
372 }
373}
374
375double vcs_VolPhase::sendToVCS_VolPM(double* const VolPM) const
376{
377 if (!m_UpToDate_VolPM) {
378 _updateVolPM();
379 }
380 for (size_t k = 0; k < m_numSpecies; k++) {
381 size_t kglob = IndSpecies[k];
382 VolPM[kglob] = PartialMolarVol[k];
383 }
384 return m_totalVol;
385}
386
387void vcs_VolPhase::sendToVCS_GStar(double* const gstar) const
388{
389 if (!m_UpToDate_GStar) {
390 _updateGStar();
391 }
392 for (size_t k = 0; k < m_numSpecies; k++) {
393 size_t kglob = IndSpecies[k];
394 gstar[kglob] = StarChemicalPotential[k];
395 }
396}
397
399{
400 m_phi = phi;
401 TP_ptr->setElectricPotential(m_phi);
402 // We have changed the state variable. Set uptodate flags to false
403 m_UpToDate_AC = false;
404 m_UpToDate_VolStar = false;
405 m_UpToDate_VolPM = false;
406 m_UpToDate_GStar = false;
407}
408
410{
411 return m_phi;
412}
413
414void vcs_VolPhase::setState_TP(const double temp, const double pres)
415{
416 if (Temp_ == temp && Pres_ == pres) {
417 return;
418 }
419 TP_ptr->setElectricPotential(m_phi);
420 TP_ptr->setState_TP(temp, pres);
421 Temp_ = temp;
422 Pres_ = pres;
423 m_UpToDate_AC = false;
424 m_UpToDate_VolStar = false;
425 m_UpToDate_VolPM = false;
426 m_UpToDate_GStar = false;
427 m_UpToDate_G0 = false;
428}
429
430void vcs_VolPhase::setState_T(const double temp)
431{
432 setState_TP(temp, Pres_);
433}
434
436{
437 TP_ptr->getStandardVolumes(&StarMolarVol[0]);
438 m_UpToDate_VolStar = true;
439}
440
441double vcs_VolPhase::VolStar_calc_one(size_t kspec) const
442{
443 if (!m_UpToDate_VolStar) {
445 }
446 return StarMolarVol[kspec];
447}
448
450{
451 TP_ptr->getPartialMolarVolumes(&PartialMolarVol[0]);
452 m_totalVol = 0.0;
453 for (size_t k = 0; k < m_numSpecies; k++) {
455 }
457
458 if (m_totalMolesInert > 0.0) {
459 if (m_gasPhase) {
460 double volI = m_totalMolesInert * GasConstant * Temp_ / Pres_;
461 m_totalVol += volI;
462 } else {
463 throw CanteraError("vcs_VolPhase::_updateVolPM", "unknown situation");
464 }
465 }
466 m_UpToDate_VolPM = true;
467 return m_totalVol;
468}
469
471{
472 double phaseTotalMoles = v_totalMoles;
473 if (phaseTotalMoles < 1.0E-14) {
474 phaseTotalMoles = 1.0;
475 }
476
477 // Evaluate the current base activity coefficients if necessary
478 if (!m_UpToDate_AC) {
480 }
481 if (!TP_ptr) {
482 return;
483 }
484 TP_ptr->getdlnActCoeffdlnN(m_numSpecies, &np_dLnActCoeffdMolNumber(0,0));
485 for (size_t j = 0; j < m_numSpecies; j++) {
486 double moles_j_base = phaseTotalMoles * Xmol_[j];
487 double* const np_lnActCoeffCol = np_dLnActCoeffdMolNumber.ptrColumn(j);
488 if (moles_j_base < 1.0E-200) {
489 moles_j_base = 1.0E-7 * moles_j_base + 1.0E-13 * phaseTotalMoles + 1.0E-150;
490 }
491 for (size_t k = 0; k < m_numSpecies; k++) {
492 np_lnActCoeffCol[k] = np_lnActCoeffCol[k] * phaseTotalMoles / moles_j_base;
493 }
494 }
495
496 double deltaMoles_j = 0.0;
497 // Make copies of ActCoeff and Xmol_ for use in taking differences
498 vector<double> ActCoeff_Base(ActCoeff);
499 vector<double> Xmol_Base(Xmol_);
500 double TMoles_base = phaseTotalMoles;
501
502 // Loop over the columns species to be deltad
503 for (size_t j = 0; j < m_numSpecies; j++) {
504 // Calculate a value for the delta moles of species j. Note Xmol_[] and
505 // Tmoles are always positive or zero quantities.
506 double moles_j_base = phaseTotalMoles * Xmol_Base[j];
507 deltaMoles_j = 1.0E-7 * moles_j_base + 1.0E-13 * phaseTotalMoles + 1.0E-150;
508
509 // Now, update the total moles in the phase and all of the mole
510 // fractions based on this.
511 phaseTotalMoles = TMoles_base + deltaMoles_j;
512 for (size_t k = 0; k < m_numSpecies; k++) {
513 Xmol_[k] = Xmol_Base[k] * TMoles_base / phaseTotalMoles;
514 }
515 Xmol_[j] = (moles_j_base + deltaMoles_j) / phaseTotalMoles;
516
517 // Go get new values for the activity coefficients. Note this calls
518 // setState_PX();
521
522 // Revert to the base case Xmol_, v_totalMoles
523 v_totalMoles = TMoles_base;
524 Xmol_ = Xmol_Base;
525 }
526
527 // Go get base values for the activity coefficients. Note this calls
528 // setState_TPX() again; Just wanted to make sure that cantera is in sync
529 // with VolPhase after this call.
530 setMoleFractions(&Xmol_Base[0]);
533}
534
536{
537 // update the Ln Act Coeff Jacobian entries with respect to the mole number
538 // of species in the phase -> we always assume that they are out of date.
540
541 // Now copy over the values
542 for (size_t j = 0; j < m_numSpecies; j++) {
543 size_t jglob = IndSpecies[j];
544 for (size_t k = 0; k < m_numSpecies; k++) {
545 size_t kglob = IndSpecies[k];
546 np_LnACJac_VCS(kglob,jglob) = np_dLnActCoeffdMolNumber(k,j);
547 }
548 }
549}
550
552{
553 TP_ptr = tp_ptr;
554 Temp_ = TP_ptr->temperature();
555 Pres_ = TP_ptr->pressure();
557 m_phi = TP_ptr->electricPotential();
558 size_t nsp = TP_ptr->nSpecies();
559 size_t nelem = TP_ptr->nElements();
560 if (nsp != m_numSpecies) {
561 if (m_numSpecies != 0) {
562 warn_user("vcs_VolPhase::setPtrThermoPhase",
563 "Nsp != NVolSpeces: {} {}", nsp, m_numSpecies);
564 }
565 resize(VP_ID_, nsp, nelem, PhaseName.c_str());
566 }
567 TP_ptr->getMoleFractions(&Xmol_[0]);
570
571 // figure out ideal solution tag
572 if (nsp == 1) {
573 m_isIdealSoln = true;
574 } else {
575 m_isIdealSoln = TP_ptr->isIdeal();
576 }
577}
578
580{
581 return TP_ptr;
582}
583
585{
586 return v_totalMoles;
587}
588
589double vcs_VolPhase::molefraction(size_t k) const
590{
591 return Xmol_[k];
592}
593
594void vcs_VolPhase::setCreationMoleNumbers(const double* const n_k,
595 const vector<size_t> &creationGlobalRxnNumbers)
596{
597 creationMoleNumbers_.assign(n_k, n_k+m_numSpecies);
598 for (size_t k = 0; k < m_numSpecies; k++) {
599 creationGlobalRxnNumbers_[k] = creationGlobalRxnNumbers[k];
600 }
601}
602
604 vector<size_t> &creationGlobalRxnNumbers) const
605{
606 creationGlobalRxnNumbers = creationGlobalRxnNumbers_;
608}
609
610void vcs_VolPhase::setTotalMoles(const double totalMols)
611{
612 v_totalMoles = totalMols;
613 if (m_totalMolesInert > 0.0) {
616 "vcs_VolPhase::setTotalMoles",
617 "totalMoles less than inert moles: {} {}",
618 totalMols, m_totalMolesInert);
619 } else {
620 if (m_singleSpecies && (m_phiVarIndex == 0)) {
622 } else {
623 if (totalMols > 0.0) {
625 } else {
627 }
628 }
629 }
630}
631
633{
634 m_UpToDate = false;
635 if (stateCalc != -1) {
636 m_vcsStateStatus = stateCalc;
637 }
638}
639
641{
642 m_UpToDate = true;
643 m_vcsStateStatus = stateCalc;
644}
645
647{
648 return m_isIdealSoln;
649}
650
652{
653 return m_phiVarIndex;
654}
655
656void vcs_VolPhase::setPhiVarIndex(size_t phiVarIndex)
657{
658 m_phiVarIndex = phiVarIndex;
659 m_speciesUnknownType[m_phiVarIndex] = VCS_SPECIES_TYPE_INTERFACIALVOLTAGE;
660 if (m_singleSpecies && m_phiVarIndex == 0) {
661 m_existence = VCS_PHASE_EXIST_ALWAYS;
662 }
663}
664
666{
667 return ListSpeciesPtr[kindex];
668}
669
671{
672 return m_existence;
673}
674
675void vcs_VolPhase::setExistence(const int existence)
676{
677 if (existence == VCS_PHASE_EXIST_NO || existence == VCS_PHASE_EXIST_ZEROEDPHASE) {
678 if (v_totalMoles != 0.0) {
679 throw CanteraError("vcs_VolPhase::setExistence",
680 "setting false existence for phase with moles");
681 }
682 } else if (m_totalMolesInert == 0.0) {
683 if (v_totalMoles == 0.0 && (!m_singleSpecies || m_phiVarIndex != 0)) {
684 throw CanteraError("vcs_VolPhase::setExistence",
685 "setting true existence for phase with no moles");
686 }
687 }
688 if (m_singleSpecies && m_phiVarIndex == 0 && (existence == VCS_PHASE_EXIST_NO || existence == VCS_PHASE_EXIST_ZEROEDPHASE)) {
689 throw CanteraError("vcs_VolPhase::setExistence",
690 "Trying to set existence of an electron phase to false");
691 }
692 m_existence = existence;
693}
694
695size_t vcs_VolPhase::spGlobalIndexVCS(const size_t spIndex) const
696{
697 return IndSpecies[spIndex];
698}
699
700void vcs_VolPhase::setSpGlobalIndexVCS(const size_t spIndex,
701 const size_t spGlobalIndex)
702{
703 IndSpecies[spIndex] = spGlobalIndex;
704 if (spGlobalIndex >= m_numElemConstraints) {
705 creationGlobalRxnNumbers_[spIndex] = spGlobalIndex - m_numElemConstraints;
706 }
707}
708
709void vcs_VolPhase::setTotalMolesInert(const double tMolesInert)
710{
711 if (m_totalMolesInert != tMolesInert) {
712 m_UpToDate = false;
713 m_UpToDate_AC = false;
714 m_UpToDate_VolStar = false;
715 m_UpToDate_VolPM = false;
716 m_UpToDate_GStar = false;
717 m_UpToDate_G0 = false;
718 v_totalMoles += (tMolesInert - m_totalMolesInert);
719 m_totalMolesInert = tMolesInert;
720 }
721 if (m_totalMolesInert > 0.0) {
723 } else if (m_singleSpecies && (m_phiVarIndex == 0)) {
725 } else {
726 if (v_totalMoles > 0.0) {
728 } else {
730 }
731 }
732}
733
735{
736 return m_totalMolesInert;
737}
738
739size_t vcs_VolPhase::elemGlobalIndex(const size_t e) const
740{
741 AssertThrow(e < m_numElemConstraints, " vcs_VolPhase::elemGlobalIndex");
742 return m_elemGlobalIndex[e];
743}
744
745void vcs_VolPhase::setElemGlobalIndex(const size_t eLocal, const size_t eGlobal)
746{
748 "vcs_VolPhase::setElemGlobalIndex");
749 m_elemGlobalIndex[eLocal] = eGlobal;
750}
751
753{
755}
756
757string vcs_VolPhase::elementName(size_t e) const
758{
759
760 if (e < m_elementNames.size()) {
761 return m_elementNames[e];
762 }
763 throw IndexError("vcs_VolPhase::elementName", "element", e, m_elementNames.size());
764}
765
766//! This function decides whether a phase has charged species or not.
767static bool hasChargedSpecies(const ThermoPhase* const tPhase)
768{
769 for (size_t k = 0; k < tPhase->nSpecies(); k++) {
770 if (tPhase->charge(k) != 0.0) {
771 return true;
772 }
773 }
774 return false;
775}
776
777//! This utility routine decides whether a Cantera ThermoPhase needs
778//! a constraint equation representing the charge neutrality of the
779//! phase. It does this by searching for charged species. If it
780//! finds one, and if the phase needs one, then it returns true.
781static bool chargeNeutralityElement(const ThermoPhase* const tPhase)
782{
783 int hasCharge = hasChargedSpecies(tPhase);
784 if (tPhase->chargeNeutralityNecessary() && hasCharge) {
785 return true;
786 }
787 return false;
788}
789
791{
792 size_t nebase = tPhase->nElements();
793 size_t ne = nebase;
794 size_t ns = tPhase->nSpecies();
795
796 // Decide whether we need an extra element constraint for charge
797 // neutrality of the phase
798 bool cne = chargeNeutralityElement(tPhase);
799 if (cne) {
801 ne++;
802 }
803
804 // Assign and resize structures
805 elemResize(ne);
806
809 }
810
811 size_t eFound = npos;
812 if (hasChargedSpecies(tPhase)) {
813 if (cne) {
814 // We need a charge neutrality constraint. We also have an Electron
815 // Element. These are duplicates of each other. To avoid trouble
816 // with possible range error conflicts, sometimes we eliminate the
817 // Electron condition. Flag that condition for elimination by
818 // toggling the ElActive variable. If we find we need it later, we
819 // will retoggle ElActive to true.
820 for (size_t eT = 0; eT < nebase; eT++) {
821 if (tPhase->elementName(eT) == "E") {
822 eFound = eT;
823 m_elementActive[eT] = 0;
825 }
826 }
827 } else {
828 for (size_t eT = 0; eT < nebase; eT++) {
829 if (tPhase->elementName(eT) == "E") {
830 eFound = eT;
832 }
833 }
834 }
835 if (eFound == npos) {
836 eFound = ne;
838 m_elementActive[ne] = 0;
839 string ename = "E";
840 m_elementNames[ne] = ename;
841 ne++;
842 elemResize(ne);
843 }
844 }
845
846 m_formulaMatrix.resize(ns, ne, 0.0);
848 elemResize(ne);
849
850 size_t e = 0;
851 for (size_t eT = 0; eT < nebase; eT++) {
852 m_elementNames[e] = tPhase->elementName(eT);
853 m_elementType[e] = tPhase->elementType(eT);
854 e++;
855 }
856
857 if (cne) {
858 string pname = tPhase->name();
859 if (pname == "") {
860 pname = fmt::format("phase{}", VP_ID_);
861 }
863 m_elementNames[e] = "cn_" + pname;
864 }
865
866 for (size_t k = 0; k < ns; k++) {
867 e = 0;
868 for (size_t eT = 0; eT < nebase; eT++) {
869 m_formulaMatrix(k,e) = tPhase->nAtoms(k, eT);
870 e++;
871 }
872 if (eFound != npos) {
873 m_formulaMatrix(k,eFound) = - tPhase->charge(k);
874 }
875 }
876
877 if (cne) {
878 for (size_t k = 0; k < ns; k++) {
880 }
881 }
882
883 // Here, we figure out what is the species types are The logic isn't set in
884 // stone, and is just for a particular type of problem that I'm solving
885 // first.
886 if (ns == 1 && tPhase->charge(0) != 0.0) {
888 setPhiVarIndex(0);
889 }
890
891 return ne;
892}
893
894int vcs_VolPhase::elementType(const size_t e) const
895{
896 return m_elementType[e];
897}
898
900{
901 return m_formulaMatrix;
902}
903
904int vcs_VolPhase::speciesUnknownType(const size_t k) const
905{
906 return m_speciesUnknownType[k];
907}
908
909int vcs_VolPhase::elementActive(const size_t e) const
910{
911 return m_elementActive[e];
912}
913
915{
916 return m_numSpecies;
917}
918
920{
921 switch (m_eqnState) {
922 case VCS_EOS_CONSTANT:
923 return "Constant";
924 case VCS_EOS_IDEAL_GAS:
925 return "Ideal Gas";
926 case VCS_EOS_STOICH_SUB:
927 return "Stoich Sub";
928 case VCS_EOS_IDEAL_SOLN:
929 return "Ideal Soln";
930 case VCS_EOS_DEBEYE_HUCKEL:
931 return "Debeye Huckel";
932 case VCS_EOS_REDLICH_KWONG:
933 return "Redlich_Kwong";
934 case VCS_EOS_REGULAR_SOLN:
935 return "Regular Soln";
936 default:
937 return fmt::format("UnkType: {:7d}", m_eqnState);
938 break;
939 }
940}
941
942}
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:270
int elementType(size_t m) const
Return the element constraint type Possible types include:
Definition Phase.cpp:109
double nAtoms(size_t k, size_t m) const
Number of atoms of element m in species k.
Definition Phase.cpp:121
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:52
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:605
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.
This is the main structure used to hold the internal data used in vcs_solve_TP(), and to solve TP sys...
Definition vcs_solve.h:45
Properties of a single species.
vector< double > StarChemicalPotential
Vector of calculated Star chemical potentials for the current Temperature and pressure.
void setCreationMoleNumbers(const double *const n_k, const vector< size_t > &creationGlobalRxnNumbers)
Sets the creationMoleNum's within the phase object.
void setElectricPotential(const double phi)
set the electric potential of the phase
bool m_UpToDate_GStar
Boolean indicating whether GStar is up to date.
double electricPotential() const
Returns the electric field of the phase.
size_t m_phiVarIndex
If the potential is a solution variable in VCS, it acts as a species.
vector< double > Xmol_
Vector of the current mole fractions for species in the phase.
double Temp_
Current value of the temperature for this object, and underlying objects.
int speciesUnknownType(const size_t k) const
Returns the type of the species unknown.
size_t elemGlobalIndex(const size_t e) const
Returns the global index of the local element index for the phase.
void _updateGStar() const
Gibbs free energy calculation for standard states.
vector< double > SS0ChemicalPotential
Vector of calculated SS0 chemical potentials for the current Temperature.
size_t nSpecies() const
Return the number of species in the phase.
string PhaseName
String name for the phase.
void setMolesFromVCSCheck(const int vcsStateStatus, const double *molesSpeciesVCS, const double *const TPhMoles)
Set the moles within the phase.
void setMoleFractionsState(const double molNum, const double *const moleFracVec, const int vcsStateStatus)
Set the moles and/or mole fractions within the phase.
void setMolesOutOfDate(int stateCalc=-1)
Sets the mole flag within the object to out of date.
vector< string > m_elementNames
vector of strings containing the element constraint names
vector< int > m_elementActive
boolean indicating whether an element constraint is active for the current problem
vector< double > creationMoleNumbers_
Vector of current creationMoleNumbers_.
double _updateVolPM() const
Calculate the partial molar volumes of all species and return the total volume.
double m_totalVol
Total Volume of the phase. Units are m**3.
size_t ChargeNeutralityElement
This is the element number for the charge neutrality condition of the phase.
vector< size_t > m_elemGlobalIndex
Index of the element number in the global list of elements stored in VCS_SOLVE.
size_t VP_ID_
Original ID of the phase in the problem.
Array2D np_dLnActCoeffdMolNumber
Vector of the derivatives of the ln activity coefficient wrt to the current mole number multiplied by...
const Array2D & getFormulaMatrix() const
Get a constant form of the Species Formula Matrix.
vector< int > m_speciesUnknownType
Type of the species unknown.
bool m_UpToDate_VolStar
Boolean indicating whether Star volumes are up to date.
void updateFromVCS_MoleNumbers(const int stateCalc)
Update the moles within the phase, if necessary.
vcs_SpeciesProperties * speciesProperty(const size_t kindex)
Retrieve the kth Species structure for the species belonging to this phase.
size_t phiVarIndex() const
Return the index of the species that represents the the voltage of the phase.
double totalMolesInert() const
Returns the value of the total kmol of inert in the phase.
double v_totalMoles
Total mols in the phase. units are kmol.
bool m_singleSpecies
If true, this phase consists of a single species.
void setMolesFromVCS(const int stateCalc, const double *molesSpeciesVCS=0)
Set the moles within the phase.
void setElemGlobalIndex(const size_t eLocal, const size_t eGlobal)
sets a local phase element to a global index value
int m_MFStartIndex
This is always equal to zero.
void resize(const size_t phaseNum, const size_t numSpecies, const size_t numElem, const char *const phaseName, const double molesInert=0.0)
The resize() function fills in all of the initial information if it is not given in the constructor.
bool m_isIdealSoln
Boolean indicating whether the phase is an ideal solution and therefore its molar-based activity coef...
vector< double > ActCoeff
Vector of calculated activity coefficients for the current state.
int exists() const
int indicating whether the phase exists or not
vector< size_t > IndSpecies
Index into the species vectors.
const vector< double > & moleFractions() const
Return a const reference to the mole fractions stored in the object.
size_t nElemConstraints() const
Returns the number of element constraints.
double m_phi
Value of the potential for the phase (Volts)
size_t spGlobalIndexVCS(const size_t spIndex) const
Return the Global VCS index of the kth species in the phase.
size_t transferElementsFM(const ThermoPhase *const tPhase)
Transfer all of the element information from the ThermoPhase object to the vcs_VolPhase object.
vector< int > m_elementType
Type of the element constraint.
double VolStar_calc_one(size_t kspec) const
Molar volume calculation for standard state of one species.
double GStar_calc_one(size_t kspec) const
Gibbs free energy calculation for standard state of one species.
bool m_gasPhase
If true, this phase is a gas-phase like phase.
void _updateVolStar() const
Molar volume calculation for standard states.
int m_vcsStateStatus
Status.
bool m_UpToDate_G0
Boolean indicating whether G0 is up to date.
void _updateG0() const
Gibbs free energy calculation at a temperature for the reference state of each species.
double Pres_
Current value of the pressure for this object, and underlying objects.
vector< double > PartialMolarVol
Vector of the Partial molar Volumes of the species. units m3 / kmol.
double molefraction(size_t kspec) const
Returns the mole fraction of the kspec species.
Array2D m_formulaMatrix
Formula Matrix for the phase.
void _updateLnActCoeffJac()
Evaluation of Activity Coefficient Jacobians.
string elementName(size_t e) const
Name of the element constraint with index e.
void setTotalMoles(const double totalMols)
Sets the total moles in the phase.
void setMolesCurrent(int vcsStateStatus)
Sets the mole flag within the object to be current.
void sendToVCS_GStar(double *const gstar) const
Fill in the standard state Gibbs free energy vector for VCS.
void sendToVCS_ActCoeff(const int stateCalc, double *const AC)
Fill in an activity coefficients vector within a VCS_SOLVE object.
size_t m_numElemConstraints
Number of element constraints within the problem.
vector< double > StarMolarVol
Vector of the Star molar Volumes of the species. units m3 / kmol.
double sendToVCS_VolPM(double *const VolPM) const
Fill in the partial molar volume vector for VCS.
void sendToVCS_LnActCoeffJac(Array2D &LnACJac_VCS)
Downloads the ln ActCoeff Jacobian into the VCS version of the ln ActCoeff Jacobian.
int m_eqnState
Type of the equation of state.
void setTotalMolesInert(const double tMolesInert)
Sets the total moles of inert in the phase.
double G0_calc_one(size_t kspec) const
Gibbs free energy calculation at a temperature for the reference state of a species,...
void _updateActCoeff() const
Evaluate the activity coefficients at the current conditions.
void setState_TP(const double temperature_Kelvin, const double pressure_PA)
Sets the temperature and pressure in this object and underlying ThermoPhase objects.
int m_existence
Current state of existence:
string eos_name() const
Return the name corresponding to the equation of state.
bool m_UpToDate_AC
Boolean indicating whether activity coefficients are up to date.
VCS_SOLVE * m_owningSolverObject
Backtrack value of VCS_SOLVE *.
ThermoPhase * TP_ptr
If we are using Cantera, this is the pointer to the ThermoPhase object.
void _updateMoleFractionDependencies()
Updates the mole fraction dependencies.
double m_totalMolesInert
Total moles of inert in the phase.
int elementType(const size_t e) const
Type of the element constraint with index e.
const ThermoPhase * ptrThermoPhase() const
Return a const ThermoPhase pointer corresponding to this phase.
vector< size_t > creationGlobalRxnNumbers_
Vector of creation global reaction numbers for the phase stability problem.
void setSpGlobalIndexVCS(const size_t spIndex, const size_t spGlobalIndex)
set the Global VCS index of the kth species in the phase
size_t m_numSpecies
Number of species in the phase.
void setExistence(const int existence)
Set the existence flag in the object.
void setMoleFractions(const double *const xmol)
Set the mole fractions from a conventional mole fraction vector.
bool m_UpToDate
Boolean indicating whether the object has an up-to-date mole number vector and potential with respect...
vector< vcs_SpeciesProperties * > ListSpeciesPtr
Vector of Species structures for the species belonging to this phase.
void setState_T(const double temperature_Kelvin)
Sets the temperature in this object and underlying ThermoPhase objects.
const vector< double > & creationMoleNumbers(vector< size_t > &creationGlobalRxnNumbers) const
Return a const reference to the creationMoleNumbers stored in the object.
bool m_UpToDate_VolPM
Boolean indicating whether partial molar volumes are up to date.
bool isIdealSoln() const
Returns whether the phase is an ideal solution phase.
double totalMoles() const
Return the total moles in the phase.
void setPtrThermoPhase(ThermoPhase *tp_ptr)
Set the pointer for Cantera's ThermoPhase parameter.
#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.
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition global.h:268
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:286
#define VCS_STATECALC_OLD
State Calculation based on the old or base mole numbers.
Definition vcs_defs.h:297
#define VCS_STATECALC_TMP
State Calculation based on a temporary set of mole numbers.
Definition vcs_defs.h:308
#define VCS_PHASE_EXIST_NO
Phase doesn't currently exist in the mixture.
Definition vcs_defs.h:200
#define VCS_PHASE_EXIST_ALWAYS
These defines are valid values for the phase existence flag.
Definition vcs_defs.h:187
#define VCS_ELEM_TYPE_ABSPOS
Normal element constraint consisting of positive coefficients for the formula matrix.
Definition vcs_defs.h:226
#define VCS_STATECALC_NEW
State Calculation based on the new or tentative mole numbers.
Definition vcs_defs.h:300
#define VCS_SPECIES_TYPE_MOLNUM
Unknown refers to mole number of a single species.
Definition vcs_defs.h:278
#define VCS_ELEM_TYPE_ELECTRONCHARGE
This refers to conservation of electrons.
Definition vcs_defs.h:232
#define VCS_ELEM_TYPE_CHARGENEUTRALITY
This refers to a charge neutrality of a single phase.
Definition vcs_defs.h:238
#define VCS_PHASE_EXIST_YES
Phase is a normal phase that currently exists.
Definition vcs_defs.h:190
#define VCS_PHASE_EXIST_ZEROEDPHASE
Phase currently is zeroed due to a programmatic issue.
Definition vcs_defs.h:206
Header file for the internal object that holds the vcs equilibrium problem (see Class VCS_SOLVE and C...