Cantera  4.0.0a1
Loading...
Searching...
No Matches
PlasmaPhase.cpp
Go to the documentation of this file.
1//! @file PlasmaPhase.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
8#include <boost/math/special_functions/gamma.hpp>
10#include "cantera/base/global.h"
11#include "cantera/numerics/eigen_dense.h"
15#include <boost/polymorphic_pointer_cast.hpp>
17
18namespace Cantera {
19
20namespace {
21 const double gamma = sqrt(2 * ElectronCharge / ElectronMass);
22}
23
24PlasmaPhase::PlasmaPhase(const string& inputFile, const string& id_)
25{
26 initThermoFile(inputFile, id_);
27
28 // initial electron temperature
30
31 // Initialize the Boltzmann Solver
32 m_eedfSolver = make_unique<EEDFTwoTermApproximation>(this);
33
34 // Set Energy Grid (Hardcoded Defaults for Now)
35 double kTe_max = 60;
36 size_t nGridCells = 301;
37 m_nPoints = nGridCells + 1;
38 m_eedfSolver->setLinearGrid(kTe_max, nGridCells);
41}
42
43PlasmaPhase::~PlasmaPhase()
44{
45 if (shared_ptr<Solution> soln = m_soln.lock()) {
46 soln->removeChangedCallback(this);
47 soln->kinetics()->removeReactionAddedCallback(this);
48 }
49 for (size_t k = 0; k < nCollisions(); k++) {
50 // remove callback
51 m_collisions[k]->removeSetRateCallback(this);
52 }
53}
54
56{
57 if (m_distributionType == "discretized") {
58 throw CanteraError("PlasmaPhase::updateElectronEnergyDistribution",
59 "Invalid for discretized electron energy distribution.");
60 } else if (m_distributionType == "isotropic") {
62 } else if (m_distributionType == "Boltzmann-two-term") {
63 auto ierr = m_eedfSolver->calculateDistributionFunction();
64 if (ierr == 0) {
65 auto y = m_eedfSolver->getEEDFEdge();
66 m_electronEnergyDist = Eigen::Map<const Eigen::ArrayXd>(y.data(), m_nPoints);
67 } else {
68 throw CanteraError("PlasmaPhase::updateElectronEnergyDistribution",
69 "Call to calculateDistributionFunction failed.");
70 }
71 bool validEEDF = (
72 static_cast<size_t>(m_electronEnergyDist.size()) == m_nPoints &&
73 m_electronEnergyDist.allFinite() &&
74 m_electronEnergyDist.maxCoeff() > 0.0 &&
75 m_electronEnergyDist.sum() > 0.0
76 );
77
78 if (validEEDF) {
80 } else {
81 writelog("Skipping Te update: EEDF is empty, non-finite, or unnormalized.\n");
82 }
83 } else {
84 throw CanteraError("PlasmaPhase::updateElectronEnergyDistribution",
85 "Unknown method '{}' for determining EEDF", m_distributionType);
86 }
89}
90
92 Eigen::ArrayXd eps32 = m_electronEnergyLevels.pow(3./2.);
93 double norm = 2./3. * numericalQuadrature(m_quadratureMethod,
95 if (norm < 0.0) {
96 throw CanteraError("PlasmaPhase::normalizeElectronEnergyDistribution",
97 "The norm is negative. This might be caused by bad "
98 "electron energy distribution");
99 }
100 m_electronEnergyDist /= norm;
101}
102
104{
105 if (type == "discretized" ||
106 type == "isotropic" ||
107 type == "Boltzmann-two-term") {
109 } else {
110 throw CanteraError("PlasmaPhase::setElectronEnergyDistributionType",
111 "Unknown type for electron energy distribution.");
112 }
113}
114
116{
118 double x = m_isotropicShapeFactor;
119 double gamma1 = boost::math::tgamma(3.0 / 2.0 / x);
120 double gamma2 = boost::math::tgamma(5.0 / 2.0 / x);
121 double c1 = x * std::pow(gamma2, 1.5) / std::pow(gamma1, 2.5);
122 double c2 = std::pow(gamma2 / gamma1, x);
124 c1 / std::pow(meanElectronEnergy(), 1.5) *
125 (-c2 * (m_electronEnergyLevels /
126 meanElectronEnergy()).pow(x)).exp();
128}
129
131 m_electronTemp = Te;
133}
134
136{
138
139 if (!m_inEquilibrate) {
140 m_inEquilibrate = true;
141 // Remember current Te and lock Te -> T for the duration
144 }
145}
146
148{
149 if (m_inEquilibrate) {
150 // Restore Te to the pre-equilibrate value
152 m_inEquilibrate = false;
153 }
154
156}
157
159 m_electronTemp = 2.0 / 3.0 * energy * ElectronCharge / Boltzmann;
161}
162
163void PlasmaPhase::setElectronEnergyLevels(span<const double> levels)
164{
165 m_nPoints = levels.size();
166 m_electronEnergyLevels = Eigen::Map<const Eigen::ArrayXd>(levels.data(), m_nPoints);
170}
171
173{
174 m_distNum++;
175}
176
178{
179 m_levelNum++;
180 // Cross sections are interpolated on the energy levels
181 if (m_collisions.size() > 0) {
182 for (shared_ptr<Reaction> collision : m_collisions) {
183 const auto& rate = boost::polymorphic_pointer_downcast
186 }
187 }
188}
189
191{
192 Eigen::ArrayXd h = m_electronEnergyLevels.tail(m_nPoints - 1) -
194 if (m_electronEnergyLevels[0] < 0.0 || (h <= 0.0).any()) {
195 throw CanteraError("PlasmaPhase::checkElectronEnergyLevels",
196 "Values of electron energy levels need to be positive and "
197 "monotonically increasing.");
198 }
199}
200
202{
203 Eigen::ArrayXd h = m_electronEnergyLevels.tail(m_nPoints - 1) -
205 if ((m_electronEnergyDist < 0.0).any()) {
206 throw CanteraError("PlasmaPhase::checkElectronEnergyDistribution",
207 "Values of electron energy distribution cannot be negative.");
208 }
209 if (m_electronEnergyDist[m_nPoints - 1] > 0.01) {
210 warn_user("PlasmaPhase::checkElectronEnergyDistribution",
211 "The value of the last element of electron energy distribution exceed 0.01. "
212 "This indicates that the value of electron energy level is not high enough "
213 "to contain the isotropic distribution at mean electron energy of "
214 "{} eV", meanElectronEnergy());
215 }
216}
217
219 span<const double> dist)
220{
221 m_distributionType = "discretized";
222 m_nPoints = levels.size();
228 }
234}
235
237{
238 // calculate mean electron energy and electron temperature
239 Eigen::ArrayXd eps52 = m_electronEnergyLevels.pow(5./2.);
240 double epsilon_m = 2.0 / 5.0 * numericalQuadrature(m_quadratureMethod,
241 m_electronEnergyDist, eps52);
242 if (epsilon_m < 0.0 && m_quadratureMethod == "simpson") {
243 // try trapezoidal method
244 epsilon_m = 2.0 / 5.0 * numericalQuadrature(
245 "trapezoidal", m_electronEnergyDist, eps52);
246 }
247
248 if (epsilon_m < 0.0) {
249 throw CanteraError("PlasmaPhase::updateElectronTemperatureFromEnergyDist",
250 "The electron energy distribution produces negative electron temperature.");
251 }
252
253 m_electronTemp = 2.0 / 3.0 * epsilon_m * ElectronCharge / Boltzmann;
254}
255
257 m_isotropicShapeFactor = x;
259}
260
261void PlasmaPhase::getParameters(AnyMap& phaseNode) const
262{
264 AnyMap eedf;
265 eedf["type"] = m_distributionType;
266 vector<double> levels(m_nPoints);
267 Eigen::Map<Eigen::ArrayXd>(levels.data(), m_nPoints) = m_electronEnergyLevels;
268 eedf["energy-levels"] = levels;
269 if (m_distributionType == "isotropic") {
270 eedf["shape-factor"] = m_isotropicShapeFactor;
271 eedf["mean-electron-energy"].setQuantity(meanElectronEnergy(), "eV");
272 } else if (m_distributionType == "discretized") {
273 vector<double> dist(m_nPoints);
274 Eigen::Map<Eigen::ArrayXd>(dist.data(), m_nPoints) = m_electronEnergyDist;
275 eedf["distribution"] = dist;
276 eedf["normalize"] = m_do_normalizeElectronEnergyDist;
277 }
278 phaseNode["electron-energy-distribution"] = std::move(eedf);
279}
280
281void PlasmaPhase::setParameters(const AnyMap& phaseNode, const AnyMap& rootNode)
282{
283 IdealGasPhase::setParameters(phaseNode, rootNode);
284 if (phaseNode.hasKey("electron-energy-distribution")) {
285 const AnyMap eedf = phaseNode["electron-energy-distribution"].as<AnyMap>();
286 m_distributionType = eedf["type"].asString();
287 if (m_distributionType == "isotropic") {
288 if (eedf.hasKey("shape-factor")) {
289 setIsotropicShapeFactor(eedf["shape-factor"].asDouble());
290 } else {
291 throw CanteraError("PlasmaPhase::setParameters",
292 "isotropic type requires shape-factor key.");
293 }
294 if (eedf.hasKey("mean-electron-energy")) {
295 double energy = eedf.convert("mean-electron-energy", "eV");
296 setMeanElectronEnergy(energy);
297 } else {
298 throw CanteraError("PlasmaPhase::setParameters",
299 "isotropic type requires electron-temperature key.");
300 }
301 if (eedf.hasKey("energy-levels")) {
302 auto levels = eedf["energy-levels"].asVector<double>();
304 }
306 } else if (m_distributionType == "discretized") {
307 if (!eedf.hasKey("energy-levels")) {
308 throw CanteraError("PlasmaPhase::setParameters",
309 "Cannot find key energy-levels.");
310 }
311 if (!eedf.hasKey("distribution")) {
312 throw CanteraError("PlasmaPhase::setParameters",
313 "Cannot find key distribution.");
314 }
315 if (eedf.hasKey("normalize")) {
316 enableNormalizeElectronEnergyDist(eedf["normalize"].asBool());
317 }
318 auto levels = eedf["energy-levels"].asVector<double>();
319 auto distribution = eedf["distribution"].asVector<double>(levels.size());
320 setDiscretizedElectronEnergyDist(levels, distribution);
321 }
322 }
323
324 if (rootNode.hasKey("electron-collisions")) {
325 for (const auto& item : rootNode["electron-collisions"].asVector<AnyMap>()) {
326 auto rate = make_shared<ElectronCollisionPlasmaRate>(item);
327 Composition reactants, products;
328 reactants[item["target"].asString()] = 1;
329 reactants[electronSpeciesName()] = 1;
330 if (item.hasKey("product")) {
331 products[item["product"].asString()] = 1;
332 } else {
333 products[item["target"].asString()] = 1;
334 }
335 products[electronSpeciesName()] = 1;
336 if (rate->kind() == "ionization") {
337 products[electronSpeciesName()] += 1;
338 } else if (rate->kind() == "attachment") {
339 products[electronSpeciesName()] -= 1;
340 }
341 auto R = make_shared<Reaction>(reactants, products, rate);
342 addCollision(R);
343 }
344 }
345}
346
347bool PlasmaPhase::addSpecies(shared_ptr<Species> spec)
348{
349 bool added = IdealGasPhase::addSpecies(spec);
350 size_t k = m_kk - 1;
351
352 if ((spec->name == "e" || spec->name == "Electron") ||
353 (spec->composition.find("E") != spec->composition.end() &&
354 spec->composition.size() == 1 &&
355 spec->composition["E"] == 1)) {
358 } else {
359 throw CanteraError("PlasmaPhase::addSpecies",
360 "Cannot add species, {}. "
361 "Only one electron species is allowed.", spec->name);
362 }
363 }
364 return added;
365}
366
368{
370
371 // Check electron species
373 throw CanteraError("PlasmaPhase::initThermo",
374 "No electron species found.");
375 }
376}
377
378void PlasmaPhase::setSolution(std::weak_ptr<Solution> soln) {
380 // register callback function to be executed
381 // when the thermo or kinetics object changed
382 if (shared_ptr<Solution> soln = m_soln.lock()) {
383 soln->registerChangedCallback(this, [&]() {
385 });
386 }
387}
388
390{
391 if (shared_ptr<Solution> soln = m_soln.lock()) {
392 shared_ptr<Kinetics> kin = soln->kinetics();
393 if (!kin) {
394 return;
395 }
396
397 // add collision from the initial list of reactions. Only add reactions we
398 // haven't seen before
399 set<Reaction*> existing;
400 for (auto& R : m_collisions) {
401 existing.insert(R.get());
402 }
403 for (size_t i = 0; i < kin->nReactions(); i++) {
404 shared_ptr<Reaction> R = kin->reaction(i);
405 if (R->rate()->type() != "electron-collision-plasma"
406 || existing.count(R.get())) {
407 continue;
408 }
409 addCollision(R);
410 }
411
412 // register callback when reaction is added later
413 // Modifying collision reactions is not supported
414 kin->registerReactionAddedCallback(this, [this, kin]() {
415 size_t i = kin->nReactions() - 1;
416 if (kin->reaction(i)->type() == "electron-collision-plasma") {
417 addCollision(kin->reaction(i));
418 }
419 });
420 }
421}
422
423void PlasmaPhase::addCollision(shared_ptr<Reaction> collision)
424{
425 size_t i = nCollisions();
426
427 // setup callback to signal updating the cross-section-related
428 // parameters
429 collision->registerSetRateCallback(this, [this, i, collision]() {
430 m_interp_cs_ready[i] = false;
432 std::dynamic_pointer_cast<ElectronCollisionPlasmaRate>(collision->rate());
433 });
434
435 // Identify target species for electron-collision reactions
436 string target;
437 for (const auto& [name, _] : collision->reactants) {
438 // Reactants are expected to be electrons and the target species
439 if (name != electronSpeciesName()) {
440 m_targetSpeciesIndices.emplace_back(speciesIndex(name, true));
441 target = name;
442 break;
443 }
444 }
445 if (target.empty()) {
446 throw CanteraError("PlasmaPhase::addCollision", "Error identifying target for"
447 " collision with equation '{}'", collision->equation());
448 }
449
450 m_collisions.emplace_back(collision);
451 m_collisionRates.emplace_back(
452 std::dynamic_pointer_cast<ElectronCollisionPlasmaRate>(collision->rate()));
453 m_interp_cs_ready.emplace_back(false);
454
455 // resize parameters
458
459 // Set up data used by Boltzmann solver
460 auto& rate = *m_collisionRates.back();
461 string kind = m_collisionRates.back()->kind();
462
463 if ((kind == "effective" || kind == "elastic")) {
464 for (size_t k = 0; k < m_collisions.size() - 1; k++) {
465 if (m_collisions[k]->reactants == collision->reactants &&
466 (m_collisionRates[k]->kind() == "elastic" ||
467 m_collisionRates[k]->kind() == "effective") && !collision->duplicate)
468 {
469 throw CanteraError("PlasmaPhase::addCollision", "Phase already contains"
470 " an effective/elastic cross section for '{}'.", target);
471 }
472 }
473 m_kElastic.push_back(i);
474 } else {
475 m_kInelastic.push_back(i);
476 }
477
478 auto levels = rate.energyLevels();
479 m_energyLevels.emplace_back(levels.begin(), levels.end());
480 auto sections = rate.crossSections();
481 m_crossSections.emplace_back(sections.begin(), sections.end());
482 m_eedfSolver->setGridCache();
483}
484
486{
487 if (m_interp_cs_ready[i]) {
488 return false;
489 }
490 vector<double> levels(m_nPoints);
491 Eigen::Map<Eigen::ArrayXd>(levels.data(), m_nPoints) = m_electronEnergyLevels;
492 m_collisionRates[i]->updateInterpolatedCrossSection(levels);
493 m_interp_cs_ready[i] = true;
494 return true;
495}
496
498{
500 // Forward difference for the first point
504
505 // Central difference for the middle points
506 for (size_t i = 1; i < m_nPoints - 1; i++) {
510 (h1 * h1 - h0 * h0) * m_electronEnergyDist[i] -
511 h1 * h1 * m_electronEnergyDist[i-1]) /
512 (h1 * h0) / (h1 + h0);
513 }
514
515 // Backward difference for the last point
521}
522
524{
525 // cache of cross section plus distribution plus energy-level number
526 static const int cacheId = m_cache.getId();
527 CachedScalar last_stateNum = m_cache.getScalar(cacheId);
528
529 // combine the distribution and energy level number
530 int stateNum = m_distNum + m_levelNum;
531
532 vector<bool> interpChanged(m_collisions.size());
533 for (size_t i = 0; i < m_collisions.size(); i++) {
534 interpChanged[i] = updateInterpolatedCrossSection(i);
535 }
536
537 if (last_stateNum.validate(temperature(), stateNum)) {
538 // check each cross section, and only update coefficients that
539 // the interpolated cross sections change
540 for (size_t i = 0; i < m_collisions.size(); i++) {
541 if (interpChanged[i]) {
543 }
544 }
545 } else {
546 // update every coefficient if distribution, temperature,
547 // or energy levels change.
548 for (size_t i = 0; i < m_collisions.size(); i++) {
550 }
551 }
552}
553
555{
556 // @todo exclude attachment collisions
557 size_t k = m_targetSpeciesIndices[i];
558
559 // Map cross sections to Eigen::ArrayXd
560 auto cs_array = Eigen::Map<const Eigen::ArrayXd>(
561 m_collisionRates[i]->crossSectionInterpolated().data(),
562 m_collisionRates[i]->crossSectionInterpolated().size()
563 );
564
565 // Mass ratio calculation
566 double mass_ratio = ElectronMass / molecularWeight(k) * Avogadro;
567
568 // Calculate the rate using Simpson's rule or trapezoidal rule
569 Eigen::ArrayXd f0_plus = m_electronEnergyDist + Boltzmann * temperature() /
571 m_elasticElectronEnergyLossCoefficients[i] = 2.0 * mass_ratio * gamma *
573 m_quadratureMethod, 1.0 / 3.0 * f0_plus.cwiseProduct(cs_array),
574 m_electronEnergyLevels.pow(3.0));
575}
576
578{
579 if (m_electronEnergyDist.size() != m_nPoints
580 || m_electronEnergyDistDiff.size() != m_nPoints) {
581 throw CanteraError("PlasmaPhase::elasticPowerLoss:",
582 "EEDF not initialized");
583 }
584
586 // The elastic power loss includes the contributions from inelastic
587 // collisions (inelastic recoil effects).
588 double rate = 0.0;
589 for (size_t i = 0; i < nCollisions(); i++) {
592 }
593 const double q_elastic = Avogadro * Avogadro * ElectronCharge *
595
596 if (!std::isfinite(q_elastic)) {
597 throw CanteraError("PlasmaPhase::elasticPowerLoss:",
598 "Non-finite elastic power loss");
599 }
600
601 return q_elastic;
602}
603
605{
607 static const int cacheId = m_cache.getId();
608 CachedScalar cached = m_cache.getScalar(cacheId);
609 double tempNow = temperature();
610 double electronTempNow = electronTemperature();
611 size_t k = m_electronSpeciesIndex;
612 // If the electron temperature has changed since the last time these
613 // properties were computed, recompute them.
614 if (cached.state1 != tempNow || cached.state2 != electronTempNow) {
616 m_cp0_R[k], m_h0_RT[k], m_s0_R[k]);
617 cached.state1 = tempNow;
618 cached.state2 = electronTempNow;
619 }
620 // update the species Gibbs functions
621 m_g0_RT[k] = m_h0_RT[k] - m_s0_R[k];
622}
623
625 double value = IdealGasPhase::enthalpy_mole();
626 value += GasConstant * (electronTemperature() - temperature()) *
629 return value;
630}
631
633{
634 m_work.resize(m_kk);
636 double u = 0.0;
637 for (size_t k = 0; k < m_kk; ++k) {
638 u += moleFraction(k) * m_work[k];
639 }
640 return u;
641}
642
644{
645 m_work.resize(m_kk);
647 double s = 0.0;
648 for (size_t k = 0; k < m_kk; ++k) {
649 s += moleFraction(k) * m_work[k];
650 }
651 return s;
652}
653
655{
656 m_work.resize(m_kk);
658 double g = 0.0;
659 for (size_t k = 0; k < m_kk; ++k) {
660 g += moleFraction(k) * m_work[k];
661 }
662 return g;
663}
664
665void PlasmaPhase::getGibbs_ref(span<double> g) const
666{
669}
670
671void PlasmaPhase::getStandardVolumes_ref(span<double> vol) const
672{
675}
676
677void PlasmaPhase::getPartialMolarEnthalpies(span<double> hbar) const
678{
681}
682
683void PlasmaPhase::getPartialMolarEntropies(span<double> sbar) const
684{
686 double logp = log(pressure());
687 double logpe = log(electronPressure());
688 sbar[m_electronSpeciesIndex] += GasConstant * (logp - logpe);
689}
690
691void PlasmaPhase::getPartialMolarIntEnergies(span<double> ubar) const
692{
693 checkArraySize("PlasmaPhase::getPartialMolarIntEnergies", ubar.size(), m_kk);
694 auto _h = enthalpy_RT_ref();
695 for (size_t k = 0; k < m_kk; k++) {
696 ubar[k] = RT() * (_h[k] - 1.0);
697 }
698 size_t k = m_electronSpeciesIndex;
699 ubar[k] = RTe() * (_h[k] - 1.0);
700}
701
702void PlasmaPhase::getChemPotentials(span<double> mu) const
703{
705 size_t k = m_electronSpeciesIndex;
706 double xx = std::max(SmallNumber, moleFraction(k));
707 mu[k] += (RTe() - RT()) * log(xx);
708}
709
710void PlasmaPhase::getStandardChemPotentials(span<double> muStar) const
711{
713 size_t k = m_electronSpeciesIndex;
714 muStar[k] -= log(pressure() / refPressure()) * RT();
715 muStar[k] += log(electronPressure() / refPressure()) * RTe();
716}
717
718void PlasmaPhase::getEntropy_R(span<double> sr) const
719{
720 checkArraySize("PlasmaPhase::getEntropy_R", sr.size(), m_kk);
721 auto _s = entropy_R_ref();
722 copy(_s.begin(), _s.end(), sr.begin());
723 double tmp = log(pressure() / refPressure());
724 for (size_t k = 0; k < m_kk; k++) {
725 if (k != m_electronSpeciesIndex) {
726 sr[k] -= tmp;
727 } else {
728 sr[k] -= log(electronPressure() / refPressure());
729 }
730 }
731}
732
733void PlasmaPhase::getGibbs_RT(span<double> grt) const
734{
735 checkArraySize("PlasmaPhase::getGibbs_RT", grt.size(), m_kk);
736 auto gibbsrt = gibbs_RT_ref();
737 copy(gibbsrt.begin(), gibbsrt.end(), grt.begin());
738 double tmp = log(pressure() / refPressure());
739 for (size_t k = 0; k < m_kk; k++) {
740 if (k != m_electronSpeciesIndex) {
741 grt[k] += tmp;
742 } else {
743 grt[k] += log(electronPressure() / refPressure());
744 }
745 }
746}
747
749{
750 // Only implemented when using the Boltzmann two-term EEDF
751 if (m_distributionType == "Boltzmann-two-term") {
752 return m_eedfSolver->getElectronMobility();
753 } else {
754 throw NotImplementedError("PlasmaPhase::electronMobility",
755 "Electron mobility is only available for 'Boltzmann-two-term' "
756 "electron energy distributions.");
757 }
758}
759
760
762{
763 // sigma = e * n_e * mu_e [S/m]; q_J = sigma * E^2 [W/m^3]
764 const double mu_e = electronMobility(); // m^2 / (V·s)
765 if (mu_e <= 0.0) {
766 return 0.0;
767 }
768 const double ne = concentration(m_electronSpeciesIndex) * Avogadro; // m^-3
769 if (ne <= 0.0) {
770 return 0.0;
771 }
772 const double E = electricField(); // V/m
773 if (E <= 0.0) {
774 return 0.0;
775 }
776 const double sigma = ElectronCharge * ne * mu_e; // S/m
777 return sigma * E * E; // W/m^3
778}
779
781{
782 // Joule heating: sigma * E^2 [W/m^3]
783 const double qJ = jouleHeatingPower();
784 checkFinite(qJ);
785
786 // Elastic + inelastic recoil power loss [W/m^3]
787 double qElastic = elasticPowerLoss();
788 checkFinite(qElastic);
789
790 return qJ + qElastic;
791}
792
793
794}
EEDF Two-Term approximation solver.
Header for plasma reaction rates parameterized by electron collision cross section and electron energ...
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
Header file for class PlasmaPhase.
Declaration for class Cantera::Species.
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
Definition AnyMap.cpp:1477
double convert(const string &key, const string &units) const
Convert the item stored by the given key to the units specified in units.
Definition AnyMap.cpp:1595
Base class for exceptions thrown by Cantera classes.
Electron collision plasma reaction rate type.
void updateInterpolatedCrossSection(span< const double >)
Update the value of m_crossSectionsInterpolated [m2].
void getGibbs_ref(span< double > g) const override
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
double enthalpy_mole() const override
Return the Molar enthalpy. Units: J/kmol.
void getPartialMolarEnthalpies(span< double > hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
double pressure() const override
Pressure.
vector< double > m_g0_RT
Temporary storage for dimensionless reference state Gibbs energies.
vector< double > m_h0_RT
Temporary storage for dimensionless reference state enthalpies.
span< const double > enthalpy_RT_ref() const
Returns a reference to the dimensionless reference state enthalpy vector.
void getStandardChemPotentials(span< double > mu) const override
Get the array of chemical potentials at unit activity for the species at their standard states at the...
virtual void updateThermo() const
Update the species reference state thermodynamic functions.
vector< double > m_s0_R
Temporary storage for dimensionless reference state entropies.
span< const double > entropy_R_ref() const
Returns a reference to the dimensionless reference state Entropy vector.
void getPartialMolarEntropies(span< double > sbar) const override
Returns an array of partial molar entropies of the species in the solution.
vector< double > m_cp0_R
Temporary storage for dimensionless reference state heat capacities.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
span< const double > gibbs_RT_ref() const
Returns a reference to the dimensionless reference state Gibbs free energy vector.
void getChemPotentials(span< double > mu) const override
Get the species chemical potentials. Units: J/kmol.
void getStandardVolumes_ref(span< double > vol) const override
Get the molar volumes of the species reference states at the current T and P_ref of the solution.
virtual void update_single(size_t k, double T, double &cp_R, double &h_RT, double &s_R) const
Get reference-state properties for a single species.
An error indicating that an unimplemented function has been called.
ValueCache m_cache
Cached for saved calculations within each ThermoPhase.
Definition Phase.h:862
size_t m_kk
Number of species in the phase.
Definition Phase.h:882
size_t speciesIndex(const string &name, bool raise=true) const
Returns the index of a species named 'name' within the Phase object.
Definition Phase.cpp:127
double temperature() const
Temperature (K).
Definition Phase.h:586
virtual double concentration(const size_t k) const
Concentration of species k.
Definition Phase.cpp:495
double moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition Phase.cpp:457
double molecularWeight(size_t k) const
Molecular weight of species k.
Definition Phase.cpp:398
string name() const
Return the name of the phase.
Definition Phase.cpp:20
void checkElectronEnergyDistribution() const
Check the electron energy distribution.
void getStandardChemPotentials(span< double > muStar) const override
Get the array of chemical potentials at unit activity for the species at their standard states at the...
vector< vector< double > > m_energyLevels
Electron energy levels corresponding to the cross section data.
void setCollisions()
Set collisions.
double meanElectronEnergy() const
Mean electron energy [eV].
void getGibbs_ref(span< double > g) const override
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
double m_electronTempEquil
Saved electron temperature during an equilibrium solve.
double enthalpy_mole() const override
Return the Molar enthalpy. Units: J/kmol.
size_t m_nPoints
Number of points of electron energy levels.
void addCollision(shared_ptr< Reaction > collision)
Add a collision and record the target species.
virtual void setSolution(std::weak_ptr< Solution > soln) override
Set the link to the Solution object that owns this ThermoPhase.
void normalizeElectronEnergyDistribution()
Electron energy distribution norm.
void updateThermo() const override
Update the species reference state thermodynamic functions.
void getPartialMolarEnthalpies(span< double > hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
vector< size_t > m_targetSpeciesIndices
The collision-target species indices of m_collisions.
void setElectronTemperature(double Te) override
Set the internally stored electron temperature of the phase (K).
void electronEnergyLevelChanged()
When electron energy level changed, plasma properties such as electron-collision reaction rates need ...
double elasticPowerLoss()
The elastic power loss (J/s/m³)
int m_levelNum
Electron energy level change variable.
bool updateInterpolatedCrossSection(size_t k)
Update interpolated cross section of a collision.
bool m_inEquilibrate
Lock flag (default off)
void electronEnergyDistributionChanged()
When electron energy distribution changed, plasma properties such as electron-collision reaction rate...
size_t nElectronEnergyLevels() const
Number of electron levels.
size_t nCollisions() const
Number of electron collision cross sections.
void endEquilibrate() override
Hook called at the end of an equilibrium calculation on this phase.
Eigen::ArrayXd m_electronEnergyDist
Normalized electron energy distribution vector [-] Length: m_nPoints.
double electricField() const
Get the applied electric field strength [V/m].
Eigen::ArrayXd m_electronEnergyLevels
electron energy levels [ev]. Length: m_nPoints
void getGibbs_RT(span< double > grt) const override
Get the nondimensional Gibbs functions for the species in their standard states at the current T and ...
double intrinsicHeating() override
Intrinsic volumetric heating rate [W/m³].
double electronMobility() const
The electron mobility (m²/V/s)
void getParameters(AnyMap &phaseNode) const override
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
string type() const override
String indicating the thermodynamic model implemented.
void checkElectronEnergyLevels() const
Check the electron energy levels.
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
void updateElasticElectronEnergyLossCoefficients()
Update elastic electron energy loss coefficients.
void updateElectronTemperatureFromEnergyDist()
Update electron temperature (K) From energy distribution.
string m_distributionType
Electron energy distribution type.
void updateElectronEnergyDistribution()
Update the electron energy distribution.
vector< double > m_elasticElectronEnergyLossCoefficients
Elastic electron energy loss coefficients (eV m3/s)
string m_quadratureMethod
Numerical quadrature method for electron energy distribution.
PlasmaPhase(const string &inputFile="", const string &id="")
Construct and initialize a PlasmaPhase object directly from an input file.
void beginEquilibrate() override
Hook called at the beginning of an equilibrium calculation on this phase.
void getEntropy_R(span< double > sr) const override
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
void setDiscretizedElectronEnergyDist(span< const double > levels, span< const double > distrb)
Set discretized electron energy distribution.
double m_electronTemp
Electron temperature [K].
double RTe() const
Return the Gas Constant multiplied by the current electron temperature.
double intEnergy_mole() const override
Return the molar internal energy. Units: J/kmol.
double entropy_mole() const override
Return the molar entropy. Units: J/kmol/K.
bool m_do_normalizeElectronEnergyDist
Flag of normalizing electron energy distribution.
void updateElectronEnergyDistDifference()
Update electron energy distribution difference.
void updateElasticElectronEnergyLossCoefficient(size_t i)
Updates the elastic electron energy loss coefficient for collision index i.
vector< size_t > m_kElastic
Indices of elastic collisions in m_crossSections.
unique_ptr< EEDFTwoTermApproximation > m_eedfSolver
Solver used to calculate the EEDF based on electron collision rates.
string electronSpeciesName() const
Electron species name.
void setIsotropicElectronEnergyDistribution()
Set isotropic electron energy distribution.
Eigen::ArrayXd m_electronEnergyDistDiff
ionization degree for the electron-electron collisions (tmp is the previous one)
void getPartialMolarEntropies(span< double > sbar) const override
Returns an array of partial molar entropies of the species in the solution.
double gibbs_mole() const override
Return the molar Gibbs free energy. Units: J/kmol.
std::vector< double > m_work
Work array.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
const shared_ptr< Reaction > collision(size_t i) const
Get the Reaction object associated with electron collision i.
vector< bool > m_interp_cs_ready
The list of whether the interpolated cross sections is ready.
vector< shared_ptr< ElectronCollisionPlasmaRate > > m_collisionRates
The list of shared pointers of collision rates.
void getChemPotentials(span< double > mu) const override
Get the species chemical potentials. Units: J/kmol.
void setElectronEnergyLevels(span< const double > levels)
Set electron energy levels.
vector< shared_ptr< Reaction > > m_collisions
The list of shared pointers of plasma collision reactions.
void setParameters(const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap()) override
Set equation of state parameters from an AnyMap phase description.
void setMeanElectronEnergy(double energy)
Set mean electron energy [eV].
void getStandardVolumes_ref(span< double > vol) const override
Get the molar volumes of the species reference states at the current T and P_ref of the solution.
size_t m_electronSpeciesIndex
Index of electron species.
vector< vector< double > > m_crossSections
Cross section data.
void setElectronEnergyDistributionType(const string &type)
Set electron energy distribution type.
virtual double electronPressure() const
Electron pressure.
double jouleHeatingPower() const
The joule heating power (W/m³)
vector< size_t > m_kInelastic
Indices of inelastic collisions in m_crossSections.
double electronTemperature() const override
Electron Temperature (K)
void setIsotropicShapeFactor(double x)
Set the shape factor of isotropic electron energy distribution.
void enableNormalizeElectronEnergyDist(bool enable)
Set flag of automatically normalize electron energy distribution Flag: m_do_normalizeElectronEnergyDi...
void getPartialMolarIntEnergies(span< double > ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
int m_distNum
Electron energy distribution change variable.
virtual void endEquilibrate()
Hook called at the end of an equilibrium calculation on this phase.
virtual void setParameters(const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap())
Set equation of state parameters from an AnyMap phase description.
virtual void getParameters(AnyMap &phaseNode) const
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
double RT() const
Return the Gas Constant multiplied by the current temperature.
virtual void setSolution(std::weak_ptr< Solution > soln)
Set the link to the Solution object that owns this ThermoPhase.
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
void initThermoFile(const string &inputFile, const string &id)
Initialize a ThermoPhase object using an input file.
std::weak_ptr< Solution > m_soln
reference to Solution
virtual void beginEquilibrate()
Hook called at the beginning of an equilibrium calculation on this phase.
MultiSpeciesThermo m_spthermo
Pointer to the calculation manager for species reference-state thermodynamic properties.
virtual double refPressure() const
Returns the reference pressure in Pa.
CachedScalar getScalar(int id)
Get a reference to a CachedValue object representing a scalar (double) with the given id.
Definition ValueCache.h:161
int getId()
Get a unique id for a cached value.
Header for a file containing miscellaneous numerical functions.
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition global.h:171
double numericalQuadrature(const string &method, const Eigen::ArrayXd &f, const Eigen::ArrayXd &x)
Numerical integration of a function.
Definition funcs.cpp:116
const double Boltzmann
Boltzmann constant [J/K].
Definition ct_defs.h:87
const double Avogadro
Avogadro's Number [number/kmol].
Definition ct_defs.h:84
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
const double ElectronCharge
Elementary charge [C].
Definition ct_defs.h:93
const double ElectronMass
Electron Mass [kg].
Definition ct_defs.h:114
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition global.h:263
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:183
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
MappedVector asVectorXd(vector< double > &v)
Convenience wrapper for accessing std::vector as an Eigen VectorXd.
Definition eigen_dense.h:60
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
Definition eigen_dense.h:46
const double SmallNumber
smallest number to compare to zero.
Definition ct_defs.h:161
map< string, double > Composition
Map from string names to doubles.
Definition ct_defs.h:180
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
A cached property value and the state at which it was evaluated.
Definition ValueCache.h:33
double state2
Value of the second state variable for the state at which value was evaluated, for example density or...
Definition ValueCache.h:106
bool validate(double state1New)
Check whether the currently cached value is valid based on a single state variable.
Definition ValueCache.h:39
double state1
Value of the first state variable for the state at which value was evaluated, for example temperature...
Definition ValueCache.h:102