8#include <boost/math/special_functions/gamma.hpp>
11#include "cantera/numerics/eigen_dense.h"
15#include <boost/polymorphic_pointer_cast.hpp>
32 m_eedfSolver = make_unique<EEDFTwoTermApproximation>(
this);
36 size_t nGridCells = 301;
43PlasmaPhase::~PlasmaPhase()
45 if (shared_ptr<Solution> soln =
m_soln.lock()) {
46 soln->removeChangedCallback(
this);
47 soln->kinetics()->removeReactionAddedCallback(
this);
58 throw CanteraError(
"PlasmaPhase::updateElectronEnergyDistribution",
59 "Invalid for discretized electron energy distribution.");
63 auto ierr =
m_eedfSolver->calculateDistributionFunction();
68 throw CanteraError(
"PlasmaPhase::updateElectronEnergyDistribution",
69 "Call to calculateDistributionFunction failed.");
81 writelog(
"Skipping Te update: EEDF is empty, non-finite, or unnormalized.\n");
84 throw CanteraError(
"PlasmaPhase::updateElectronEnergyDistribution",
96 throw CanteraError(
"PlasmaPhase::normalizeElectronEnergyDistribution",
97 "The norm is negative. This might be caused by bad "
98 "electron energy distribution");
105 if (
type ==
"discretized" ||
106 type ==
"isotropic" ||
107 type ==
"Boltzmann-two-term") {
110 throw CanteraError(
"PlasmaPhase::setElectronEnergyDistributionType",
111 "Unknown type for electron energy distribution.");
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);
185 const auto& rate = boost::polymorphic_pointer_downcast
197 throw CanteraError(
"PlasmaPhase::checkElectronEnergyLevels",
198 "Values of electron energy levels need to be positive and "
199 "monotonically increasing.");
208 throw CanteraError(
"PlasmaPhase::checkElectronEnergyDistribution",
209 "Values of electron energy distribution cannot be negative.");
212 warn_user(
"PlasmaPhase::checkElectronEnergyDistribution",
213 "The value of the last element of electron energy distribution exceed 0.01. "
214 "This indicates that the value of electron energy level is not high enough "
215 "to contain the isotropic distribution at mean electron energy of "
221 span<const double> dist)
226 Eigen::Map<const Eigen::ArrayXd>(levels.data(),
m_nPoints);
228 Eigen::Map<const Eigen::ArrayXd>(dist.data(),
m_nPoints);
252 if (epsilon_m < 0.0) {
253 throw CanteraError(
"PlasmaPhase::updateElectronTemperatureFromEnergyDist",
254 "The electron energy distribution produces negative electron temperature.");
261 m_isotropicShapeFactor = x;
272 eedf[
"energy-levels"] = levels;
274 eedf[
"shape-factor"] = m_isotropicShapeFactor;
279 eedf[
"distribution"] = dist;
282 phaseNode[
"electron-energy-distribution"] = std::move(eedf);
288 if (phaseNode.
hasKey(
"electron-energy-distribution")) {
289 const AnyMap eedf = phaseNode[
"electron-energy-distribution"].as<
AnyMap>();
292 if (eedf.
hasKey(
"shape-factor")) {
296 "isotropic type requires shape-factor key.");
298 if (eedf.
hasKey(
"mean-electron-energy")) {
299 double energy = eedf.
convert(
"mean-electron-energy",
"eV");
303 "isotropic type requires electron-temperature key.");
305 if (eedf.
hasKey(
"energy-levels")) {
306 auto levels = eedf[
"energy-levels"].asVector<
double>();
311 if (!eedf.
hasKey(
"energy-levels")) {
313 "Cannot find key energy-levels.");
315 if (!eedf.
hasKey(
"distribution")) {
317 "Cannot find key distribution.");
319 if (eedf.
hasKey(
"normalize")) {
322 auto levels = eedf[
"energy-levels"].asVector<
double>();
323 auto distribution = eedf[
"distribution"].asVector<
double>(levels.size());
328 if (rootNode.
hasKey(
"electron-collisions")) {
329 for (
const auto& item : rootNode[
"electron-collisions"].asVector<
AnyMap>()) {
330 auto rate = make_shared<ElectronCollisionPlasmaRate>(item);
332 reactants[item[
"target"].asString()] = 1;
334 if (item.hasKey(
"product")) {
335 products[item[
"product"].asString()] = 1;
337 products[item[
"target"].asString()] = 1;
340 if (rate->kind() ==
"ionization") {
342 }
else if (rate->kind() ==
"attachment") {
345 auto R = make_shared<Reaction>(reactants, products, rate);
356 if ((spec->name ==
"e" || spec->name ==
"Electron") ||
357 (spec->composition.find(
"E") != spec->composition.end() &&
358 spec->composition.size() == 1 &&
359 spec->composition[
"E"] == 1)) {
364 "Cannot add species, {}. "
365 "Only one electron species is allowed.", spec->name);
378 "No electron species found.");
386 if (shared_ptr<Solution> soln =
m_soln.lock()) {
387 soln->registerChangedCallback(
this, [&]() {
395 if (shared_ptr<Solution> soln =
m_soln.lock()) {
396 shared_ptr<Kinetics> kin = soln->kinetics();
403 set<Reaction*> existing;
405 existing.insert(R.get());
407 for (
size_t i = 0; i < kin->nReactions(); i++) {
408 shared_ptr<Reaction> R = kin->reaction(i);
409 if (R->rate()->type() !=
"electron-collision-plasma"
410 || existing.count(R.get())) {
418 kin->registerReactionAddedCallback(
this, [
this, kin]() {
419 size_t i = kin->nReactions() - 1;
420 if (kin->reaction(i)->type() ==
"electron-collision-plasma") {
421 addCollision(kin->reaction(i));
436 std::dynamic_pointer_cast<ElectronCollisionPlasmaRate>(
collision->rate());
449 if (target.empty()) {
450 throw CanteraError(
"PlasmaPhase::addCollision",
"Error identifying target for"
451 " collision with equation '{}'",
collision->equation());
456 std::dynamic_pointer_cast<ElectronCollisionPlasmaRate>(
collision->rate()));
467 if ((kind ==
"effective" || kind ==
"elastic")) {
473 throw CanteraError(
"PlasmaPhase::addCollision",
"Phase already contains"
474 " an effective/elastic cross section for '{}'.", target);
482 auto levels = rate.energyLevels();
484 auto sections = rate.crossSections();
510 for (
size_t i = 1; i <
m_nPoints - 1; i++) {
516 (h1 * h0) / (h1 + h0);
545 if (interpChanged[i]) {
564 auto cs_array = Eigen::Map<const Eigen::ArrayXd>(
586 "EEDF not initialized");
600 if (!std::isfinite(q_elastic)) {
602 "Non-finite elastic power loss");
618 if (cached.
state1 != tempNow || cached.
state2 != electronTempNow) {
622 cached.
state2 = electronTempNow;
641 for (
size_t k = 0; k <
m_kk; ++k) {
652 for (
size_t k = 0; k <
m_kk; ++k) {
663 for (
size_t k = 0; k <
m_kk; ++k) {
699 for (
size_t k = 0; k <
m_kk; k++) {
700 ubar[k] =
RT() * (_h[k] - 1.0);
703 ubar[k] =
RTe() * (_h[k] - 1.0);
711 mu[k] += (
RTe() -
RT()) * log(xx);
726 copy(_s.begin(), _s.end(), sr.begin());
728 for (
size_t k = 0; k <
m_kk; k++) {
741 copy(gibbsrt.begin(), gibbsrt.end(), grt.begin());
743 for (
size_t k = 0; k <
m_kk; k++) {
759 "Electron mobility is only available for 'Boltzmann-two-term' "
760 "electron energy distributions.");
781 return sigma * E * E;
794 return qJ + qElastic;
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.
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
double convert(const string &key, const string &units) const
Convert the item stored by the given key to the units specified in units.
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.
size_t m_kk
Number of species in the phase.
size_t speciesIndex(const string &name, bool raise=true) const
Returns the index of a species named 'name' within the Phase object.
double temperature() const
Temperature (K).
virtual double concentration(const size_t k) const
Concentration of species k.
double moleFraction(size_t k) const
Return the mole fraction of a single species.
double molecularWeight(size_t k) const
Molecular weight of species k.
string name() const
Return the name of the phase.
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.
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.
double numericalQuadrature(const string &method, const Eigen::ArrayXd &f, const Eigen::ArrayXd &x)
Numerical integration of a function.
const double Boltzmann
Boltzmann constant [J/K].
const double Avogadro
Avogadro's Number [number/kmol].
const double GasConstant
Universal Gas Constant [J/kmol/K].
const double ElectronCharge
Elementary charge [C].
const double ElectronMass
Electron Mass [kg].
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Namespace for the Cantera kernel.
const size_t npos
index returned by functions to indicate "no position"
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
const double SmallNumber
smallest number to compare to zero.
map< string, double > Composition
Map from string names to doubles.
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.
double state2
Value of the second state variable for the state at which value was evaluated, for example density or...
bool validate(double state1New)
Check whether the currently cached value is valid based on a single state variable.
double state1
Value of the first state variable for the state at which value was evaluated, for example temperature...