Cantera  3.1.0a2
Loading...
Searching...
No Matches
PlasmaPhase.h
Go to the documentation of this file.
1/**
2 * @file PlasmaPhase.h
3 * Header file for class PlasmaPhase.
4 */
5
6// This file is part of Cantera. See License.txt in the top-level directory or
7// at https://cantera.org/license.txt for license and copyright information.
8
9#ifndef CT_PLASMAPHASE_H
10#define CT_PLASMAPHASE_H
11
13#include "cantera/numerics/eigen_sparse.h"
14
15namespace Cantera
16{
17
18/**
19 * Base class for a phase with plasma properties. This class manages the
20 * plasma properties such as electron energy distribution function (EEDF).
21 * There are two ways to define the electron distribution and electron
22 * temperature. The first method uses setElectronTemperature() to set
23 * the electron temperature which is used to calculate the electron energy
24 * distribution with isotropic-velocity model. The generalized electron
25 * energy distribution for isotropic-velocity distribution can be
26 * expressed as [1,2],
27 * @f[
28 * f(\epsilon) = c_1 \frac{\sqrt{\epsilon}}{\epsilon_m^{3/2}}
29 * \exp(-c_2 (\frac{\epsilon}{\epsilon_m})^x),
30 * @f]
31 * where @f$ x = 1 @f$ and @f$ x = 2 @f$ correspond to the Maxwellian and
32 * Druyvesteyn (default) electron energy distribution, respectively.
33 * @f$ \epsilon_m = 3/2 T_e @f$ [eV] (mean electron energy). The second
34 * method uses setDiscretizedElectronEnergyDist() to manually set electron
35 * energy distribution and calculate electron temperature from mean electron
36 * energy, which is calculated as [3],
37 * @f[
38 * \epsilon_m = \int_0^{\infty} \epsilon^{3/2} f(\epsilon) d\epsilon,
39 * @f]
40 * which can be calculated using trapezoidal rule,
41 * @f[
42 * \epsilon_m = \sum_i (\epsilon^{5/2}_{i+1} - \epsilon^{5/2}_i)
43 * (f(\epsilon_{i+1}) + f(\epsilon_i)) / 2,
44 * @f]
45 * where @f$ i @f$ is the index of energy levels.
46 *
47 * For references, see Gudmundsson @cite gudmundsson2001; Khalilpour and Foroutan
48 * @cite khalilpour2020; Hagelaar and Pitchford @cite hagelaar2005, and BOLOS
49 * @cite BOLOS.
50 *
51 * @warning This class is an experimental part of %Cantera and may be
52 * changed or removed without notice.
53 * @todo Implement electron Boltzmann equation solver to solve EEDF.
54 * https://github.com/Cantera/enhancements/issues/127
55 * @ingroup phase
56 */
58{
59public:
60 //! Construct and initialize a PlasmaPhase object
61 //! directly from an input file. The constructor initializes the electron
62 //! energy distribution to be Druyvesteyn distribution (m_x = 2.0). The initial
63 //! electron energy grid is set to a linear space which starts at 0.01 eV and ends
64 //! at 1 eV with 1000 points.
65 /*!
66 * @param inputFile Name of the input file containing the phase definition
67 * to set up the object. If blank, an empty phase will be
68 * created.
69 * @param id ID of the phase in the input file. Defaults to the
70 * empty string.
71 */
72 explicit PlasmaPhase(const string& inputFile="", const string& id="");
73
74 string type() const override {
75 return "plasma";
76 }
77
78 void initThermo() override;
79
80 //! Set electron energy levels.
81 //! @param levels The vector of electron energy levels (eV).
82 //! Length: #m_nPoints.
83 //! @param length The length of the @c levels.
84 void setElectronEnergyLevels(const double* levels, size_t length);
85
86 //! Get electron energy levels.
87 //! @param levels The vector of electron energy levels (eV). Length: #m_nPoints
88 void getElectronEnergyLevels(double* levels) const {
89 Eigen::Map<Eigen::ArrayXd>(levels, m_nPoints) = m_electronEnergyLevels;
90 }
91
92 //! Set discretized electron energy distribution.
93 //! @param levels The vector of electron energy levels (eV).
94 //! Length: #m_nPoints.
95 //! @param distrb The vector of electron energy distribution.
96 //! Length: #m_nPoints.
97 //! @param length The length of the vectors, which equals #m_nPoints.
98 void setDiscretizedElectronEnergyDist(const double* levels,
99 const double* distrb,
100 size_t length);
101
102 //! Set discretized electron energy distribution.
103 //! @param distrb The vector of electron energy distribution.
104 //! Length: #m_nPoints.
105 //! @param length The length of the vectors, which equals #m_nPoints.
106 void setDiscretizedElectronEnergyDist(const double* distrb,
107 size_t length);
108
109 //! Get electron energy distribution.
110 //! @param distrb The vector of electron energy distribution.
111 //! Length: #m_nPoints.
112 void getElectronEnergyDistribution(double* distrb) const {
113 Eigen::Map<Eigen::ArrayXd>(distrb, m_nPoints) = m_electronEnergyDist;
114 }
115
116 //! Set the shape factor of isotropic electron energy distribution.
117 //! Note that @f$ x = 1 @f$ and @f$ x = 2 @f$ correspond to the
118 //! Maxwellian and Druyvesteyn distribution, respectively.
119 //! @param x The shape factor
120 void setIsotropicShapeFactor(double x);
121
122 //! The shape factor of isotropic electron energy distribution
123 double isotropicShapeFactor() const {
124 return m_isotropicShapeFactor;
125 }
126
127 //! Set the internally stored electron temperature of the phase (K).
128 //! @param Te Electron temperature in Kelvin
129 void setElectronTemperature(double Te) override;
130
131 //! Set mean electron energy [eV]. This method also sets electron temperature
132 //! accordingly.
133 void setMeanElectronEnergy(double energy);
134
135 //! Get electron energy distribution type
137 return m_distributionType;
138 }
139
140 //! Set electron energy distribution type
141 void setElectronEnergyDistributionType(const string& type);
142
143 //! Numerical quadrature method. Method: #m_quadratureMethod
144 string quadratureMethod() const {
145 return m_quadratureMethod;
146 }
147
148 //! Set numerical quadrature method for integrating electron
149 //! energy distribution function. Method: #m_quadratureMethod
150 void setQuadratureMethod(const string& method) {
151 m_quadratureMethod = method;
152 }
153
154 //! Mean electron energy [eV]
155 double meanElectronEnergy() const {
156 return 3.0 / 2.0 * electronTemperature() * Boltzmann / ElectronCharge;
157 }
158
159 //! Set flag of automatically normalize electron energy distribution
160 //! Flag: #m_do_normalizeElectronEnergyDist
163 }
164
165 //! Flag of automatically normalize electron energy distribution.
166 //! Flag: #m_do_normalizeElectronEnergyDist
169 }
170
171 bool addSpecies(shared_ptr<Species> spec) override;
172
173 //! Electron Temperature (K)
174 //! @return The electron temperature of the phase
175 double electronTemperature() const override {
176 return m_electronTemp;
177 }
178
179 //! Return the Gas Constant multiplied by the current electron temperature
180 /*!
181 * The units are Joules kmol-1
182 */
183 double RTe() const {
185 }
186
187 /**
188 * Electron pressure. Units: Pa.
189 * @f[P = n_{k_e} R T_e @f]
190 */
191 virtual double electronPressure() const {
194 }
195
196 //! Number of electron levels
197 size_t nElectronEnergyLevels() const {
198 return m_nPoints;
199 }
200
201 //! Electron Species Index
202 size_t electronSpeciesIndex() const {
204 }
205
206 //! Return the Molar enthalpy. Units: J/kmol.
207 /*!
208 * For an ideal gas mixture with additional electron,
209 * @f[
210 * \hat h(T) = \sum_{k \neq k_e} X_k \hat h^0_k(T) + X_{k_e} \hat h^0_{k_e}(T_e),
211 * @f]
212 * and is a function only of temperature. The standard-state pure-species
213 * enthalpies @f$ \hat h^0_k(T) @f$ are computed by the species
214 * thermodynamic property manager.
215 *
216 * \see MultiSpeciesThermo
217 */
218 double enthalpy_mole() const override;
219
220 double cp_mole() const override {
221 throw NotImplementedError("PlasmaPhase::cp_mole");
222 }
223
224 double entropy_mole() const override {
225 throw NotImplementedError("PlasmaPhase::entropy_mole");
226 }
227
228 double gibbs_mole() const override {
229 throw NotImplementedError("PlasmaPhase::gibbs_mole");
230 }
231
232 double intEnergy_mole() const override {
233 throw NotImplementedError("PlasmaPhase::intEnergy_mole");
234 }
235
236 void getEntropy_R(double* sr) const override;
237
238 void getGibbs_RT(double* grt) const override;
239
240 void getGibbs_ref(double* g) const override;
241
242 void getStandardVolumes_ref(double* vol) const override;
243
244 void getChemPotentials(double* mu) const override;
245
246 void getStandardChemPotentials(double* muStar) const override;
247
248 void getPartialMolarEnthalpies(double* hbar) const override;
249
250 void getPartialMolarEntropies(double* sbar) const override;
251
252 void getPartialMolarIntEnergies(double* ubar) const override;
253
254 void getParameters(AnyMap& phaseNode) const override;
255
256 void setParameters(const AnyMap& phaseNode,
257 const AnyMap& rootNode=AnyMap()) override;
258
259 //! Electron species name
260 string electronSpeciesName() const {
262 }
263
264 //! Return the distribution Number #m_distNum
265 int distributionNumber() const {
266 return m_distNum;
267 }
268
269 //! Return the electron energy level Number #m_levelNum
270 int levelNumber() const {
271 return m_levelNum;
272 }
273
274protected:
275 void updateThermo() const override;
276
277 //! When electron energy distribution changed, plasma properties such as
278 //! electron-collision reaction rates need to be re-evaluated.
280
281 //! When electron energy level changed, plasma properties such as
282 //! electron-collision reaction rates need to be re-evaluate.
283 //! In addition, the cross-sections need to be interpolated at
284 //! the new level.
286
287 //! Check the electron energy levels
288 /*!
289 * The values of electron energy levels need to be positive and
290 * monotonically increasing.
291 */
292 void checkElectronEnergyLevels() const;
293
294 //! Check the electron energy distribution
295 /*!
296 * This method check the electron energy distribution for the criteria
297 * below.
298 *
299 * 1. The values of electron energy distribution cannot be negative.
300 *
301 * 2. If the last value of electron energy distribution is larger
302 * than 0.01, it will raise a warning to suggest using a higher electron
303 * energy levels.
304 */
306
307 //! Update electron energy distribution.
309
310 //! Set isotropic electron energy distribution
312
313 //! Update electron temperature (K) From energy distribution.
314 //! #m_electronTemp
316
317 //! Electron energy distribution norm
319
320 // Electron energy order in the exponential term
321 double m_isotropicShapeFactor = 2.0;
322
323 //! Number of points of electron energy levels
324 size_t m_nPoints = 1001;
325
326 //! electron energy levels [ev]. Length: #m_nPoints
328
329 //! Normalized electron energy distribution vector [-]
330 //! Length: #m_nPoints
331 Eigen::ArrayXd m_electronEnergyDist;
332
333 //! Index of electron species
335
336 //! Electron temperature [K]
338
339 //! Electron energy distribution type
340 string m_distributionType = "isotropic";
341
342 //! Numerical quadrature method for electron energy distribution
343 string m_quadratureMethod = "simpson";
344
345 //! Flag of normalizing electron energy distribution
347
348private:
349 //! Electron energy distribution change variable. Whenever
350 //! #m_electronEnergyDist changes, this int is incremented.
351 int m_distNum = -1;
352
353 //! Electron energy level change variable. Whenever
354 //! #m_electronEnergyLevels changes, this int is incremented.
355 int m_levelNum = -1;
356};
357
358}
359
360#endif
ThermoPhase object for the ideal gas equation of state - workhorse for Cantera (see Thermodynamic Pro...
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:427
Class IdealGasPhase represents low-density gases that obey the ideal gas equation of state.
An error indicating that an unimplemented function has been called.
virtual double concentration(const size_t k) const
Concentration of species k.
Definition Phase.cpp:476
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:142
Base class for a phase with plasma properties.
Definition PlasmaPhase.h:58
void checkElectronEnergyDistribution() const
Check the electron energy distribution.
double meanElectronEnergy() const
Mean electron energy [eV].
int distributionNumber() const
Return the distribution Number m_distNum.
double enthalpy_mole() const override
Return the Molar enthalpy. Units: J/kmol.
void setQuadratureMethod(const string &method)
Set numerical quadrature method for integrating electron energy distribution function.
size_t m_nPoints
Number of points of electron energy levels.
int levelNumber() const
Return the electron energy level Number m_levelNum.
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 ...
string electronEnergyDistributionType() const
Get electron energy distribution type.
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...
string quadratureMethod() const
Numerical quadrature method. Method: m_quadratureMethod.
size_t nElectronEnergyLevels() const
Number of electron levels.
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.
size_t electronSpeciesIndex() const
Electron Species Index.
double m_electronTemp
Electron temperature [K].
void getElectronEnergyDistribution(double *distrb) const
Get electron energy distribution.
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 ...
double intEnergy_mole() const override
Molar internal energy. Units: J/kmol.
double entropy_mole() const override
Molar entropy.
bool m_do_normalizeElectronEnergyDist
Flag of normalizing electron energy distribution.
bool normalizeElectronEnergyDistEnabled() const
Flag of automatically normalize electron energy distribution.
void getPartialMolarIntEnergies(double *ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
double cp_mole() const override
Molar heat capacity at constant pressure.
string electronSpeciesName() const
Electron species name.
void setIsotropicElectronEnergyDistribution()
Set isotropic electron energy distribution.
double isotropicShapeFactor() const
The shape factor of isotropic electron energy distribution.
double gibbs_mole() const override
Molar Gibbs function. Units: J/kmol.
void getElectronEnergyLevels(double *levels) const
Get electron energy levels.
Definition PlasmaPhase.h:88
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.
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
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