Cantera  3.1.0a2
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
7#include <boost/math/special_functions/gamma.hpp>
11
12namespace Cantera {
13
14PlasmaPhase::PlasmaPhase(const string& inputFile, const string& id_)
15{
16 initThermoFile(inputFile, id_);
17
18 // initial grid
19 m_electronEnergyLevels = Eigen::ArrayXd::LinSpaced(m_nPoints, 0.0, 1.0);
20
21 // initial electron temperature
23}
24
26{
27 if (m_distributionType == "discretized") {
28 throw CanteraError("PlasmaPhase::updateElectronEnergyDistribution",
29 "Invalid for discretized electron energy distribution.");
30 } else if (m_distributionType == "isotropic") {
32 }
34}
35
37 Eigen::ArrayXd eps32 = m_electronEnergyLevels.pow(3./2.);
38 double norm = 2./3. * numericalQuadrature(m_quadratureMethod,
40 if (norm < 0.0) {
41 throw CanteraError("PlasmaPhase::normalizeElectronEnergyDistribution",
42 "The norm is negative. This might be caused by bad "
43 "electron energy distribution");
44 }
46}
47
49{
50 if (type == "discretized" ||
51 type == "isotropic") {
53 } else {
54 throw CanteraError("PlasmaPhase::setElectronEnergyDistributionType",
55 "Unknown type for electron energy distribution.");
56 }
57}
58
60{
62 double x = m_isotropicShapeFactor;
63 double gamma1 = boost::math::tgamma(3.0 / 2.0 * x);
64 double gamma2 = boost::math::tgamma(5.0 / 2.0 * x);
65 double c1 = x * std::pow(gamma2, 1.5) / std::pow(gamma1, 2.5);
66 double c2 = x * std::pow(gamma2 / gamma1, x);
68 c1 * m_electronEnergyLevels.sqrt() /
69 std::pow(meanElectronEnergy(), 1.5) *
71 meanElectronEnergy()).pow(x)).exp();
73}
74
76 m_electronTemp = Te;
78}
79
81 m_electronTemp = 2.0 / 3.0 * energy * ElectronCharge / Boltzmann;
83}
84
85void PlasmaPhase::setElectronEnergyLevels(const double* levels, size_t length)
86{
87 m_nPoints = length;
88 m_electronEnergyLevels = Eigen::Map<const Eigen::ArrayXd>(levels, length);
92}
93
95{
96 m_distNum++;
97}
98
100{
101 m_levelNum++;
102}
103
105{
106 Eigen::ArrayXd h = m_electronEnergyLevels.tail(m_nPoints - 1) -
108 if (m_electronEnergyLevels[0] < 0.0 || (h <= 0.0).any()) {
109 throw CanteraError("PlasmaPhase::checkElectronEnergyLevels",
110 "Values of electron energy levels need to be positive and "
111 "monotonically increasing.");
112 }
113}
114
116{
117 Eigen::ArrayXd h = m_electronEnergyLevels.tail(m_nPoints - 1) -
119 if ((m_electronEnergyDist < 0.0).any()) {
120 throw CanteraError("PlasmaPhase::checkElectronEnergyDistribution",
121 "Values of electron energy distribution cannot be negative.");
122 }
123 if (m_electronEnergyDist[m_nPoints - 1] > 0.01) {
124 warn_user("PlasmaPhase::checkElectronEnergyDistribution",
125 "The value of the last element of electron energy distribution exceed 0.01. "
126 "This indicates that the value of electron energy level is not high enough "
127 "to contain the isotropic distribution at mean electron energy of "
128 "{} eV", meanElectronEnergy());
129 }
130}
131
133 const double* dist,
134 size_t length)
135{
136 m_distributionType = "discretized";
137 m_nPoints = length;
139 Eigen::Map<const Eigen::ArrayXd>(levels, length);
141 Eigen::Map<const Eigen::ArrayXd>(dist, length);
145 }
150}
151
153 size_t length)
154{
155 m_distributionType = "discretized";
156 m_nPoints = length;
158 Eigen::Map<const Eigen::ArrayXd>(dist, length);
162 }
166}
167
169{
170 // calculate mean electron energy and electron temperature
171 Eigen::ArrayXd eps52 = m_electronEnergyLevels.pow(5./2.);
172 double epsilon_m = 2.0 / 5.0 * numericalQuadrature(m_quadratureMethod,
173 m_electronEnergyDist, eps52);
174 m_electronTemp = 2.0 / 3.0 * epsilon_m * ElectronCharge / Boltzmann;
175}
176
178 m_isotropicShapeFactor = x;
180}
181
182void PlasmaPhase::getParameters(AnyMap& phaseNode) const
183{
185 AnyMap eedf;
186 eedf["type"] = m_distributionType;
187 vector<double> levels(m_nPoints);
188 Eigen::Map<Eigen::ArrayXd>(levels.data(), m_nPoints) = m_electronEnergyLevels;
189 eedf["energy-levels"] = levels;
190 if (m_distributionType == "isotropic") {
191 eedf["shape-factor"] = m_isotropicShapeFactor;
192 eedf["mean-electron-energy"].setQuantity(meanElectronEnergy(), "eV");
193 } else if (m_distributionType == "discretized") {
194 vector<double> dist(m_nPoints);
195 Eigen::Map<Eigen::ArrayXd>(dist.data(), m_nPoints) = m_electronEnergyDist;
196 eedf["distribution"] = dist;
197 eedf["normalize"] = m_do_normalizeElectronEnergyDist;
198 }
199 phaseNode["electron-energy-distribution"] = std::move(eedf);
200}
201
202void PlasmaPhase::setParameters(const AnyMap& phaseNode, const AnyMap& rootNode)
203{
204 IdealGasPhase::setParameters(phaseNode, rootNode);
205 if (phaseNode.hasKey("electron-energy-distribution")) {
206 const AnyMap eedf = phaseNode["electron-energy-distribution"].as<AnyMap>();
207 m_distributionType = eedf["type"].asString();
208 if (m_distributionType == "isotropic") {
209 if (eedf.hasKey("shape-factor")) {
210 setIsotropicShapeFactor(eedf["shape-factor"].asDouble());
211 } else {
212 throw CanteraError("PlasmaPhase::setParameters",
213 "isotropic type requires shape-factor key.");
214 }
215 if (eedf.hasKey("mean-electron-energy")) {
216 double energy = eedf.convert("mean-electron-energy", "eV");
217 setMeanElectronEnergy(energy);
218 } else {
219 throw CanteraError("PlasmaPhase::setParameters",
220 "isotropic type requires electron-temperature key.");
221 }
222 if (eedf.hasKey("energy-levels")) {
223 setElectronEnergyLevels(eedf["energy-levels"].asVector<double>().data(),
224 eedf["energy-levels"].asVector<double>().size());
225 }
227 } else if (m_distributionType == "discretized") {
228 if (!eedf.hasKey("energy-levels")) {
229 throw CanteraError("PlasmaPhase::setParameters",
230 "Cannot find key energy-levels.");
231 }
232 if (!eedf.hasKey("distribution")) {
233 throw CanteraError("PlasmaPhase::setParameters",
234 "Cannot find key distribution.");
235 }
236 if (eedf.hasKey("normalize")) {
237 enableNormalizeElectronEnergyDist(eedf["normalize"].asBool());
238 }
239 setDiscretizedElectronEnergyDist(eedf["energy-levels"].asVector<double>().data(),
240 eedf["distribution"].asVector<double>().data(),
241 eedf["energy-levels"].asVector<double>().size());
242 }
243 }
244}
245
246bool PlasmaPhase::addSpecies(shared_ptr<Species> spec)
247{
248 bool added = IdealGasPhase::addSpecies(spec);
249 size_t k = m_kk - 1;
250
251 if (spec->composition.find("E") != spec->composition.end() &&
252 spec->composition.size() == 1 &&
253 spec->composition["E"] == 1) {
256 } else {
257 throw CanteraError("PlasmaPhase::addSpecies",
258 "Cannot add species, {}. "
259 "Only one electron species is allowed.", spec->name);
260 }
261 }
262 return added;
263}
264
266{
268 // check electron species
270 throw CanteraError("PlasmaPhase::initThermo",
271 "No electron species found.");
272 }
273}
274
276{
278 static const int cacheId = m_cache.getId();
279 CachedScalar cached = m_cache.getScalar(cacheId);
280 double tempNow = temperature();
281 double electronTempNow = electronTemperature();
282 size_t k = m_electronSpeciesIndex;
283 // If the electron temperature has changed since the last time these
284 // properties were computed, recompute them.
285 if (cached.state1 != tempNow || cached.state2 != electronTempNow) {
287 &m_cp0_R[k], &m_h0_RT[k], &m_s0_R[k]);
288 cached.state1 = tempNow;
289 cached.state2 = electronTempNow;
290 }
291 // update the species Gibbs functions
292 m_g0_RT[k] = m_h0_RT[k] - m_s0_R[k];
293}
294
296 double value = IdealGasPhase::enthalpy_mole();
297 value += GasConstant * (electronTemperature() - temperature()) *
300 return value;
301}
302
303void PlasmaPhase::getGibbs_ref(double* g) const
304{
307}
308
310{
313}
314
316{
319}
320
322{
324 double logp = log(pressure());
325 double logpe = log(electronPressure());
326 sbar[m_electronSpeciesIndex] += GasConstant * (logp - logpe);
327}
328
330{
331 const vector<double>& _h = enthalpy_RT_ref();
332 for (size_t k = 0; k < m_kk; k++) {
333 ubar[k] = RT() * (_h[k] - 1.0);
334 }
335 size_t k = m_electronSpeciesIndex;
336 ubar[k] = RTe() * (_h[k] - 1.0);
337}
338
339void PlasmaPhase::getChemPotentials(double* mu) const
340{
342 size_t k = m_electronSpeciesIndex;
343 double xx = std::max(SmallNumber, moleFraction(k));
344 mu[k] += (RTe() - RT()) * log(xx);
345}
346
348{
350 size_t k = m_electronSpeciesIndex;
351 muStar[k] -= log(pressure() / refPressure()) * RT();
352 muStar[k] += log(electronPressure() / refPressure()) * RTe();
353}
354
355void PlasmaPhase::getEntropy_R(double* sr) const
356{
357 const vector<double>& _s = entropy_R_ref();
358 copy(_s.begin(), _s.end(), sr);
359 double tmp = log(pressure() / refPressure());
360 for (size_t k = 0; k < m_kk; k++) {
361 if (k != m_electronSpeciesIndex) {
362 sr[k] -= tmp;
363 } else {
364 sr[k] -= log(electronPressure() / refPressure());
365 }
366 }
367}
368
369void PlasmaPhase::getGibbs_RT(double* grt) const
370{
371 const vector<double>& gibbsrt = gibbs_RT_ref();
372 copy(gibbsrt.begin(), gibbsrt.end(), grt);
373 double tmp = log(pressure() / refPressure());
374 for (size_t k = 0; k < m_kk; k++) {
375 if (k != m_electronSpeciesIndex) {
376 grt[k] += tmp;
377 } else {
378 grt[k] += log(electronPressure() / refPressure());
379 }
380 }
381}
382
383}
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:427
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
Definition AnyMap.cpp:1423
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:1535
Base class for exceptions thrown by Cantera classes.
const vector< double > & entropy_R_ref() const
Returns a reference to the dimensionless reference state Entropy vector.
double enthalpy_mole() const override
Return the Molar enthalpy. Units: J/kmol.
void getPartialMolarEnthalpies(double *hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
void getChemPotentials(double *mu) const override
Get the species chemical potentials. Units: J/kmol.
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.
void getGibbs_ref(double *g) const override
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
void getStandardChemPotentials(double *mu) const override
Get the array of chemical potentials at unit activity for the species at their standard states at the...
const vector< double > & gibbs_RT_ref() const
Returns a reference to the dimensionless reference state Gibbs free energy vector.
void getStandardVolumes_ref(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 updateThermo() const
Update the species reference state thermodynamic functions.
vector< double > m_s0_R
Temporary storage for dimensionless reference state entropies.
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.
void getPartialMolarEntropies(double *sbar) const override
Returns an array of partial molar entropies of the species in the solution.
const vector< double > & enthalpy_RT_ref() const
Returns a reference to the dimensionless reference state enthalpy vector.
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.
ValueCache m_cache
Cached for saved calculations within each ThermoPhase.
Definition Phase.h:834
size_t m_kk
Number of species in the phase.
Definition Phase.h:854
double temperature() const
Temperature (K).
Definition Phase.h:562
double moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition Phase.cpp:439
void checkElectronEnergyDistribution() const
Check the electron energy distribution.
double meanElectronEnergy() const
Mean electron energy [eV].
double enthalpy_mole() const override
Return the Molar enthalpy. Units: J/kmol.
size_t m_nPoints
Number of points of electron energy levels.
void getPartialMolarEnthalpies(double *hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
void getChemPotentials(double *mu) const override
Get the species chemical potentials. Units: J/kmol.
void normalizeElectronEnergyDistribution()
Electron energy distribution norm.
void getStandardChemPotentials(double *muStar) const override
Get the array of chemical potentials at unit activity for the species at their standard states at the...
void updateThermo() const override
Update the species reference state thermodynamic functions.
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 ...
int m_levelNum
Electron energy level change variable.
void electronEnergyDistributionChanged()
When electron energy distribution changed, plasma properties such as electron-collision reaction rate...
void getEntropy_R(double *sr) const override
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
void getGibbs_ref(double *g) const override
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
Eigen::ArrayXd m_electronEnergyDist
Normalized electron energy distribution vector [-] Length: m_nPoints.
Eigen::ArrayXd m_electronEnergyLevels
electron energy levels [ev]. Length: m_nPoints
void setDiscretizedElectronEnergyDist(const double *levels, const double *distrb, size_t length)
Set discretized electron energy distribution.
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.
Definition PlasmaPhase.h:74
void checkElectronEnergyLevels() const
Check the electron energy levels.
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
void updateElectronTemperatureFromEnergyDist()
Update electron temperature (K) From energy distribution.
string m_distributionType
Electron energy distribution type.
void getStandardVolumes_ref(double *vol) const override
Get the molar volumes of the species reference states at the current T and P_ref of the solution.
void updateElectronEnergyDistribution()
Update electron energy distribution.
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.
double m_electronTemp
Electron temperature [K].
double RTe() const
Return the Gas Constant multiplied by the current electron temperature.
void setElectronEnergyLevels(const double *levels, size_t length)
Set electron energy levels.
void getGibbs_RT(double *grt) const override
Get the nondimensional Gibbs functions for the species in their standard states at the current T and ...
bool m_do_normalizeElectronEnergyDist
Flag of normalizing electron energy distribution.
void getPartialMolarIntEnergies(double *ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
void setIsotropicElectronEnergyDistribution()
Set isotropic electron energy distribution.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
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].
size_t m_electronSpeciesIndex
Index of electron species.
void setElectronEnergyDistributionType(const string &type)
Set electron energy distribution type.
void getPartialMolarEntropies(double *sbar) const override
Returns an array of partial molar entropies of the species in the solution.
virtual double electronPressure() const
Electron pressure.
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...
int m_distNum
Electron energy distribution change variable.
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 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.
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,...
double numericalQuadrature(const string &method, const Eigen::ArrayXd &f, const Eigen::ArrayXd &x)
Numerical integration of a function.
Definition funcs.cpp:112
const double Boltzmann
Boltzmann constant [J/K].
Definition ct_defs.h:84
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
const double ElectronCharge
Elementary charge [C].
Definition ct_defs.h:90
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition global.h:267
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
const double SmallNumber
smallest number to compare to zero.
Definition ct_defs.h:158
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
double state1
Value of the first state variable for the state at which value was evaluated, for example temperature...
Definition ValueCache.h:102