Boltzmann equation solver for the electron energy distribution function based on the two-term approximation. More...
#include <EEDFTwoTermApproximation.h>
Boltzmann equation solver for the electron energy distribution function based on the two-term approximation.
This class implements a solver for the electron energy distribution function based on a steady-state solution to the Boltzmann equation using the classical two-term expansion, applicable to weakly ionized plasmas. The numerical approach and theory are primarily derived from the work of Hagelaar and Pitchford [14].
The two-term approximation assumes that the EEDF can be represented as:
\[ f(\epsilon, \mu) = f_0(\epsilon) + \mu f_1(\epsilon), \]
where \( \epsilon \) is the electron energy and \( \mu \) is the cosine of the angle between the electron velocity vector and the electric field. The Boltzmann equation is projected onto the zeroth moment over mu to obtain an equation for \( f_0(\epsilon) \) , the isotropic part of the distribution. The first-order anisotropic term \( f_1(\epsilon) \) is not solved directly, but is approximated and substituted into the drift and collision terms. This results in a second-order differential equation for \( f_0(\epsilon) \) alone.
The governing equation for \( f_0(\epsilon) \) is discretized on an energy grid using a finite difference method and solved using a tridiagonal matrix algorithm.
Definition at line 48 of file EEDFTwoTermApproximation.h.
Public Member Functions | |
EEDFTwoTermApproximation (PlasmaPhase *s) | |
Constructor combined with the initialization function. | |
int | calculateDistributionFunction () |
void | setLinearGrid (double &kTe_max, size_t &ncell) |
void | setGridCache () |
vector< double > | getGridEdge () const |
vector< double > | getEEDFEdge () const |
double | getElectronMobility () const |
Protected Member Functions | |
void | converge (Eigen::VectorXd &f0) |
Iterate f0 (EEDF) until convergence. | |
Eigen::VectorXd | iterate (const Eigen::VectorXd &f0, double delta) |
An iteration of solving electron energy distribution function. | |
double | integralPQ (double a, double b, double u0, double u1, double g, double x0) |
The integral in [a, b] of \(x u(x) \exp[g (x_0 - x)]\) assuming that u is linear with u(a) = u0 and u(b) = u1. | |
vector< double > | vector_g (const Eigen::VectorXd &f0) |
Vector g is used by matrix_P() and matrix_Q(). | |
Eigen::SparseMatrix< double > | matrix_P (const vector< double > &g, size_t k) |
The matrix of scattering-out. | |
Eigen::SparseMatrix< double > | matrix_Q (const vector< double > &g, size_t k) |
The matrix of scattering-in. | |
Eigen::SparseMatrix< double > | matrix_A (const Eigen::VectorXd &f0) |
Matrix A (Ax = b) of the equation of EEDF, which is discretized by the exponential scheme of Scharfetter and Gummel,. | |
double | netProductionFrequency (const Eigen::VectorXd &f0) |
Reduced net production frequency. | |
double | electronDiffusivity (const Eigen::VectorXd &f0) |
Diffusivity. | |
double | electronMobility (const Eigen::VectorXd &f0) |
Mobility. | |
void | initSpeciesIndexCrossSections () |
Initialize species indices associated with cross-section data. | |
void | updateCrossSections () |
Update the total cross sections based on the current state. | |
void | updateMoleFractions () |
Update the vector of species mole fractions. | |
void | calculateTotalElasticCrossSection () |
Compute the total elastic collision cross section. | |
void | calculateTotalCrossSection () |
Compute the total (elastic + inelastic) cross section. | |
double | norm (const Eigen::VectorXd &f, const Eigen::VectorXd &grid) |
Compute the L1 norm of a function f defined over a given energy grid. | |
Protected Attributes | |
double | m_delta0 = 1e14 |
Formerly options for the EEDF solver. | |
size_t | m_maxn = 200 |
Maximum number of iterations. | |
double | m_factorM = 4.0 |
The factor for step size change. | |
size_t | m_points = 150 |
The number of points in the EEDF grid. | |
double | m_rtol = 1e-5 |
Error tolerance for convergence. | |
std::string | m_growth = "temporal" |
The growth model of EEDF. | |
double | m_moleFractionThreshold = 0.01 |
The threshold for species mole fractions. | |
std::string | m_firstguess = "maxwell" |
The first guess for the EEDF. | |
double | m_init_kTe = 2.0 |
The initial electron temperature [eV]. | |
PlasmaPhase * | m_phase |
Pointer to the PlasmaPhase object used to initialize this object. | |
double | m_electronMobility |
Electron mobility [m²/V·s]. | |
Eigen::VectorXd | m_gridCenter |
Grid of electron energy (cell center) [eV]. | |
vector< double > | m_gridEdge |
Grid of electron energy (cell boundary i-1/2) [eV]. | |
vector< vector< size_t > > | m_j |
Location of cell j for grid cache. | |
vector< vector< size_t > > | m_i |
Location of cell i for grid cache. | |
vector< vector< vector< double > > > | m_sigma |
Cross section at the boundaries of the overlap of cell i and j. | |
vector< vector< vector< double > > > | m_eps |
The energy boundaries of the overlap of cell i and j. | |
Eigen::VectorXd | m_f0 |
normalized electron energy distribution function | |
vector< double > | m_f0_edge |
EEDF at grid edges (cell boundaries) | |
vector< double > | m_totalCrossSectionCenter |
Total electron cross section on the cell center of energy grid. | |
vector< double > | m_totalCrossSectionEdge |
Total electron cross section on the cell boundary (i-1/2) of energy grid. | |
vector< double > | m_sigmaElastic |
vector of total elastic cross section weighted with mass ratio | |
vector< size_t > | m_kTargets |
list of target species indices in global Cantera numbering (1 index per cs) | |
vector< size_t > | m_klocTargets |
list of target species indices in local X EEDF numbering (1 index per cs) | |
vector< size_t > | m_kOthers |
Indices of species which has no cross-section data. | |
vector< size_t > | m_k_lg_Targets |
Local to global indices. | |
vector< double > | m_X_targets |
Mole fraction of targets. | |
vector< double > | m_X_targets_prev |
Previous mole fraction of targets used to compute eedf. | |
vector< int > | m_inFactor |
in factor. | |
double | m_gamma |
bool | m_has_EEDF |
flag of having an EEDF | |
bool | m_first_call |
First call to calculateDistributionFunction. | |
Constructor combined with the initialization function.
This constructor initializes the EEDFTwoTermApproximation object with everything it needs to start solving EEDF.
s | PlasmaPhase object that will be used in the solver calls. |
Definition at line 22 of file EEDFTwoTermApproximation.cpp.
int calculateDistributionFunction | ( | ) |
Definition at line 46 of file EEDFTwoTermApproximation.cpp.
void setLinearGrid | ( | double & | kTe_max, |
size_t & | ncell | ||
) |
Definition at line 31 of file EEDFTwoTermApproximation.cpp.
void setGridCache | ( | ) |
Definition at line 519 of file EEDFTwoTermApproximation.cpp.
|
inline |
Definition at line 75 of file EEDFTwoTermApproximation.h.
|
inline |
Definition at line 79 of file EEDFTwoTermApproximation.h.
|
inline |
Definition at line 83 of file EEDFTwoTermApproximation.h.
|
protected |
Iterate f0 (EEDF) until convergence.
Definition at line 86 of file EEDFTwoTermApproximation.cpp.
|
protected |
An iteration of solving electron energy distribution function.
Definition at line 125 of file EEDFTwoTermApproximation.cpp.
|
protected |
The integral in [a, b] of \(x u(x) \exp[g (x_0 - x)]\) assuming that u is linear with u(a) = u0 and u(b) = u1.
Definition at line 176 of file EEDFTwoTermApproximation.cpp.
|
protected |
Vector g is used by matrix_P() and matrix_Q().
\[ g_i = \frac{1}{\epsilon_{i+1} - \epsilon_{i-1}} \ln(\frac{F_{0, i+1}}{F_{0, i-1}}) \]
Definition at line 203 of file EEDFTwoTermApproximation.cpp.
|
protected |
The matrix of scattering-out.
\[ P_{i,k} = \gamma \int_{\epsilon_i - 1/2}^{\epsilon_i + 1/2} \epsilon \sigma_k exp[(\epsilon_i - \epsilon)g_i] d \epsilon \]
Definition at line 228 of file EEDFTwoTermApproximation.cpp.
|
protected |
The matrix of scattering-in.
\[ Q_{i,j,k} = \gamma \int_{\epsilon_1}^{\epsilon_2} \epsilon \sigma_k exp[(\epsilon_j - \epsilon)g_j] d \epsilon \]
where the interval \([\epsilon_1, \epsilon_2]\) is the overlap of cell j, and cell i shifted by the threshold energy:
\[ \epsilon_1 = \min(\max(\epsilon_{i-1/2}+u_k, \epsilon_{j-1/2}),\epsilon_{j+1/2}), \]
\[ \epsilon_2 = \min(\max(\epsilon_{i+1/2}+u_k, \epsilon_{j-1/2}),\epsilon_{j+1/2}) \]
Definition at line 247 of file EEDFTwoTermApproximation.cpp.
|
protected |
Matrix A (Ax = b) of the equation of EEDF, which is discretized by the exponential scheme of Scharfetter and Gummel,.
\[ \left[ \tilde{W} F_0 - \tilde{D} \frac{d F_0}{\epsilon} \right]_{i+1/2} = \frac{\tilde{W}_{i+1/2} F_{0,i}}{1 - \exp[-z_{i+1/2}]} + \frac{\tilde{W}_{i+1/2} F_{0,i+1}}{1 - \exp[z_{i+1/2}]} \]
where \( z_{i+1/2} = \tilde{w}_{i+1/2} / \tilde{D}_{i+1/2} \) (Peclet number).
Definition at line 267 of file EEDFTwoTermApproximation.cpp.
|
protected |
Reduced net production frequency.
Equation (10) of ref. [1] divided by N.
f0 | EEDF |
Definition at line 362 of file EEDFTwoTermApproximation.cpp.
|
protected |
Diffusivity.
Definition at line 381 of file EEDFTwoTermApproximation.cpp.
|
protected |
Mobility.
Definition at line 397 of file EEDFTwoTermApproximation.cpp.
|
protected |
Initialize species indices associated with cross-section data.
Definition at line 415 of file EEDFTwoTermApproximation.cpp.
|
protected |
Update the total cross sections based on the current state.
Definition at line 461 of file EEDFTwoTermApproximation.cpp.
|
protected |
Update the vector of species mole fractions.
Definition at line 469 of file EEDFTwoTermApproximation.cpp.
|
protected |
Compute the total elastic collision cross section.
Definition at line 502 of file EEDFTwoTermApproximation.cpp.
|
protected |
Compute the total (elastic + inelastic) cross section.
Definition at line 483 of file EEDFTwoTermApproximation.cpp.
|
protected |
Compute the L1 norm of a function f defined over a given energy grid.
f | Vector representing the function values (EEDF) |
grid | Vector representing the energy grid corresponding to f |
Definition at line 590 of file EEDFTwoTermApproximation.cpp.
|
protected |
Formerly options for the EEDF solver.
The first step size
Definition at line 92 of file EEDFTwoTermApproximation.h.
|
protected |
Maximum number of iterations.
Definition at line 95 of file EEDFTwoTermApproximation.h.
|
protected |
The factor for step size change.
Definition at line 98 of file EEDFTwoTermApproximation.h.
|
protected |
The number of points in the EEDF grid.
Definition at line 101 of file EEDFTwoTermApproximation.h.
|
protected |
Error tolerance for convergence.
Definition at line 104 of file EEDFTwoTermApproximation.h.
|
protected |
The growth model of EEDF.
Definition at line 107 of file EEDFTwoTermApproximation.h.
|
protected |
The threshold for species mole fractions.
Definition at line 110 of file EEDFTwoTermApproximation.h.
|
protected |
The first guess for the EEDF.
Definition at line 113 of file EEDFTwoTermApproximation.h.
|
protected |
The initial electron temperature [eV].
Definition at line 116 of file EEDFTwoTermApproximation.h.
|
protected |
Pointer to the PlasmaPhase object used to initialize this object.
This PlasmaPhase object provides species, element, and cross-section data used by the EEDF solver. It is set during construction and is not modified afterwards. All subsequent calls to compute functions must use the same PlasmaPhase context.
Definition at line 125 of file EEDFTwoTermApproximation.h.
|
protected |
Electron mobility [m²/V·s].
Definition at line 219 of file EEDFTwoTermApproximation.h.
|
protected |
Grid of electron energy (cell center) [eV].
Definition at line 222 of file EEDFTwoTermApproximation.h.
|
protected |
Grid of electron energy (cell boundary i-1/2) [eV].
Definition at line 225 of file EEDFTwoTermApproximation.h.
|
protected |
Location of cell j for grid cache.
Definition at line 228 of file EEDFTwoTermApproximation.h.
|
protected |
Location of cell i for grid cache.
Definition at line 231 of file EEDFTwoTermApproximation.h.
|
protected |
Cross section at the boundaries of the overlap of cell i and j.
Definition at line 234 of file EEDFTwoTermApproximation.h.
|
protected |
The energy boundaries of the overlap of cell i and j.
Definition at line 237 of file EEDFTwoTermApproximation.h.
|
protected |
normalized electron energy distribution function
Definition at line 240 of file EEDFTwoTermApproximation.h.
|
protected |
EEDF at grid edges (cell boundaries)
Definition at line 243 of file EEDFTwoTermApproximation.h.
|
protected |
Total electron cross section on the cell center of energy grid.
Definition at line 246 of file EEDFTwoTermApproximation.h.
|
protected |
Total electron cross section on the cell boundary (i-1/2) of energy grid.
Definition at line 250 of file EEDFTwoTermApproximation.h.
|
protected |
vector of total elastic cross section weighted with mass ratio
Definition at line 253 of file EEDFTwoTermApproximation.h.
|
protected |
list of target species indices in global Cantera numbering (1 index per cs)
Definition at line 256 of file EEDFTwoTermApproximation.h.
|
protected |
list of target species indices in local X EEDF numbering (1 index per cs)
Definition at line 259 of file EEDFTwoTermApproximation.h.
|
protected |
Indices of species which has no cross-section data.
Definition at line 262 of file EEDFTwoTermApproximation.h.
|
protected |
Local to global indices.
Definition at line 265 of file EEDFTwoTermApproximation.h.
|
protected |
Mole fraction of targets.
Definition at line 268 of file EEDFTwoTermApproximation.h.
|
protected |
Previous mole fraction of targets used to compute eedf.
Definition at line 271 of file EEDFTwoTermApproximation.h.
|
protected |
in factor.
This is used for calculating the Q matrix of scattering-in processes.
Definition at line 275 of file EEDFTwoTermApproximation.h.
|
protected |
Definition at line 277 of file EEDFTwoTermApproximation.h.
|
protected |
flag of having an EEDF
Definition at line 280 of file EEDFTwoTermApproximation.h.
|
protected |
First call to calculateDistributionFunction.
Definition at line 283 of file EEDFTwoTermApproximation.h.