Cantera 2.6.0
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 * References:
48 *
49 * [1] J. T. Gudmundsson. On the effect of the electron energy distribution on the
50 * plasma parameters of an argon discharge: a global (volume-averaged) model study.
51 * Plasma Sources Science and Technology, 10.1 (2001): 76.
52 * doi: https://doi.org/10.1088/0963-0252/10/1/310
53 *
54 * [2] H. Khalilpour and G. Foroutan. The effects of electron energy distribution
55 * function on the plasma sheath structure in the presence of charged nanoparticles
56 * Journal of Plasma Physics 86.2 (2020).
57 * doi: https://doi.org/10.1017/S0022377820000161
58 *
59 * [3] G. J. M. Hagelaar and L. C. Pitchford
60 * "Solving the Boltzmann equation to obtain electron transport
61 * coefficients and rate coefficients for fluid models."
62 * Plasma Sources Science and Technology 14.4 (2005): 722.
63 * doi: https://doi.org/10.1088/0963-0252/14/4/011
64 *
65 * [4] A. Luque, "BOLOS: An open source solver for the Boltzmann equation,"
66 * https://github.com/aluque/bolos.
67 *
68 * @warning This class is an experimental part of %Cantera and may be
69 * changed or removed without notice.
70 * @todo Implement electron Boltzmann equation solver to solve EEDF.
71 * https://github.com/Cantera/enhancements/issues/127
72 * @ingroup phase
73 */
75{
76public:
77 //! Construct and initialize a PlasmaPhase object
78 //! directly from an input file. The constructor initializes the electron
79 //! energy distribution to be Druyvesteyn distribution (m_x = 2.0). The initial
80 //! electron energy grid is set to a linear space which starts at 0.01 eV and ends
81 //! at 1 eV with 1000 points.
82 /*!
83 * @param inputFile Name of the input file containing the phase definition
84 * to set up the object. If blank, an empty phase will be
85 * created.
86 * @param id ID of the phase in the input file. Defaults to the
87 * empty string.
88 */
89 explicit PlasmaPhase(const std::string& inputFile="",
90 const std::string& id="");
91
92 virtual std::string type() const {
93 return "plasma";
94 }
95
96 virtual void initThermo();
97
98 //! Set electron energy levels.
99 //! @param levels The vector of electron energy levels (eV).
100 //! Length: #m_nPoints.
101 //! @param length The length of the \p levels.
102 void setElectronEnergyLevels(const double* levels, size_t length);
103
104 //! Get electron energy levels.
105 //! @param levels The vector of electron energy levels (eV). Length: #m_nPoints
106 void getElectronEnergyLevels(double* levels) const {
107 Eigen::Map<Eigen::ArrayXd>(levels, m_nPoints) = m_electronEnergyLevels;
108 }
109
110 //! Set discretized electron energy distribution.
111 //! @param levels The vector of electron energy levels (eV).
112 //! Length: #m_nPoints.
113 //! @param distrb The vector of electron energy distribution.
114 //! Length: #m_nPoints.
115 //! @param length The length of the vectors, which equals #m_nPoints.
116 void setDiscretizedElectronEnergyDist(const double* levels,
117 const double* distrb,
118 size_t length);
119
120 //! Get electron energy distribution.
121 //! @param distrb The vector of electron energy distribution.
122 //! Length: #m_nPoints.
123 void getElectronEnergyDistribution(double* distrb) const {
124 Eigen::Map<Eigen::ArrayXd>(distrb, m_nPoints) = m_electronEnergyDist;
125 }
126
127 //! Set the shape factor of isotropic electron energy distribution.
128 //! Note that \f$ x = 1 \f$ and \f$ x = 2 \f$ correspond to the
129 //! Maxwellian and Druyvesteyn distribution, respectively.
130 //! @param x The shape factor
131 void setIsotropicShapeFactor(double x);
132
133 //! The shape factor of isotropic electron energy distribution
134 double isotropicShapeFactor() const {
135 return m_isotropicShapeFactor;
136 }
137
138 //! Set the internally stored electron temperature of the phase (K).
139 //! @param Te Electron temperature in Kelvin
140 virtual void setElectronTemperature(double Te);
141
142 //! Set mean electron energy [eV]. This method also sets electron temperature
143 //! accordingly.
144 void setMeanElectronEnergy(double energy);
145
146 //! Get electron energy distribution type
147 std::string electronEnergyDistributionType() const {
148 return m_distributionType;
149 }
150
151 //! Set electron energy distribution type
152 void setElectronEnergyDistributionType(const std::string& type);
153
154 //! Numerical quadrature method. Method: #m_quadratureMethod
155 std::string quadratureMethod() const {
156 return m_quadratureMethod;
157 }
158
159 //! Set numerical quadrature method for intergating electron
160 //! energy distribution function. Method: #m_quadratureMethod
161 void setQuadratureMethod(const std::string& method) {
162 m_quadratureMethod = method;
163 }
164
165 //! Mean electron energy [eV]
166 double meanElectronEnergy() const {
167 return 3.0 / 2.0 * electronTemperature() * Boltzmann / ElectronCharge;
168 }
169
170 //! Set flag of automatically normalize electron energy distribution
171 //! Flag: #m_do_normalizeElectronEnergyDist
174 }
175
176 //! Flag of automatically normalize electron energy distribution.
177 //! Flag: #m_do_normalizeElectronEnergyDist
180 }
181
182 virtual bool addSpecies(shared_ptr<Species> spec);
183
184 //! Electron Temperature (K)
185 //! @return The electron temperature of the phase
186 virtual double electronTemperature() const {
187 return m_electronTemp;
188 }
189
190 //! Number of electron levels
191 size_t nElectronEnergyLevels() const {
192 return m_nPoints;
193 }
194
195 virtual void getParameters(AnyMap& phaseNode) const;
196
197 virtual void setParameters(const AnyMap& phaseNode,
198 const AnyMap& rootNode=AnyMap());
199
200protected:
201 virtual void updateThermo() const;
202
203 //! Check the electron energy levels
204 /*!
205 * The values of electron energy levels need to be positive and
206 * monotonically increasing.
207 */
208 void checkElectronEnergyLevels() const;
209
210 //! Check the electron energy distribution
211 /*!
212 * This method check the electron energy distribution for the criteria
213 * below.
214 *
215 * 1. The values of electron energy distribution cannot be negative.
216 *
217 * 2. If the last value of electron energy distribution is larger
218 * than 0.01, it will raise a warning to suggest using a higher electron
219 * energy levels.
220 */
222
223 //! Update electron energy distribution.
225
226 //! Set isotropic electron energy distribution
228
229 //! Update electron temperature (K) From energy distribution.
230 //! #m_electronTemp
232
233 //! Electron energy distribution norm
235
236 // Electron energy order in the exponential term
237 double m_isotropicShapeFactor;
238
239 //! Number of points of electron energy levels
240 size_t m_nPoints;
241
242 //! electron energy levels [ev]. Length: #m_nPoints
244
245 //! Normalized electron energy distribution vector [-]
246 //! Length: #m_nPoints
247 Eigen::ArrayXd m_electronEnergyDist;
248
249 //! Index of electron species
251
252 //! Electron temperature [K]
254
255 //! Electron energy distribution type
257
258 //! Numerical quadrature method for electron energy distribution
260
261 //! Flag of normalizing electron energy distribution
263};
264
265}
266
267#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:399
Class IdealGasPhase represents low-density gases that obey the ideal gas equation of state.
Base class for a phase with plasma properties.
Definition: PlasmaPhase.h:75
void checkElectronEnergyDistribution() const
Check the electron energy distribution.
double meanElectronEnergy() const
Mean electron energy [eV].
Definition: PlasmaPhase.h:166
size_t m_nPoints
Number of points of electron energy levels.
Definition: PlasmaPhase.h:240
virtual bool addSpecies(shared_ptr< Species > spec)
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 ...
void normalizeElectronEnergyDistribution()
Electron energy distribution norm.
Definition: PlasmaPhase.cpp:41
void setQuadratureMethod(const std::string &method)
Set numerical quadrature method for intergating electron energy distribution function.
Definition: PlasmaPhase.h:161
std::string m_quadratureMethod
Numerical quadrature method for electron energy distribution.
Definition: PlasmaPhase.h:259
virtual void setElectronTemperature(double Te)
Set the internally stored electron temperature of the phase (K).
Definition: PlasmaPhase.cpp:80
PlasmaPhase(const std::string &inputFile="", const std::string &id="")
Construct and initialize a PlasmaPhase object directly from an input file.
Definition: PlasmaPhase.cpp:14
size_t nElectronEnergyLevels() const
Number of electron levels.
Definition: PlasmaPhase.h:191
Eigen::ArrayXd m_electronEnergyDist
Normalized electron energy distribution vector [-] Length: m_nPoints.
Definition: PlasmaPhase.h:247
Eigen::ArrayXd m_electronEnergyLevels
electron energy levels [ev]. Length: m_nPoints
Definition: PlasmaPhase.h:243
void setDiscretizedElectronEnergyDist(const double *levels, const double *distrb, size_t length)
Set discretized electron energy distribution.
void checkElectronEnergyLevels() const
Check the electron energy levels.
Definition: PlasmaPhase.cpp:98
void updateElectronTemperatureFromEnergyDist()
Update electron temperature (K) From energy distribution.
void setElectronEnergyDistributionType(const std::string &type)
Set electron energy distribution type.
Definition: PlasmaPhase.cpp:53
void updateElectronEnergyDistribution()
Update electron energy distribution.
Definition: PlasmaPhase.cpp:31
std::string quadratureMethod() const
Numerical quadrature method. Method: m_quadratureMethod.
Definition: PlasmaPhase.h:155
virtual double electronTemperature() const
Electron Temperature (K)
Definition: PlasmaPhase.h:186
virtual void updateThermo() const
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
double m_electronTemp
Electron temperature [K].
Definition: PlasmaPhase.h:253
void getElectronEnergyDistribution(double *distrb) const
Get electron energy distribution.
Definition: PlasmaPhase.h:123
void setElectronEnergyLevels(const double *levels, size_t length)
Set electron energy levels.
Definition: PlasmaPhase.cpp:90
bool m_do_normalizeElectronEnergyDist
Flag of normalizing electron energy distribution.
Definition: PlasmaPhase.h:262
virtual std::string type() const
String indicating the thermodynamic model implemented.
Definition: PlasmaPhase.h:92
bool normalizeElectronEnergyDistEnabled() const
Flag of automatically normalize electron energy distribution.
Definition: PlasmaPhase.h:178
std::string m_distributionType
Electron energy distribution type.
Definition: PlasmaPhase.h:256
void setIsotropicElectronEnergyDistribution()
Set isotropic electron energy distribution.
Definition: PlasmaPhase.cpp:64
double isotropicShapeFactor() const
The shape factor of isotropic electron energy distribution.
Definition: PlasmaPhase.h:134
void getElectronEnergyLevels(double *levels) const
Get electron energy levels.
Definition: PlasmaPhase.h:106
std::string electronEnergyDistributionType() const
Get electron energy distribution type.
Definition: PlasmaPhase.h:147
void setMeanElectronEnergy(double energy)
Set mean electron energy [eV].
Definition: PlasmaPhase.cpp:85
size_t m_electronSpeciesIndex
Index of electron species.
Definition: PlasmaPhase.h:250
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...
Definition: PlasmaPhase.h:172
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
const double Boltzmann
Boltzmann constant [J/K].
Definition: ct_defs.h:69
const double ElectronCharge
Elementary charge [C].
Definition: ct_defs.h:75