Cantera  4.0.0a1
Loading...
Searching...
No Matches
EEDFTwoTermApproximation.cpp
Go to the documentation of this file.
1/**
2 * @file EEDFTwoTermApproximation.cpp
3 * EEDF Two-Term approximation solver. Implementation file for class
4 * EEDFTwoTermApproximation.
5 */
6
7// This file is part of Cantera. See License.txt in the top-level directory or
8// at https://cantera.org/license.txt for license and copyright information.
9
12#include "cantera/numerics/eigen_dense.h"
16
17namespace Cantera
18{
19
20typedef Eigen::SparseMatrix<double> SparseMat;
21
22EEDFTwoTermApproximation::EEDFTwoTermApproximation(PlasmaPhase* s)
23{
24 // store a pointer to s.
25 m_phase = s;
26 m_first_call = true;
27 m_has_EEDF = false;
28 m_gamma = pow(2.0 * ElectronCharge / ElectronMass, 0.5);
29}
30
31void EEDFTwoTermApproximation::setLinearGrid(double& kTe_max, size_t& ncell)
32{
33 m_points = ncell;
34 m_gridCenter.resize(m_points);
35 m_gridEdge.resize(m_points + 1);
36 m_f0.resize(m_points);
37 m_f0_edge.resize(m_points + 1);
38 for (size_t j = 0; j < m_points; j++) {
39 m_gridCenter[j] = kTe_max * ( j + 0.5 ) / m_points;
40 m_gridEdge[j] = kTe_max * j / m_points;
41 }
42 m_gridEdge[m_points] = kTe_max;
43 setGridCache();
44}
45
46int EEDFTwoTermApproximation::calculateDistributionFunction()
47{
48 if (m_first_call) {
50 m_first_call = false;
51 }
52
55
56 if (!m_has_EEDF) {
57 if (m_firstguess == "maxwell") {
58 for (size_t j = 0; j < m_points; j++) {
59 m_f0(j) = 2.0 * pow(1.0 / Pi, 0.5) * pow(m_init_kTe, -3. / 2.) *
60 exp(-m_gridCenter[j] / m_init_kTe);
61 }
62 } else {
63 throw CanteraError("EEDFTwoTermApproximation::calculateDistributionFunction",
64 " unknown EEDF first guess");
65 }
66 }
67
69
70 // write the EEDF at grid edges
71 vector<double> f(m_f0.data(), m_f0.data() + m_f0.rows() * m_f0.cols());
72 vector<double> x(m_gridCenter.data(), m_gridCenter.data() + m_gridCenter.rows() * m_gridCenter.cols());
73 for (size_t i = 0; i < m_points + 1; i++) {
74 m_f0_edge[i] = linearInterp(m_gridEdge[i], x, f);
75 }
76
77 m_has_EEDF = true;
78
79 // update electron mobility
81 return 0;
82}
83
84void EEDFTwoTermApproximation::converge(Eigen::VectorXd& f0)
85{
86 double err0 = 0.0;
87 double err1 = 0.0;
88 double delta = m_delta0;
89
90 if (m_maxn == 0) {
91 throw CanteraError("EEDFTwoTermApproximation::converge",
92 "m_maxn is zero; no iterations will occur.");
93 }
94 if (m_points == 0) {
95 throw CanteraError("EEDFTwoTermApproximation::converge",
96 "m_points is zero; the EEDF grid is empty.");
97 }
98 if (isnan(delta) || delta == 0.0) {
99 throw CanteraError("EEDFTwoTermApproximation::converge",
100 "m_delta0 is NaN or zero; solver cannot update.");
101 }
102
103 for (size_t n = 0; n < m_maxn; n++) {
104 if (0.0 < err1 && err1 < err0) {
105 delta *= log(m_factorM) / (log(err0) - log(err1));
106 }
107
108 Eigen::VectorXd f0_old = f0;
109 f0 = iterate(f0_old, delta);
110 checkFinite("EEDFTwoTermApproximation::converge: f0", asSpan(f0));
111
112 err0 = err1;
113 Eigen::VectorXd Df0 = (f0_old - f0).cwiseAbs();
114 err1 = norm(Df0, m_gridCenter);
115 if (err1 < m_rtol) {
116 break;
117 } else if (n == m_maxn - 1) {
118 throw CanteraError("WeaklyIonizedGas::converge", "Convergence failed");
119 }
120 }
121}
122
123Eigen::VectorXd EEDFTwoTermApproximation::iterate(const Eigen::VectorXd& f0, double delta)
124{
125 // CQM multiple call to vector_* and matrix_*
126 // probably extremely ineficient
127 // must be refactored!!
128
129 SparseMat PQ(m_points, m_points);
130 vector<double> g = vector_g(f0);
131
132 for (size_t k : m_phase->kInelastic()) {
133 SparseMat Q_k = matrix_Q(g, k);
134 SparseMat P_k = matrix_P(g, k);
135 PQ += (matrix_Q(g, k) - matrix_P(g, k)) * m_X_targets[m_klocTargets[k]];
136 }
137
138 SparseMat A = matrix_A(f0);
139 SparseMat I(m_points, m_points);
140 for (size_t i = 0; i < m_points; i++) {
141 I.insert(i,i) = 1.0;
142 }
143 A -= PQ;
144 A *= delta;
145 A += I;
146
147 // SparseLU :
148 Eigen::SparseLU<SparseMat> solver(A);
149 if (solver.info() == Eigen::NumericalIssue) {
150 throw CanteraError("EEDFTwoTermApproximation::iterate",
151 "Error SparseLU solver: NumericalIssue");
152 } else if (solver.info() == Eigen::InvalidInput) {
153 throw CanteraError("EEDFTwoTermApproximation::iterate",
154 "Error SparseLU solver: InvalidInput");
155 }
156 if (solver.info() != Eigen::Success) {
157 throw CanteraError("EEDFTwoTermApproximation::iterate",
158 "Error SparseLU solver", "Decomposition failed");
159 return f0;
160 }
161
162 // solve f0
163 Eigen::VectorXd f1 = solver.solve(f0);
164 if(solver.info() != Eigen::Success) {
165 throw CanteraError("EEDFTwoTermApproximation::iterate", "Solving failed");
166 return f0;
167 }
168
169 checkFinite("EEDFTwoTermApproximation::converge: f0", asSpan(f1));
170 f1 /= norm(f1, m_gridCenter);
171 return f1;
172}
173
174double EEDFTwoTermApproximation::integralPQ(double a, double b, double u0, double u1,
175 double g, double x0)
176{
177 double A1;
178 double A2;
179 if (g != 0.0) {
180 double expm1a = expm1(g * (-a + x0));
181 double expm1b = expm1(g * (-b + x0));
182 double ag = a * g;
183 double ag1 = ag + 1;
184 double bg = b * g;
185 double bg1 = bg + 1;
186 A1 = (expm1a * ag1 + ag - expm1b * bg1 - bg) / (g*g);
187 A2 = (expm1a * (2 * ag1 + ag * ag) + ag * (ag + 2) -
188 expm1b * (2 * bg1 + bg * bg) - bg * (bg + 2)) / (g*g*g);
189 } else {
190 A1 = 0.5 * (b*b - a*a);
191 A2 = 1.0 / 3.0 * (b*b*b - a*a*a);
192 }
193
194 // The interpolation formula of u(x) = c0 + c1 * x
195 double c0 = (a * u1 - b * u0) / (a - b);
196 double c1 = (u0 - u1) / (a - b);
197
198 return c0 * A1 + c1 * A2;
199}
200
201vector<double> EEDFTwoTermApproximation::vector_g(const Eigen::VectorXd& f0)
202{
203 vector<double> g(m_points, 0.0);
204 const double f_min = 1e-300; // Smallest safe floating-point value
205
206 // Handle first point (i = 0)
207 double f1 = std::max(f0(1), f_min);
208 double f0_ = std::max(f0(0), f_min);
209 g[0] = log(f1 / f0_) / (m_gridCenter[1] - m_gridCenter[0]);
210
211 // Handle last point (i = N)
212 size_t N = m_points - 1;
213 double fN = std::max(f0(N), f_min);
214 double fNm1 = std::max(f0(N - 1), f_min);
215 g[N] = log(fN / fNm1) / (m_gridCenter[N] - m_gridCenter[N - 1]);
216
217 // Handle interior points
218 for (size_t i = 1; i < N; ++i) {
219 double f_up = std::max(f0(i + 1), f_min);
220 double f_down = std::max(f0(i - 1), f_min);
221 g[i] = log(f_up / f_down) / (m_gridCenter[i + 1] - m_gridCenter[i - 1]);
222 }
223 return g;
224}
225
226SparseMat EEDFTwoTermApproximation::matrix_P(span<const double> g, size_t k)
227{
228 SparseTriplets tripletList;
229 for (size_t n = 0; n < m_eps[k].size(); n++) {
230 double eps_a = m_eps[k][n][0];
231 double eps_b = m_eps[k][n][1];
232 double sigma_a = m_sigma[k][n][0];
233 double sigma_b = m_sigma[k][n][1];
234 auto j = static_cast<SparseMat::StorageIndex>(m_j[k][n]);
235 double r = integralPQ(eps_a, eps_b, sigma_a, sigma_b, g[j], m_gridCenter[j]);
236 double p = m_gamma * r;
237
238 tripletList.emplace_back(j, j, p);
239 }
240 SparseMat P(m_points, m_points);
241 P.setFromTriplets(tripletList.begin(), tripletList.end());
242 return P;
243}
244
245SparseMat EEDFTwoTermApproximation::matrix_Q(span<const double> g, size_t k)
246{
247 SparseTriplets tripletList;
248 for (size_t n = 0; n < m_eps[k].size(); n++) {
249 double eps_a = m_eps[k][n][0];
250 double eps_b = m_eps[k][n][1];
251 double sigma_a = m_sigma[k][n][0];
252 double sigma_b = m_sigma[k][n][1];
253 auto i = static_cast<SparseMat::StorageIndex>(m_i[k][n]);
254 auto j = static_cast<SparseMat::StorageIndex>(m_j[k][n]);
255 double r = integralPQ(eps_a, eps_b, sigma_a, sigma_b, g[j], m_gridCenter[j]);
256 double q = m_inFactor[k] * m_gamma * r;
257
258 tripletList.emplace_back(i, j, q);
259 }
260 SparseMat Q(m_points, m_points);
261 Q.setFromTriplets(tripletList.begin(), tripletList.end());
262 return Q;
263}
264
265SparseMat EEDFTwoTermApproximation::matrix_A(const Eigen::VectorXd& f0)
266{
267 vector<double> a0(m_points + 1);
268 vector<double> a1(m_points + 1);
269 size_t N = m_points - 1;
270 // Scharfetter-Gummel scheme
271 double nu = netProductionFrequency(f0);
272 a0[0] = NAN;
273 a1[0] = NAN;
274 a0[N+1] = NAN;
275 a1[N+1] = NAN;
276
277 double nDensity = m_phase->molarDensity() * Avogadro;
278 double alpha;
279 double E = m_phase->electricField();
280 if (m_growth == "spatial") {
281 double mu = electronMobility(f0);
282 double D = electronDiffusivity(f0);
283 alpha = (mu * E - sqrt(pow(mu * E, 2) - 4 * D * nu * nDensity)) / 2.0 / D / nDensity;
284 } else {
285 alpha = 0.0;
286 }
287
288 double sigma_tilde;
289 double omega = 2 * Pi * m_phase->electricFieldFrequency();
290 for (size_t j = 1; j < m_points; j++) {
291 if (m_growth == "temporal") {
292 sigma_tilde = m_totalCrossSectionEdge[j] + nu / pow(m_gridEdge[j], 0.5) / m_gamma;
293 } else {
294 sigma_tilde = m_totalCrossSectionEdge[j];
295 }
296 double q = omega / (nDensity * m_gamma * pow(m_gridEdge[j], 0.5));
297 double W = -m_gamma * m_gridEdge[j] * m_gridEdge[j] * m_sigmaElastic[j];
298 double F = sigma_tilde * sigma_tilde / (sigma_tilde * sigma_tilde + q * q);
299 double DA = m_gamma / 3.0 * pow(E / nDensity, 2.0) * m_gridEdge[j];
300 double DB = m_gamma * m_phase->temperature() * Boltzmann / ElectronCharge * m_gridEdge[j] * m_gridEdge[j] * m_sigmaElastic[j];
301 double D = DA / sigma_tilde * F + DB;
302 if (m_growth == "spatial") {
303 W -= m_gamma / 3.0 * 2 * alpha * E / nDensity * m_gridEdge[j] / sigma_tilde;
304 }
305
306 double z = W * (m_gridCenter[j] - m_gridCenter[j-1]) / D;
307 if (!std::isfinite(z)) {
308 throw CanteraError("matrix_A", "Non-finite Peclet number encountered");
309 }
310 if (std::abs(z) > 500) {
311 warn_user("EEDFTwoTermApproximation::matrix_A",
312 "Large Peclet number z = {:.3e} at j = {}. "
313 "W = {:.3e}, D = {:.3e}, E/N = {:.3e}\n",
314 z, j, W, D, E / nDensity);
315 }
316 a0[j] = W / (1 - std::exp(-z));
317 a1[j] = W / (1 - std::exp(z));
318 }
319
320 SparseTriplets tripletList;
321 // center diagonal
322 // zero flux b.c. at energy = 0
323 tripletList.emplace_back(0, 0, a0[1]);
324
325 for (size_t j = 1; j < m_points - 1; j++) {
326 tripletList.emplace_back(j, j, a0[j+1] - a1[j]);
327 }
328
329 // upper diagonal
330 for (size_t j = 0; j < m_points - 1; j++) {
331 tripletList.emplace_back(j, j+1, a1[j+1]);
332 }
333
334 // lower diagonal
335 for (size_t j = 1; j < m_points; j++) {
336 tripletList.emplace_back(j, j-1, -a0[j]);
337 }
338
339 // zero flux b.c.
340 tripletList.emplace_back(N, N, -a1[N]);
341
342 SparseMat A(m_points, m_points);
343 A.setFromTriplets(tripletList.begin(), tripletList.end());
344
345 //plus G
346 SparseMat G(m_points, m_points);
347 if (m_growth == "temporal") {
348 for (size_t i = 0; i < m_points; i++) {
349 G.insert(i, i) = 2.0 / 3.0 * (pow(m_gridEdge[i+1], 1.5) - pow(m_gridEdge[i], 1.5)) * nu;
350 }
351 } else if (m_growth == "spatial") {
352 double nDensity = m_phase->molarDensity() * Avogadro;
353 for (size_t i = 0; i < m_points; i++) {
354 double sigma_c = 0.5 * (m_totalCrossSectionEdge[i] + m_totalCrossSectionEdge[i + 1]);
355 G.insert(i, i) = - alpha * m_gamma / 3 * (alpha * (pow(m_gridEdge[i + 1], 2) - pow(m_gridEdge[i], 2)) / sigma_c / 2
356 - E / nDensity * (m_gridEdge[i + 1] / m_totalCrossSectionEdge[i + 1] - m_gridEdge[i] / m_totalCrossSectionEdge[i]));
357 }
358 }
359 return A + G;
360}
361
363{
364 double nu = 0.0;
365 vector<double> g = vector_g(f0);
366
367 for (size_t k = 0; k < m_phase->nCollisions(); k++) {
368 if (m_phase->collisionRate(k)->kind() == "ionization" ||
369 m_phase->collisionRate(k)->kind() == "attachment") {
370 SparseMat PQ = (matrix_Q(g, k) - matrix_P(g, k)) *
372 Eigen::VectorXd s = PQ * f0;
373 checkFinite("EEDFTwoTermApproximation::netProductionFrequency: s",
374 asSpan(s));
375 nu += s.sum();
376 }
377 }
378 return nu;
379}
380
381double EEDFTwoTermApproximation::electronDiffusivity(const Eigen::VectorXd& f0)
382{
383 vector<double> y(m_points, 0.0);
384 double nu = netProductionFrequency(f0);
385 for (size_t i = 0; i < m_points; i++) {
386 if (m_gridCenter[i] != 0.0) {
387 y[i] = m_gridCenter[i] * f0(i) /
388 (m_totalCrossSectionCenter[i] + nu / m_gamma / pow(m_gridCenter[i], 0.5));
389 }
390 }
391 double nDensity = m_phase->molarDensity() * Avogadro;
392 auto f = Eigen::Map<const Eigen::ArrayXd>(y.data(), y.size());
393 auto x = Eigen::Map<const Eigen::ArrayXd>(m_gridCenter.data(), m_gridCenter.size());
394 return 1./3. * m_gamma * simpson(f, x) / nDensity;
395}
396
397double EEDFTwoTermApproximation::electronMobility(const Eigen::VectorXd& f0)
398{
399 double nu = netProductionFrequency(f0);
400 vector<double> y(m_points + 1, 0.0);
401 for (size_t i = 1; i < m_points; i++) {
402 // calculate df0 at i-1/2
403 double df0 = (f0(i) - f0(i-1)) / (m_gridCenter[i] - m_gridCenter[i-1]);
404 if (m_gridEdge[i] != 0.0) {
405 y[i] = m_gridEdge[i] * df0 /
406 (m_totalCrossSectionEdge[i] + nu / m_gamma / pow(m_gridEdge[i], 0.5));
407 }
408 }
409 double nDensity = m_phase->molarDensity() * Avogadro;
410 return -1./3. * m_gamma * simpson(asVectorXd(y), asVectorXd(m_gridEdge)) / nDensity;
411}
412
414{
415 // set up target index
416 m_kTargets.resize(m_phase->nCollisions());
418 m_inFactor.resize(m_phase->nCollisions());
419 for (size_t k = 0; k < m_phase->nCollisions(); k++) {
421 // Check if it is a new target or not :
422 auto it = find(m_k_lg_Targets.begin(), m_k_lg_Targets.end(), m_kTargets[k]);
423
424 if (it == m_k_lg_Targets.end()){
425 m_k_lg_Targets.push_back(m_kTargets[k]);
426 m_klocTargets[k] = m_k_lg_Targets.size() - 1;
427 } else {
428 m_klocTargets[k] = distance(m_k_lg_Targets.begin(), it);
429 }
430
431 const auto& kind = m_phase->collisionRate(k)->kind();
432
433 if (kind == "ionization") {
434 m_inFactor[k] = 2;
435 } else if (kind == "attachment") {
436 m_inFactor[k] = 0;
437 } else {
438 m_inFactor[k] = 1;
439 }
440 }
441
442 m_X_targets.resize(m_k_lg_Targets.size());
443 m_X_targets_prev.resize(m_k_lg_Targets.size());
444 for (size_t k = 0; k < m_X_targets.size(); k++) {
445 size_t k_glob = m_k_lg_Targets[k];
446 m_X_targets[k] = m_phase->moleFraction(k_glob);
448 }
449
450 // set up indices of species which has no cross-section data
451 for (size_t k = 0; k < m_phase->nSpecies(); k++) {
452 auto it = std::find(m_kTargets.begin(), m_kTargets.end(), k);
453 if (it == m_kTargets.end()) {
454 m_kOthers.push_back(k);
455 }
456 }
457}
458
460{
461 // Compute sigma_m and sigma_\epsilon
464}
465
466// Update the species mole fractions used for EEDF computation
468{
469 double tmp_sum = 0.0;
470 for (size_t k = 0; k < m_X_targets.size(); k++) {
472 tmp_sum = tmp_sum + m_phase->moleFraction(m_k_lg_Targets[k]);
473 }
474
475 // Normalize the mole fractions to unity:
476 for (size_t k = 0; k < m_X_targets.size(); k++) {
477 m_X_targets[k] = m_X_targets[k] / tmp_sum;
478 }
479}
480
482{
484 m_totalCrossSectionEdge.assign(m_points + 1, 0.0);
485 for (size_t k = 0; k < m_phase->nCollisions(); k++) {
486 auto x = m_phase->collisionRate(k)->energyLevels();
487 auto y = m_phase->collisionRate(k)->crossSections();
488
489 for (size_t i = 0; i < m_points; i++) {
491 linearInterp(m_gridCenter[i], x, y);
492 }
493 for (size_t i = 0; i < m_points + 1; i++) {
495 linearInterp(m_gridEdge[i], x, y);
496 }
497 }
498}
499
501{
502 m_sigmaElastic.clear();
503 m_sigmaElastic.resize(m_points, 0.0);
504 for (size_t k : m_phase->kElastic()) {
505 auto x = m_phase->collisionRate(k)->energyLevels();
506 auto y = m_phase->collisionRate(k)->crossSections();
507 // Note:
508 // moleFraction(m_kTargets[k]) <=> m_X_targets[m_klocTargets[k]]
509 double mass_ratio = ElectronMass / (m_phase->molecularWeight(m_kTargets[k]) / Avogadro);
510 for (size_t i = 0; i < m_points; i++) {
511 m_sigmaElastic[i] += 2.0 * mass_ratio * m_X_targets[m_klocTargets[k]] *
512 linearInterp(m_gridEdge[i], x, y);
513 }
514 }
515}
516
517void EEDFTwoTermApproximation::setGridCache()
518{
519 m_sigma.clear();
520 m_sigma.resize(m_phase->nCollisions());
521 m_eps.clear();
522 m_eps.resize(m_phase->nCollisions());
523 m_j.clear();
524 m_j.resize(m_phase->nCollisions());
525 m_i.clear();
526 m_i.resize(m_phase->nCollisions());
527 for (size_t k = 0; k < m_phase->nCollisions(); k++) {
528 auto& collision = m_phase->collisionRate(k);
529 auto x = collision->energyLevels();
530 auto y = collision->crossSections();
531 vector<double> eps1(m_points + 1);
532 int shiftFactor = (collision->kind() == "ionization") ? 2 : 1;
533
534 for (size_t i = 0; i < m_points + 1; i++) {
535 eps1[i] = clip(shiftFactor * m_gridEdge[i] + collision->threshold(),
536 m_gridEdge[0] + 1e-9, m_gridEdge[m_points] - 1e-9);
537 }
538 vector<double> nodes = eps1;
539 for (size_t i = 0; i < m_points + 1; i++) {
540 if (m_gridEdge[i] >= eps1[0] && m_gridEdge[i] <= eps1[m_points]) {
541 nodes.push_back(m_gridEdge[i]);
542 }
543 }
544 for (size_t i = 0; i < x.size(); i++) {
545 if (x[i] >= eps1[0] && x[i] <= eps1[m_points]) {
546 nodes.push_back(x[i]);
547 }
548 }
549
550 std::sort(nodes.begin(), nodes.end());
551 auto last = std::unique(nodes.begin(), nodes.end());
552 nodes.resize(std::distance(nodes.begin(), last));
553 vector<double> sigma0(nodes.size());
554 for (size_t i = 0; i < nodes.size(); i++) {
555 sigma0[i] = linearInterp(nodes[i], x, y);
556 }
557
558 // search position of cell j
559 for (size_t i = 1; i < nodes.size(); i++) {
560 auto low = std::lower_bound(m_gridEdge.begin(), m_gridEdge.end(), nodes[i]);
561 m_j[k].push_back(low - m_gridEdge.begin() - 1);
562 }
563
564 // search position of cell i
565 for (size_t i = 1; i < nodes.size(); i++) {
566 auto low = std::lower_bound(eps1.begin(), eps1.end(), nodes[i]);
567 m_i[k].push_back(low - eps1.begin() - 1);
568 }
569
570 // construct sigma
571 for (size_t i = 0; i < nodes.size() - 1; i++) {
572 m_sigma[k].push_back({sigma0[i], sigma0[i+1]});
573 }
574
575 // construct eps
576 for (size_t i = 0; i < nodes.size() - 1; i++) {
577 m_eps[k].push_back({nodes[i], nodes[i+1]});
578 }
579
580 // construct sigma_offset
581 vector<double> x_offset(collision->energyLevels().begin(),
582 collision->energyLevels().end());
583 for (auto& element : x_offset) {
584 element -= collision->threshold();
585 }
586 }
587}
588
589double EEDFTwoTermApproximation::norm(const Eigen::VectorXd& f, const Eigen::VectorXd& grid)
590{
591 string m_quadratureMethod = "simpson";
592 Eigen::VectorXd p(f.size());
593 for (int i = 0; i < f.size(); i++) {
594 p[i] = f(i) * pow(grid[i], 0.5);
595 }
596 return numericalQuadrature(m_quadratureMethod, p, grid);
597}
598
599}
EEDF Two-Term approximation solver.
Header for plasma reaction rates parameterized by electron collision cross section and electron energ...
Header file for class PlasmaPhase.
Base class for exceptions thrown by Cantera classes.
double m_rtol
Error tolerance for convergence.
Eigen::VectorXd m_f0
normalized electron energy distribution function
vector< double > m_gridEdge
Grid of electron energy (cell boundary i-1/2) [eV].
vector< vector< size_t > > m_i
Location of cell i for grid cache.
void calculateTotalCrossSection()
Compute the total (elastic + inelastic) cross section.
vector< vector< size_t > > m_j
Location of cell j for grid cache.
vector< double > m_X_targets_prev
Previous mole fraction of targets used to compute eedf.
vector< vector< vector< double > > > m_eps
The energy boundaries of the overlap of cell i and j.
vector< vector< vector< double > > > m_sigma
Cross section at the boundaries of the overlap of cell i and j.
vector< size_t > m_k_lg_Targets
Local to global indices.
bool m_first_call
First call to calculateDistributionFunction.
void converge(Eigen::VectorXd &f0)
Iterate f0 (EEDF) until convergence.
double m_delta0
Formerly options for the EEDF solver.
double electronDiffusivity(const Eigen::VectorXd &f0)
Diffusivity.
PlasmaPhase * m_phase
Pointer to the PlasmaPhase object used to initialize this object.
double norm(const Eigen::VectorXd &f, const Eigen::VectorXd &grid)
Compute the L1 norm of a function f defined over a given energy grid.
void updateMoleFractions()
Update the vector of species mole fractions.
vector< size_t > m_kTargets
list of target species indices in global Cantera numbering (1 index per cs)
double electronMobility(const Eigen::VectorXd &f0)
Mobility.
Eigen::VectorXd m_gridCenter
Grid of electron energy (cell center) [eV].
size_t m_maxn
Maximum number of iterations.
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 Scharfet...
Eigen::SparseMatrix< double > matrix_Q(span< const double > g, size_t k)
The matrix of scattering-in.
Eigen::VectorXd iterate(const Eigen::VectorXd &f0, double delta)
An iteration of solving electron energy distribution function.
double m_electronMobility
Electron mobility [m²/V·s].
size_t m_points
The number of points in the EEDF grid.
vector< size_t > m_kOthers
Indices of species which has no cross-section data.
double m_factorM
The factor for step size change.
void calculateTotalElasticCrossSection()
Compute the total elastic collision cross section.
void updateCrossSections()
Update the total cross sections based on the current state.
void initSpeciesIndexCrossSections()
Initialize species indices associated with cross-section data.
double m_init_kTe
The initial electron temperature [eV].
double netProductionFrequency(const Eigen::VectorXd &f0)
Reduced net production frequency.
std::string m_firstguess
The first guess for the EEDF.
std::string m_growth
The growth model of EEDF.
vector< double > m_totalCrossSectionEdge
Total electron cross section on the cell boundary (i-1/2) of energy grid.
Eigen::SparseMatrix< double > matrix_P(span< const double > g, size_t k)
The matrix of scattering-out.
vector< double > m_f0_edge
EEDF at grid edges (cell boundaries)
vector< size_t > m_klocTargets
list of target species indices in local X EEDF numbering (1 index per cs)
vector< double > m_X_targets
Mole fraction of targets.
double integralPQ(double a, double b, double u0, double u1, double g, double x0)
The integral in [a, b] of assuming that u is linear with u(a) = u0 and u(b) = u1.
vector< double > m_totalCrossSectionCenter
Total electron cross section on the cell center of energy grid.
vector< double > vector_g(const Eigen::VectorXd &f0)
Vector g is used by matrix_P() and matrix_Q().
vector< double > m_sigmaElastic
vector of total elastic cross section weighted with mass ratio
virtual double molarDensity() const
Molar density (kmol/m^3).
Definition Phase.cpp:597
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:247
double temperature() const
Temperature (K).
Definition Phase.h:586
double moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition Phase.cpp:457
double molecularWeight(size_t k) const
Molecular weight of species k.
Definition Phase.cpp:398
Base class for handling plasma properties, specifically focusing on the electron energy distribution.
Definition PlasmaPhase.h:84
double electricFieldFrequency() const
Get the frequency of the applied electric field [Hz].
size_t nCollisions() const
Number of electron collision cross sections.
double electricField() const
Get the applied electric field strength [V/m].
const vector< size_t > & kElastic() const
Get the indices for elastic electron collisions.
const shared_ptr< ElectronCollisionPlasmaRate > collisionRate(size_t i) const
Get the ElectronCollisionPlasmaRate object associated with electron collision i.
size_t targetIndex(size_t i) const
target of a specific process
const vector< size_t > & kInelastic() const
Get the indicies for inelastic electron collisions.
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
Header for a file containing miscellaneous numerical functions.
T clip(const T &value, const T &lower, const T &upper)
Clip value such that lower <= value <= upper.
Definition global.h:326
double numericalQuadrature(const string &method, const Eigen::ArrayXd &f, const Eigen::ArrayXd &x)
Numerical integration of a function.
Definition funcs.cpp:116
double simpson(const Eigen::ArrayXd &f, const Eigen::ArrayXd &x)
Numerical integration of a function using Simpson's rule with flexibility of taking odd and even numb...
Definition funcs.cpp:91
double linearInterp(double x, span< const double > xpts, span< const double > fpts)
Linearly interpolate a function defined on a discrete grid.
Definition funcs.cpp:13
const double Boltzmann
Boltzmann constant [J/K].
Definition ct_defs.h:87
const double Avogadro
Avogadro's Number [number/kmol].
Definition ct_defs.h:84
const double ElectronCharge
Elementary charge [C].
Definition ct_defs.h:93
const double Pi
Pi.
Definition ct_defs.h:71
const double ElectronMass
Electron Mass [kg].
Definition ct_defs.h:114
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition global.h:263
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
MappedVector asVectorXd(vector< double > &v)
Convenience wrapper for accessing std::vector as an Eigen VectorXd.
Definition eigen_dense.h:60
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
Definition eigen_dense.h:46