Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Medium.cc
Go to the documentation of this file.
1#include <algorithm>
2#include <cmath>
3#include <iomanip>
4#include <iostream>
5#include <utility>
6
9#include "Garfield/Medium.hh"
10#include "Garfield/Numerics.hh"
11#include "Garfield/Random.hh"
12
13namespace {
14
15void PrintNotImplemented(const std::string& cls, const std::string& fcn) {
16 std::cerr << cls << "::" << fcn << ": Function is not implemented.\n";
17}
18
19void PrintOutOfRange(const std::string& cls, const std::string& fcn,
20 const size_t i, const size_t j, const size_t k) {
21 std::cerr << cls << "::" << fcn << ": Index (" << i << ", " << j << ", " << k
22 << ") out of range.\n";
23}
24
25void PrintDataNotAvailable(const std::string& cls, const std::string& fcn) {
26 std::cerr << cls << "::" << fcn << ": Data not available.\n";
27}
28
29bool CheckFields(const std::vector<double>& fields, const std::string& hdr,
30 const std::string& lbl) {
31 if (fields.empty()) {
32 std::cerr << hdr << ": Number of " << lbl << " must be > 0.\n";
33 return false;
34 }
35
36 // Make sure the values are not negative.
37 if (fields.front() < 0.) {
38 std::cerr << hdr << ": " << lbl << " must be >= 0.\n";
39 return false;
40 }
41
42 const size_t nEntries = fields.size();
43 // Make sure the values are in strictly monotonic, ascending order.
44 if (nEntries > 1) {
45 for (size_t i = 1; i < nEntries; ++i) {
46 if (fields[i] <= fields[i - 1]) {
47 std::cerr << hdr << ": " << lbl << " are not in ascending order.\n";
48 return false;
49 }
50 }
51 }
52 return true;
53}
54}
55
56namespace Garfield {
57
58int Medium::m_idCounter = -1;
59
60Medium::Medium() : m_id(++m_idCounter) {
61 // Initialise the transport tables.
62 m_bFields.assign(1, 0.);
63 m_bAngles.assign(1, HalfPi);
64
65 // Set the default grid.
66 SetFieldGrid(100., 100000., 20, true, 0., 0., 1, HalfPi, HalfPi, 1);
67}
68
70
71void Medium::SetTemperature(const double t) {
72 if (t <= 0.) {
73 std::cerr << m_className << "::SetTemperature:\n"
74 << " Temperature [K] must be greater than zero.\n";
75 return;
76 }
77 m_temperature = t;
78 m_isChanged = true;
79}
80
81void Medium::SetPressure(const double p) {
82 if (p <= 0.) {
83 std::cerr << m_className << "::SetPressure:\n"
84 << " Pressure [Torr] must be greater than zero.\n";
85 return;
86 }
87 m_pressure = p;
88 m_isChanged = true;
89}
90
91void Medium::SetDielectricConstant(const double eps) {
92 if (eps < 1.) {
93 std::cerr << m_className << "::SetDielectricConstant:\n"
94 << " Dielectric constant must be >= 1.\n";
95 return;
96 }
97 m_epsilon = eps;
98 m_isChanged = true;
99}
100
102 return m_density * AtomicMassUnit * m_a;
103}
104
105void Medium::GetComponent(const unsigned int i, std::string& label, double& f) {
106 if (i >= m_nComponents) {
107 std::cerr << m_className << "::GetComponent: Index out of range.\n";
108 }
109
110 label = m_name;
111 f = 1.;
112}
113
114void Medium::SetAtomicNumber(const double z) {
115 if (z < 1.) {
116 std::cerr << m_className << "::SetAtomicNumber:\n"
117 << " Atomic number must be >= 1.\n";
118 return;
119 }
120 m_z = z;
121 m_isChanged = true;
122}
123
124void Medium::SetAtomicWeight(const double a) {
125 if (a <= 0.) {
126 std::cerr << m_className << "::SetAtomicWeight:\n"
127 << " Atomic weight must be greater than zero.\n";
128 return;
129 }
130 m_a = a;
131 m_isChanged = true;
132}
133
134void Medium::SetNumberDensity(const double n) {
135 if (n <= 0.) {
136 std::cerr << m_className << "::SetNumberDensity:\n"
137 << " Density [cm-3] must be greater than zero.\n";
138 return;
139 }
140 m_density = n;
141 m_isChanged = true;
142}
143
144void Medium::SetMassDensity(const double rho) {
145 if (rho <= 0.) {
146 std::cerr << m_className << "::SetMassDensity:\n"
147 << " Density [g/cm3] must be greater than zero.\n";
148 return;
149 }
150
151 if (m_a <= 0.) {
152 std::cerr << m_className << "::SetMassDensity:\n"
153 << " Atomic weight is not defined.\n";
154 return;
155 }
156 m_density = rho / (AtomicMassUnit * m_a);
157 m_isChanged = true;
158}
159
160bool Medium::Velocity(const double ex, const double ey, const double ez,
161 const double bx, const double by, const double bz,
162 const std::vector<std::vector<std::vector<double> > >& velE,
163 const std::vector<std::vector<std::vector<double> > >& velB,
164 const std::vector<std::vector<std::vector<double> > >& velX,
165 const double q, double& vx, double& vy, double& vz) const {
166
167 vx = vy = vz = 0.;
168 // Make sure there is at least a table of velocities along E.
169 if (velE.empty()) return false;
170
171 // Compute the magnitude of the electric field.
172 const double e = sqrt(ex * ex + ey * ey + ez * ez);
173 const double e0 = ScaleElectricField(e);
174 if (e < Small || e0 < Small) return false;
175
176 // Compute the magnitude of the magnetic field.
177 const double b = sqrt(bx * bx + by * by + bz * bz);
178 // Compute the angle between B field and E field.
179 const double ebang = GetAngle(ex, ey, ez, bx, by, bz, e, b);
180
181 // Calculate the velocity along E.
182 double ve = 0.;
183 if (!Interpolate(e0, b, ebang, velE, ve, m_intpVel, m_extrVel)) {
184 std::cerr << m_className << "::Velocity: Interpolation along E failed.\n";
185 return false;
186 }
187 if (b < Small) {
188 // No magnetic field.
189 const double mu = q * ve / e;
190 vx = mu * ex;
191 vy = mu * ey;
192 vz = mu * ez;
193 return true;
194 } else if (velX.empty() || velB.empty()) {
195 // Magnetic field, velocities along ExB, Bt not available.
196 const double mu = q * ve / e;
197 const double mu2 = mu * mu;
198 const double eb = bx * ex + by * ey + bz * ez;
199 const double f = mu / (1. + mu2 * b * b);
200 vx = f * (ex + mu * (ey * bz - ez * by) + mu2 * bx * eb);
201 vy = f * (ey + mu * (ez * bx - ex * bz) + mu2 * by * eb);
202 vz = f * (ez + mu * (ex * by - ey * bx) + mu2 * bz * eb);
203 return true;
204 }
205
206 // Magnetic field, velocities along ExB and Bt available.
207 // Compute unit vectors along E, E x B and Bt.
208 double ue[3] = {ex / e, ey / e, ez / e};
209 double uexb[3] = {ey * bz - ez * by, ez * bx - ex * bz, ex * by - ey * bx};
210 const double exb =
211 sqrt(uexb[0] * uexb[0] + uexb[1] * uexb[1] + uexb[2] * uexb[2]);
212 if (exb > 0.) {
213 uexb[0] /= exb;
214 uexb[1] /= exb;
215 uexb[2] /= exb;
216 } else {
217 uexb[0] = ue[0];
218 uexb[1] = ue[1];
219 uexb[2] = ue[2];
220 }
221
222 double ubt[3] = {uexb[1] * ez - uexb[2] * ey, uexb[2] * ex - uexb[0] * ez,
223 uexb[0] * ey - uexb[1] * ex};
224 const double bt = sqrt(ubt[0] * ubt[0] + ubt[1] * ubt[1] + ubt[2] * ubt[2]);
225 if (bt > 0.) {
226 ubt[0] /= bt;
227 ubt[1] /= bt;
228 ubt[2] /= bt;
229 } else {
230 ubt[0] = ue[0];
231 ubt[1] = ue[1];
232 ubt[2] = ue[2];
233 }
234
235 if (m_debug) {
236 std::cout << std::setprecision(5);
237 std::cout << m_className << "::Velocity:\n"
238 << " unit vector along E: (" << ue[0] << ", " << ue[1]
239 << ", " << ue[2] << ")\n";
240 std::cout << " unit vector along E x B: (" << uexb[0] << ", "
241 << uexb[1] << ", " << uexb[2] << ")\n";
242 std::cout << " unit vector along Bt: (" << ubt[0] << ", " << ubt[1]
243 << ", " << ubt[2] << ")\n";
244 }
245
246 // Calculate the velocities in all directions.
247 double vexb = 0.;
248 if (!Interpolate(e0, b, ebang, velX, vexb, m_intpVel, m_extrVel)) {
249 std::cerr << m_className << "::Velocity: Interpolation along ExB failed.\n";
250 return false;
251 }
252 double vbt = 0.;
253 if (!Interpolate(e0, b, ebang, velB, vbt, m_intpVel, m_extrVel)) {
254 std::cerr << m_className << "::Velocity: Interpolation along Bt failed.\n";
255 return false;
256 }
257 if (ex * bx + ey * by + ez * bz > 0.) {
258 vbt = fabs(vbt);
259 } else {
260 vbt = -fabs(vbt);
261 }
262 vx = q * (ve * ue[0] + q * q * vbt * ubt[0] + q * vexb * uexb[0]);
263 vy = q * (ve * ue[1] + q * q * vbt * ubt[1] + q * vexb * uexb[1]);
264 vz = q * (ve * ue[2] + q * q * vbt * ubt[2] + q * vexb * uexb[2]);
265
266 return true;
267}
268
269bool Medium::Diffusion(const double ex, const double ey, const double ez,
270 const double bx, const double by, const double bz,
271 const std::vector<std::vector<std::vector<double> > >& difL,
272 const std::vector<std::vector<std::vector<double> > >& difT,
273 double& dl, double& dt) const {
274
275 dl = dt = 0.;
276 // Compute the magnitude of the electric field.
277 const double e = sqrt(ex * ex + ey * ey + ez * ez);
278 const double e0 = ScaleElectricField(e);
279 if (e < Small || e0 < Small) return true;
280
281 // Compute the magnitude of the magnetic field.
282 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
283 // Compute the angle between B field and E field.
284 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
285
286 // Interpolate.
287 if (!difL.empty()) {
288 if (!Interpolate(e0, b, ebang, difL, dl, m_intpDif, m_extrDif)) dl = 0.;
289 }
290 if (!difT.empty()) {
291 if (!Interpolate(e0, b, ebang, difT, dt, m_intpDif, m_extrDif)) dt = 0.;
292 }
293
294 // If no data available, calculate
295 // the diffusion coefficients using the Einstein relation
296 if (difL.empty() || difT.empty()) {
297 const double d = sqrt(2. * BoltzmannConstant * m_temperature / e);
298 if (difL.empty()) dl = d;
299 if (difT.empty()) dt = d;
300 }
301 // Verify values and apply scaling.
302 dl = ScaleDiffusion(std::max(dl, 0.));
303 dt = ScaleDiffusion(std::max(dt, 0.));
304 return true;
305}
306
307bool Medium::Diffusion(const double ex, const double ey, const double ez,
308 const double bx, const double by, const double bz,
309 const std::vector<std::vector<std::vector<std::vector<double> > > >& diff,
310 double cov[3][3]) const {
311
312 // Initialise the tensor.
313 cov[0][0] = cov[0][1] = cov[0][2] = 0.;
314 cov[1][0] = cov[1][1] = cov[1][2] = 0.;
315 cov[2][0] = cov[2][1] = cov[2][2] = 0.;
316
317 if (diff.empty()) return false;
318
319 // Compute the magnitude of the electric field.
320 const double e = sqrt(ex * ex + ey * ey + ez * ez);
321 const double e0 = ScaleElectricField(e);
322 if (e < Small || e0 < Small) return true;
323
324 // Compute the magnitude of the magnetic field.
325 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
326 // Compute the angle between B field and E field.
327 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
328
329 for (int j = 0; j < 6; ++j) {
330 // Interpolate.
331 double y = 0.;
332 if (!Interpolate(e0, b, ebang, diff[j], y, m_intpDif, m_extrDif)) y = 0.;
333 // Apply scaling.
335 if (j < 3) {
336 cov[j][j] = y;
337 } else if (j == 3) {
338 cov[0][1] = cov[1][0] = y;
339 } else if (j == 4) {
340 cov[0][2] = cov[2][0] = y;
341 } else if (j == 5) {
342 cov[1][2] = cov[2][1] = y;
343 }
344 }
345 return true;
346}
347
348bool Medium::Alpha(const double ex, const double ey, const double ez,
349 const double bx, const double by, const double bz,
350 const std::vector<std::vector<std::vector<double> > >& tab,
351 unsigned int intp, const unsigned int thr,
352 const std::pair<unsigned int, unsigned int>& extr,
353 double& alpha) const {
354
355 alpha = 0.;
356 if (tab.empty()) return false;
357
358 // Compute the magnitude of the electric field.
359 const double e = sqrt(ex * ex + ey * ey + ez * ez);
360 const double e0 = ScaleElectricField(e);
361 if (e < Small || e0 < Small) return true;
362
363 // Compute the magnitude of the magnetic field.
364 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
365 // Compute the angle between B field and E field.
366 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
367
368 // Interpolate.
369 if (e0 < m_eFields[thr]) intp = 1;
370 if (!Interpolate(e0, b, ebang, tab, alpha, intp, extr)) alpha = -30.;
371 if (alpha < -20.) {
372 alpha = 0.;
373 } else {
374 alpha = exp(alpha);
375 }
376 return true;
377}
378
379bool Medium::ElectronVelocity(const double ex, const double ey, const double ez,
380 const double bx, const double by, const double bz,
381 double& vx, double& vy, double& vz) {
382
383 return Velocity(ex, ey, ez, bx, by, bz, m_eVelE, m_eVelB, m_eVelX, -1.,
384 vx, vy, vz);
385}
386
387bool Medium::ElectronDiffusion(const double ex, const double ey,
388 const double ez, const double bx,
389 const double by, const double bz, double& dl,
390 double& dt) {
391
392 return Diffusion(ex, ey, ez, bx, by, bz, m_eDifL, m_eDifT, dl, dt);
393}
394
395bool Medium::ElectronDiffusion(const double ex, const double ey,
396 const double ez, const double bx,
397 const double by, const double bz,
398 double cov[3][3]) {
399
400 return Diffusion(ex, ey, ez, bx, by, bz, m_eDifM, cov);
401}
402
403bool Medium::ElectronTownsend(const double ex, const double ey, const double ez,
404 const double bx, const double by, const double bz,
405 double& alpha) {
406
407 if (!Alpha(ex, ey, ez, bx, by, bz, m_eAlp, m_intpAlp, m_eThrAlp, m_extrAlp,
408 alpha)) {
409 return false;
410 }
411 // Apply scaling.
412 alpha = ScaleTownsend(alpha);
413 return true;
414}
415
416bool Medium::ElectronAttachment(const double ex, const double ey,
417 const double ez, const double bx,
418 const double by, const double bz, double& eta) {
419
420 if (!Alpha(ex, ey, ez, bx, by, bz, m_eAtt, m_intpAtt, m_eThrAtt, m_extrAtt,
421 eta)) {
422 return false;
423 }
424 // Apply scaling.
425 eta = ScaleAttachment(eta);
426 return true;
427}
428
429bool Medium::ElectronLorentzAngle(const double ex, const double ey,
430 const double ez, const double bx,
431 const double by, const double bz,
432 double& lor) {
433 lor = 0.;
434 if (m_eLor.empty()) return false;
435
436 // Compute the magnitude of the electric field.
437 const double e = sqrt(ex * ex + ey * ey + ez * ez);
438 const double e0 = ScaleElectricField(e);
439 if (e < Small || e0 < Small) return true;
440
441 // Compute the magnitude of the magnetic field.
442 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
443 // Compute the angle between B field and E field.
444 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
445
446 // Interpolate.
447 if (!Interpolate(e0, b, ebang, m_eLor, lor, m_intpLor, m_extrLor)) lor = 0.;
448 // Apply scaling.
449 lor = ScaleLorentzAngle(lor);
450 return true;
451}
452
454 if (m_eVelE.empty()) return -1.;
455 return m_eVelE[0][0][0] / UnScaleElectricField(m_eFields[0]);
456}
457
458double Medium::GetElectronEnergy(const double px, const double py,
459 const double pz, double& vx, double& vy,
460 double& vz, const int band) {
461 if (band > 0) {
462 std::cerr << m_className << "::GetElectronEnergy:\n";
463 std::cerr << " Unknown band index.\n";
464 }
465
466 vx = SpeedOfLight * px / ElectronMass;
467 vy = SpeedOfLight * py / ElectronMass;
468 vz = SpeedOfLight * pz / ElectronMass;
469
470 return 0.5 * (px * px + py * py + pz * pz) / ElectronMass;
471}
472
473void Medium::GetElectronMomentum(const double e, double& px, double& py,
474 double& pz, int& band) {
475 const double p = sqrt(2. * ElectronMass * e) / SpeedOfLight;
476 RndmDirection(px, py, pz, p);
477 band = -1;
478}
479
480double Medium::GetElectronNullCollisionRate(const int /*band*/) {
481 if (m_debug) PrintNotImplemented(m_className, "GetElectronNullCollisionRate");
482 return 0.;
483}
484
485double Medium::GetElectronCollisionRate(const double /*e*/,
486 const int /*band*/) {
487 if (m_debug) PrintNotImplemented(m_className, "GetElectronCollisionRate");
488 return 0.;
489}
490
492 const double e, int& type, int& level, double& e1, double& dx, double& dy,
493 double& dz, std::vector<std::pair<int, double> >& /*secondaries*/,
494 int& ndxc, int& band) {
495 type = level = -1;
496 e1 = e;
497 ndxc = band = 0;
498 RndmDirection(dx, dy, dz);
499
500 if (m_debug) PrintNotImplemented(m_className, "GetElectronCollision");
501 return false;
502}
503
504bool Medium::GetDeexcitationProduct(const unsigned int /*i*/, double& t,
505 double& s, int& type,
506 double& energy) const {
507 if (m_debug) PrintNotImplemented(m_className, "GetDeexcitationProduct");
508 t = s = energy = 0.;
509 type = 0;
510 return false;
511}
512
513bool Medium::HoleVelocity(const double ex, const double ey, const double ez,
514 const double bx, const double by, const double bz,
515 double& vx, double& vy, double& vz) {
516
517 return Velocity(ex, ey, ez, bx, by, bz, m_hVelE, m_hVelB, m_hVelX, +1.,
518 vx, vy, vz);
519}
520
521bool Medium::HoleDiffusion(const double ex, const double ey, const double ez,
522 const double bx, const double by, const double bz,
523 double& dl, double& dt) {
524 return Diffusion(ex, ey, ez, bx, by, bz, m_hDifL, m_hDifT, dl, dt);
525}
526
527bool Medium::HoleDiffusion(const double ex, const double ey, const double ez,
528 const double bx, const double by, const double bz,
529 double cov[3][3]) {
530
531 return Diffusion(ex, ey, ez, bx, by, bz, m_hDifM, cov);
532}
533
534bool Medium::HoleTownsend(const double ex, const double ey, const double ez,
535 const double bx, const double by, const double bz,
536 double& alpha) {
537
538 if (!Alpha(ex, ey, ez, bx, by, bz, m_hAlp, m_intpAlp, m_hThrAlp, m_extrAlp,
539 alpha)) {
540 return false;
541 }
542 // Apply scaling.
543 alpha = ScaleTownsend(alpha);
544 return true;
545}
546
547bool Medium::HoleAttachment(const double ex, const double ey, const double ez,
548 const double bx, const double by, const double bz,
549 double& eta) {
550
551 if (!Alpha(ex, ey, ez, bx, by, bz, m_hAtt, m_intpAtt, m_hThrAtt, m_extrAtt,
552 eta)) {
553 return false;
554 }
555 // Apply scaling.
556 eta = ScaleAttachment(eta);
557 return true;
558}
559
561 if (m_hVelE.empty()) return -1.;
562 return m_hVelE[0][0][0] / UnScaleElectricField(m_eFields[0]);
563}
564
565bool Medium::IonVelocity(const double ex, const double ey, const double ez,
566 const double bx, const double by, const double bz,
567 double& vx, double& vy, double& vz) {
568 vx = vy = vz = 0.;
569 if (m_iMob.empty()) return false;
570 // Compute the magnitude of the electric field.
571 const double e = sqrt(ex * ex + ey * ey + ez * ez);
572 const double e0 = ScaleElectricField(e);
573 if (e < Small || e0 < Small) return true;
574 // Compute the magnitude of the electric field.
575 const double b = sqrt(bx * bx + by * by + bz * bz);
576
577 // Compute the angle between B field and E field.
578 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
579 double mu = 0.;
580 if (!Interpolate(e0, b, ebang, m_iMob, mu, m_intpMob, m_extrMob)) mu = 0.;
581
582 constexpr double q = 1.;
583 mu *= q;
584 if (b < Small) {
585 vx = mu * ex;
586 vy = mu * ey;
587 vz = mu * ez;
588 } else {
589 const double eb = bx * ex + by * ey + bz * ez;
590 const double mu2 = mu * mu;
591 const double f = mu / (1. + mu2 * b * b);
592 vx = f * (ex + mu * (ey * bz - ez * by) + mu2 * bx * eb);
593 vy = f * (ey + mu * (ez * bx - ex * bz) + mu2 * by * eb);
594 vz = f * (ez + mu * (ex * by - ey * bx) + mu2 * bz * eb);
595 }
596
597 return true;
598}
599
600bool Medium::IonDiffusion(const double ex, const double ey, const double ez,
601 const double bx, const double by, const double bz,
602 double& dl, double& dt) {
603
604 return Diffusion(ex, ey, ez, bx, by, bz, m_iDifL, m_iDifT, dl, dt);
605}
606
607bool Medium::IonDissociation(const double ex, const double ey, const double ez,
608 const double bx, const double by, const double bz,
609 double& diss) {
610
611 if (!Alpha(ex, ey, ez, bx, by, bz, m_iDis, m_intpDis, m_iThrDis, m_extrDis,
612 diss)) {
613 return false;
614 }
615 // Apply scaling.
616 diss = ScaleDissociation(diss);
617 return true;
618}
619
621 return m_iMob.empty() ? -1. : m_iMob[0][0][0];
622}
623
624bool Medium::GetOpticalDataRange(double& emin, double& emax,
625 const unsigned int i) {
626 if (i >= m_nComponents) {
627 std::cerr << m_className << "::GetOpticalDataRange: Index out of range.\n";
628 return false;
629 }
630
631 if (m_debug) PrintNotImplemented(m_className, "GetOpticalDataRange");
632 emin = emax = 0.;
633 return false;
634}
635
636bool Medium::GetDielectricFunction(const double e, double& eps1, double& eps2,
637 const unsigned int i) {
638 if (i >= m_nComponents) {
639 std::cerr << m_className << "::GetDielectricFunction: Index out of range.\n";
640 return false;
641 }
642
643 if (e < 0.) {
644 std::cerr << m_className << "::GetDielectricFunction: Energy must be > 0.\n";
645 return false;
646 }
647
648 if (m_debug) PrintNotImplemented(m_className, "GetDielectricFunction");
649 eps1 = 1.;
650 eps2 = 0.;
651 return false;
652}
653
654bool Medium::GetPhotoAbsorptionCrossSection(const double e, double& sigma,
655 const unsigned int i) {
656 if (i >= m_nComponents) {
657 std::cerr << m_className << "::GetPhotoAbsorptionCrossSection:\n";
658 std::cerr << " Component " << i << " does not exist.\n";
659 return false;
660 }
661
662 if (e < 0.) {
663 std::cerr << m_className << "::GetPhotoAbsorptionCrossSection:\n";
664 std::cerr << " Energy must be > 0.\n";
665 return false;
666 }
667
668 if (m_debug) {
669 PrintNotImplemented(m_className, "GetPhotoAbsorptionCrossSection");
670 }
671 sigma = 0.;
672 return false;
673}
674
675double Medium::GetPhotonCollisionRate(const double e) {
676 double sigma = 0.;
677 if (!GetPhotoAbsorptionCrossSection(e, sigma)) return 0.;
678
679 return sigma * m_density * SpeedOfLight;
680}
681
682bool Medium::GetPhotonCollision(const double e, int& type, int& level,
683 double& e1, double& ctheta, int& nsec,
684 double& esec) {
685 type = level = -1;
686 e1 = e;
687 ctheta = 1.;
688 nsec = 0;
689 esec = 0.;
690 return false;
691}
692
693void Medium::SetFieldGrid(double emin, double emax, const size_t ne, bool logE,
694 double bmin, double bmax, const size_t nb,
695 double amin, double amax, const size_t na) {
696 // Check if the requested E-field range makes sense.
697 if (ne <= 0) {
698 std::cerr << m_className << "::SetFieldGrid:\n"
699 << " Number of E-fields must be > 0.\n";
700 return;
701 }
702
703 if (emin < 0. || emax < 0.) {
704 std::cerr << m_className << "::SetFieldGrid:\n"
705 << " Electric fields must be positive.\n";
706 return;
707 }
708
709 if (emax < emin) {
710 std::cerr << m_className << "::SetFieldGrid: Swapping min./max. E-field.\n";
711 std::swap(emin, emax);
712 }
713
714 double estep = 0.;
715 if (logE) {
716 // Logarithmic scale
717 if (emin < Small) {
718 std::cerr << m_className << "::SetFieldGrid:\n"
719 << " Min. E-field must be non-zero for log. scale.\n";
720 return;
721 }
722 if (ne == 1) {
723 std::cerr << m_className << "::SetFieldGrid:\n"
724 << " Number of E-fields must be > 1 for log. scale.\n";
725 return;
726 }
727 estep = pow(emax / emin, 1. / (ne - 1.));
728 } else {
729 // Linear scale
730 if (ne > 1) estep = (emax - emin) / (ne - 1.);
731 }
732
733 // Check if the requested B-field range makes sense.
734 if (nb <= 0) {
735 std::cerr << m_className << "::SetFieldGrid:\n"
736 << " Number of B-fields must be > 0.\n";
737 return;
738 }
739 if (bmax < 0. || bmin < 0.) {
740 std::cerr << m_className << "::SetFieldGrid:\n"
741 << " Magnetic fields must be positive.\n";
742 return;
743 }
744 if (bmax < bmin) {
745 std::cerr << m_className << "::SetFieldGrid: Swapping min./max. B-field.\n";
746 std::swap(bmin, bmax);
747 }
748
749 const double bstep = nb > 1 ? (bmax - bmin) / (nb - 1.) : 0.;
750
751 // Check if the requested angular range makes sense.
752 if (na <= 0) {
753 std::cerr << m_className << "::SetFieldGrid:\n"
754 << " Number of angles must be > 0.\n";
755 return;
756 }
757 if (amax < 0. || amin < 0.) {
758 std::cerr << m_className << "::SetFieldGrid:"
759 << " Angles must be positive.\n";
760 return;
761 }
762 if (amax < amin) {
763 std::cerr << m_className << "::SetFieldGrid: Swapping min./max. angle.\n";
764 std::swap(amin, amax);
765 }
766 const double astep = na > 1 ? (amax - amin) / (na - 1.) : 0;
767
768 // Setup the field grids.
769 std::vector<double> eFields(ne);
770 std::vector<double> bFields(nb);
771 std::vector<double> bAngles(na);
772 for (size_t i = 0; i < ne; ++i) {
773 eFields[i] = logE ? emin * pow(estep, i) : emin + i * estep;
774 }
775 for (size_t i = 0; i < nb; ++i) {
776 bFields[i] = bmin + i * bstep;
777 }
778 for (size_t i = 0; i < na; ++i) {
779 bAngles[i] = amin + i * astep;
780 }
781 SetFieldGrid(eFields, bFields, bAngles);
782}
783
784void Medium::SetFieldGrid(const std::vector<double>& efields,
785 const std::vector<double>& bfields,
786 const std::vector<double>& angles) {
787 const std::string hdr = m_className + "::SetFieldGrid";
788 if (!CheckFields(efields, hdr, "E-fields")) return;
789 if (!CheckFields(bfields, hdr, "B-fields")) return;
790 if (!CheckFields(angles, hdr, "angles")) return;
791
792 if (m_debug) {
793 std::cout << m_className << "::SetFieldGrid:\n E-fields:\n";
794 for (const auto efield : efields) std::cout << " " << efield << "\n";
795 std::cout << " B-fields:\n";
796 for (const auto bfield : bfields) std::cout << " " << bfield << "\n";
797 std::cout << " Angles:\n";
798 for (const auto angle : angles) std::cout << " " << angle << "\n";
799 }
800
801 // Clone the existing tables.
802 // Electrons
803 Clone(m_eVelE, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
804 "electron velocity along E");
805 Clone(m_eVelB, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
806 "electron velocity along Bt");
807 Clone(m_eVelX, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
808 "electron velocity along ExB");
809 Clone(m_eDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
810 "electron longitudinal diffusion");
811 Clone(m_eDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
812 "electron transverse diffusion");
813 Clone(m_eAlp, efields, bfields, angles, m_intpAlp, m_extrAlp, -30.,
814 "electron Townsend coefficient");
815 Clone(m_eAtt, efields, bfields, angles, m_intpAtt, m_extrAtt, -30.,
816 "electron attachment coefficient");
817 Clone(m_eLor, efields, bfields, angles, m_intpLor, m_extrLor, 0.,
818 "electron Lorentz angle");
819 if (!m_eDifM.empty()) {
820 Clone(m_eDifM, 6, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
821 "electron diffusion tensor");
822 }
823
824 // Holes
825 Clone(m_hVelE, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
826 "hole velocity along E");
827 Clone(m_hVelB, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
828 "hole velocity along Bt");
829 Clone(m_hVelX, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
830 "hole velocity along ExB");
831 Clone(m_hDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
832 "hole longitudinal diffusion");
833 Clone(m_hDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
834 "hole transverse diffusion");
835 Clone(m_hAlp, efields, bfields, angles, m_intpAlp, m_extrAlp, -30.,
836 "hole Townsend coefficient");
837 Clone(m_hAtt, efields, bfields, angles, m_intpAtt, m_extrAtt, -30.,
838 "hole attachment coefficient");
839 if (!m_hDifM.empty()) {
840 Clone(m_hDifM, 6, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
841 "hole diffusion tensor");
842 }
843
844 // Ions
845 Clone(m_iMob, efields, bfields, angles, m_intpMob, m_extrMob, 0.,
846 "ion mobility");
847 Clone(m_iDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
848 "ion longitudinal diffusion");
849 Clone(m_iDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
850 "ion transverse diffusion");
851 Clone(m_iDis, efields, bfields, angles, m_intpDis, m_extrDis, -30.,
852 "ion dissociation");
853
854 if (bfields.size() > 1 || angles.size() > 1) m_tab2d = true;
855 m_eFields = efields;
856 m_bFields = bfields;
857 m_bAngles = angles;
858}
859
860void Medium::GetFieldGrid(std::vector<double>& efields,
861 std::vector<double>& bfields,
862 std::vector<double>& angles) {
863 efields = m_eFields;
864 bfields = m_bFields;
865 angles = m_bAngles;
866}
867
868bool Medium::SetEntry(const size_t i, const size_t j, const size_t k,
869 const std::string& fcn,
870 std::vector<std::vector<std::vector<double> > >& tab,
871 const double val) {
872
873 if (i >= m_eFields.size() || j >= m_bFields.size() || k >= m_bAngles.size()) {
874 PrintOutOfRange(m_className, "Set" + fcn, i, j, k);
875 return false;
876 }
877 if (tab.empty()) {
878 Init(m_eFields.size(), m_bFields.size(), m_bAngles.size(), tab, val);
879 }
880 tab[k][j][i] = val;
881 return true;
882}
883
884bool Medium::GetEntry(const size_t i, const size_t j, const size_t k,
885 const std::string& fcn,
886 const std::vector<std::vector<std::vector<double> > >& tab,
887 double& val) const {
888 val = 0.;
889 if (i >= m_eFields.size() || j >= m_bFields.size() || k >= m_bAngles.size()) {
890 PrintOutOfRange(m_className, "Get" + fcn, i, j, k);
891 return false;
892 }
893 if (tab.empty()) {
894 if (m_debug) PrintDataNotAvailable(m_className, "Get" + fcn);
895 return false;
896 }
897 val = tab[k][j][i];
898 return true;
899}
900
907
912
916}
917
918void Medium::Clone(std::vector<std::vector<std::vector<double> > >& tab,
919 const std::vector<double>& efields,
920 const std::vector<double>& bfields,
921 const std::vector<double>& angles,
922 const unsigned int intp,
923 const std::pair<unsigned int, unsigned int>& extr,
924 const double init, const std::string& lbl) {
925 if (m_debug) {
926 std::cout << m_className << "::Clone: Copying " << lbl << " to new grid.\n";
927 }
928
929 if (tab.empty()) {
930 if (m_debug) std::cout << m_className << "::Clone: Table is empty.\n";
931 return;
932 }
933 // Get the dimensions of the new grid.
934 const auto nE = efields.size();
935 const auto nB = bfields.size();
936 const auto nA = angles.size();
937
938 // Create a temporary table to store the values at the new grid points.
939 std::vector<std::vector<std::vector<double> > > tabClone;
940 Init(nE, nB, nA, tabClone, init);
941
942 // Fill the temporary table.
943 for (size_t i = 0; i < nE; ++i) {
944 const double e = efields[i];
945 for (size_t j = 0; j < nB; ++j) {
946 const double b = bfields[j];
947 for (size_t k = 0; k < nA; ++k) {
948 const double a = angles[k];
949 double val = 0.;
950 if (!Interpolate(e, b, a, tab, val, intp, extr)) {
951 std::cerr << m_className << "::Clone:\n"
952 << " Interpolation of " << lbl << " failed.\n"
953 << " Cannot copy value to new grid at E = " << e
954 << ", B = " << b << ", angle: " << a << "\n";
955 continue;
956 }
957 tabClone[k][j][i] = val;
958 }
959 }
960 }
961 // Copy the values to the original table.
962 tab.swap(tabClone);
963}
964
966 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
967 const size_t n, const std::vector<double>& efields,
968 const std::vector<double>& bfields, const std::vector<double>& angles,
969 const unsigned int intp, const std::pair<unsigned int, unsigned int>& extr,
970 const double init, const std::string& lbl) {
971 // If the table does not exist, do nothing.
972 if (tab.empty()) return;
973
974 // Get the dimensions of the new grid.
975 const auto nE = efields.size();
976 const auto nB = bfields.size();
977 const auto nA = angles.size();
978
979 // Create a temporary table to store the values at the new grid points.
980 std::vector<std::vector<std::vector<std::vector<double> > > > tabClone;
981 Init(nE, nB, nA, n, tabClone, init);
982
983 // Fill the temporary table.
984 for (size_t l = 0; l < n; ++l) {
985 for (size_t i = 0; i < nE; ++i) {
986 const double e = efields[i];
987 for (size_t j = 0; j < nB; ++j) {
988 const double b = bfields[j];
989 for (size_t k = 0; k < nA; ++k) {
990 const double a = angles[k];
991 double val = 0.;
992 if (!Interpolate(e, b, a, tab[l], val, intp, extr)) {
993 std::cerr << m_className << "::Clone:\n"
994 << " Interpolation of " << lbl << " failed.\n"
995 << " Cannot copy value to new grid at index " << l
996 << ", E = " << e << ", B = " << b << ", angle: " << a
997 << "\n";
998 continue;
999 }
1000 tabClone[l][k][j][i] = val;
1001 }
1002 }
1003 }
1004 }
1005 // Copy the values to the original table.
1006 tab.swap(tabClone);
1007}
1008
1009bool Medium::SetIonMobility(const size_t ie, const size_t ib,
1010 const size_t ia, const double mu) {
1011 // Check the index.
1012 if (ie >= m_eFields.size() || ib >= m_bFields.size() ||
1013 ia >= m_bAngles.size()) {
1014 PrintOutOfRange(m_className, "SetIonMobility", ie, ib, ia);
1015 return false;
1016 }
1017
1018 if (m_iMob.empty()) {
1019 std::cerr << m_className << "::SetIonMobility:\n"
1020 << " Ion mobility table not initialised.\n";
1021 return false;
1022 }
1023
1024 if (mu == 0.) {
1025 std::cerr << m_className << "::SetIonMobility: Zero value not allowed.\n";
1026 return false;
1027 }
1028
1029 m_iMob[ia][ib][ie] = mu;
1030 if (m_debug) {
1031 std::cout << m_className << "::SetIonMobility:\n Ion mobility at E = "
1032 << m_eFields[ie] << " V/cm, B = "
1033 << m_bFields[ib] << " T, angle "
1034 << m_bAngles[ia] << " set to " << mu << " cm2/(V ns).\n";
1035 }
1036 return true;
1037}
1038
1039bool Medium::SetIonMobility(const std::vector<double>& efields,
1040 const std::vector<double>& mobs) {
1041 if (efields.size() != mobs.size()) {
1042 std::cerr << m_className << "::SetIonMobility:\n"
1043 << " E-field and mobility arrays have different sizes.\n";
1044 return false;
1045 }
1046
1048 const auto nE = m_eFields.size();
1049 const auto nB = m_bFields.size();
1050 const auto nA = m_bAngles.size();
1051 Init(nE, nB, nA, m_iMob, 0.);
1052 for (size_t i = 0; i < nE; ++i) {
1053 const double e = m_eFields[i];
1054 const double mu = Interpolate1D(e, mobs, efields, m_intpMob, m_extrMob);
1055 m_iMob[0][0][i] = mu;
1056 }
1057
1058 if (m_tab2d) {
1059 for (size_t i = 0; i < nA; ++i) {
1060 for (size_t j = 0; j < nB; ++j) {
1061 for (size_t k = 0; k < nE; ++k) {
1062 m_iMob[i][j][k] = m_iMob[0][0][k];
1063 }
1064 }
1065 }
1066 }
1067 return true;
1068}
1069
1070void Medium::SetExtrapolationMethodVelocity(const std::string& low,
1071 const std::string& high) {
1072 SetExtrapolationMethod(low, high, m_extrVel, "Velocity");
1073}
1074
1076 const std::string& high) {
1077 SetExtrapolationMethod(low, high, m_extrDif, "Diffusion");
1078}
1079
1080void Medium::SetExtrapolationMethodTownsend(const std::string& low,
1081 const std::string& high) {
1082 SetExtrapolationMethod(low, high, m_extrAlp, "Townsend");
1083}
1084
1086 const std::string& high) {
1087 SetExtrapolationMethod(low, high, m_extrAtt, "Attachment");
1088}
1089
1091 const std::string& high) {
1092 SetExtrapolationMethod(low, high, m_extrMob, "IonMobility");
1093}
1094
1096 const std::string& high) {
1097 SetExtrapolationMethod(low, high, m_extrDis, "IonDissociation");
1098}
1099
1100void Medium::SetExtrapolationMethod(const std::string& low,
1101 const std::string& high,
1102 std::pair<unsigned int, unsigned int>& extr,
1103 const std::string& fcn) {
1104 unsigned int i = 0;
1105 if (GetExtrapolationIndex(low, i)) {
1106 extr.first = i;
1107 } else {
1108 std::cerr << m_className << "::SetExtrapolationMethod" << fcn << ":\n"
1109 << " Unknown extrapolation method (" << low << ")\n";
1110 }
1111 unsigned int j = 0;
1112 if (GetExtrapolationIndex(high, j)) {
1113 extr.second = j;
1114 } else {
1115 std::cerr << m_className << "::SetExtrapolationMethod" << fcn << ":\n"
1116 << " Unknown extrapolation method (" << high << ")\n";
1117 }
1118}
1119
1120bool Medium::GetExtrapolationIndex(std::string str, unsigned int& nb) const {
1121 // Convert to upper-case.
1122 std::transform(str.begin(), str.end(), str.begin(), toupper);
1123
1124 if (str == "CONST" || str == "CONSTANT") {
1125 nb = 0;
1126 } else if (str == "LIN" || str == "LINEAR") {
1127 nb = 1;
1128 } else if (str == "EXP" || str == "EXPONENTIAL") {
1129 nb = 2;
1130 } else {
1131 return false;
1132 }
1133
1134 return true;
1135}
1136
1138 const std::vector<std::vector<std::vector<double> > >& tab) const {
1139
1140 if (tab.empty()) return 0;
1141 const auto nE = m_eFields.size();
1142 const auto nB = m_bFields.size();
1143 const auto nA = m_bAngles.size();
1144 for (size_t i = 0; i < nE; ++i) {
1145 bool below = false;
1146 for (size_t k = 0; k < nA; ++k) {
1147 for (size_t j = 0; j < nB; ++j) {
1148 if (tab[k][j][i] < -20.) {
1149 below = true;
1150 break;
1151 }
1152 }
1153 if (below) break;
1154 }
1155 if (below) continue;
1156 return i;
1157 }
1158 return nE - 1;
1159}
1160
1161void Medium::SetInterpolationMethodVelocity(const unsigned int intrp) {
1162 if (intrp > 0) m_intpVel = intrp;
1163}
1164
1165void Medium::SetInterpolationMethodDiffusion(const unsigned int intrp) {
1166 if (intrp > 0) m_intpDif = intrp;
1167}
1168
1169void Medium::SetInterpolationMethodTownsend(const unsigned int intrp) {
1170 if (intrp > 0) m_intpAlp = intrp;
1171}
1172
1173void Medium::SetInterpolationMethodAttachment(const unsigned int intrp) {
1174 if (intrp > 0) m_intpAtt = intrp;
1175}
1176
1177void Medium::SetInterpolationMethodIonMobility(const unsigned int intrp) {
1178 if (intrp > 0) m_intpMob = intrp;
1179}
1180
1182 if (intrp > 0) m_intpDis = intrp;
1183}
1184
1185double Medium::GetAngle(const double ex, const double ey, const double ez,
1186 const double bx, const double by, const double bz,
1187 const double e, const double b) const {
1188 const double eb = e * b;
1189 if (eb <= 0.) return m_bAngles[0];
1190 const double einb = fabs(ex * bx + ey * by + ez * bz);
1191 if (einb > 0.2 * eb) {
1192 const double ebxy = ex * by - ey * bx;
1193 const double ebxz = ex * bz - ez * bx;
1194 const double ebzy = ez * by - ey * bz;
1195 return asin(
1196 std::min(1., sqrt(ebxy * ebxy + ebxz * ebxz + ebzy * ebzy) / eb));
1197 }
1198 return acos(std::min(1., einb / eb));
1199}
1200
1202 const double e, const double b, const double a,
1203 const std::vector<std::vector<std::vector<double> > >& table, double& y,
1204 const unsigned int intp,
1205 const std::pair<unsigned int, unsigned int>& extr) const {
1206 if (table.empty()) {
1207 y = 0.;
1208 return false; // TODO: true!
1209 }
1210
1211 if (m_tab2d) {
1213 m_bAngles.size(), m_bFields.size(), m_eFields.size(),
1214 a, b, e, y, intp);
1215 } else {
1216 y = Interpolate1D(e, table[0][0], m_eFields, intp, extr);
1217 }
1218 return true;
1219}
1220
1222 const double e, const std::vector<double>& table,
1223 const std::vector<double>& fields, const unsigned int intpMeth,
1224 const std::pair<unsigned int, unsigned int>& extr) const {
1225 // This function is a generalized version of the Fortran functions
1226 // GASVEL, GASVT1, GASVT2, GASLOR, GASMOB, GASDFT, and GASDFL
1227 // for the case of a 1D table. All variables are generic.
1228
1229 const auto nSizeTable = fields.size();
1230
1231 if (e < 0. || nSizeTable < 1) return 0.;
1232
1233 double result = 0.;
1234
1235 if (nSizeTable == 1) {
1236 // Only one point
1237 result = table[0];
1238 } else if (e < fields[0]) {
1239 // Extrapolation towards small fields
1240 if (fields[0] >= fields[1]) {
1241 if (m_debug) {
1242 std::cerr << m_className << "::Interpolate1D:\n";
1243 std::cerr << " First two field values coincide.\n";
1244 std::cerr << " No extrapolation to lower fields.\n";
1245 }
1246 result = table[0];
1247 } else if (extr.first == 1) {
1248 // Linear extrapolation
1249 const double extr4 = (table[1] - table[0]) / (fields[1] - fields[0]);
1250 const double extr3 = table[0] - extr4 * fields[0];
1251 result = extr3 + extr4 * e;
1252 } else if (extr.first == 2) {
1253 // Logarithmic extrapolation
1254 const double extr4 = log(table[1] / table[0]) / (fields[1] - fields[0]);
1255 const double extr3 = log(table[0] - extr4 * fields[0]);
1256 result = std::exp(std::min(50., extr3 + extr4 * e));
1257 } else {
1258 result = table[0];
1259 }
1260 } else if (e > fields[nSizeTable - 1]) {
1261 // Extrapolation towards large fields
1262 if (fields[nSizeTable - 1] <= fields[nSizeTable - 2]) {
1263 if (m_debug) {
1264 std::cerr << m_className << "::Interpolate1D:\n";
1265 std::cerr << " Last two field values coincide.\n";
1266 std::cerr << " No extrapolation to higher fields.\n";
1267 }
1268 result = table[nSizeTable - 1];
1269 } else if (extr.second == 1) {
1270 // Linear extrapolation
1271 const double extr2 = (table[nSizeTable - 1] - table[nSizeTable - 2]) /
1272 (fields[nSizeTable - 1] - fields[nSizeTable - 2]);
1273 const double extr1 =
1274 table[nSizeTable - 1] - extr2 * fields[nSizeTable - 1];
1275 result = extr1 + extr2 * e;
1276 } else if (extr.second == 2) {
1277 // Logarithmic extrapolation
1278 const double extr2 = log(table[nSizeTable - 1] / table[nSizeTable - 2]) /
1279 (fields[nSizeTable - 1] - fields[nSizeTable - 2]);
1280 const double extr1 =
1281 log(table[nSizeTable - 1]) - extr2 * fields[nSizeTable - 1];
1282 result = exp(std::min(50., extr1 + extr2 * e));
1283 } else {
1284 result = table[nSizeTable - 1];
1285 }
1286 } else {
1287 // Intermediate points, spline interpolation (not implemented).
1288 // Intermediate points, Newtonian interpolation
1289 result = Numerics::Divdif(table, fields, nSizeTable, e, intpMeth);
1290 }
1291
1292 return result;
1293}
1294
1295void Medium::Init(const size_t nE, const size_t nB, const size_t nA,
1296 std::vector<std::vector<std::vector<double> > >& tab,
1297 const double val) {
1298 if (nE == 0 || nB == 0 || nA == 0) {
1299 std::cerr << m_className << "::Init: Invalid grid.\n";
1300 return;
1301 }
1302 tab.assign(
1303 nA, std::vector<std::vector<double> >(nB, std::vector<double>(nE, val)));
1304}
1305
1307 const size_t nE, const size_t nB, const size_t nA, const size_t nT,
1308 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
1309 const double val) {
1310 if (nE == 0 || nB == 0 || nA == 0 || nT == 0) {
1311 std::cerr << m_className << "::Init: Invalid grid.\n";
1312 return;
1313 }
1314
1315 tab.assign(nT, std::vector<std::vector<std::vector<double> > >(
1316 nA, std::vector<std::vector<double> >(
1317 nB, std::vector<double>(nE, val))));
1318}
1319}
virtual double IonMobility()
Low-field mobility [cm2 V-1 ns-1].
Definition: Medium.cc:620
void ResetHoleAttachment()
Definition: Medium.hh:430
virtual void GetComponent(const unsigned int i, std::string &label, double &f)
Get the name and fraction of a given component.
Definition: Medium.cc:105
bool Interpolate(const double e, const double b, const double a, const std::vector< std::vector< std::vector< double > > > &table, double &y, const unsigned int intp, const std::pair< unsigned int, unsigned int > &extr) const
Definition: Medium.cc:1201
double m_density
Definition: Medium.hh:514
void ResetIonMobility()
Definition: Medium.hh:432
virtual bool GetPhotonCollision(const double e, int &type, int &level, double &e1, double &ctheta, int &nsec, double &esec)
Definition: Medium.cc:682
virtual double UnScaleElectricField(const double e) const
Definition: Medium.hh:464
void ResetElectronAttachment()
Definition: Medium.hh:416
virtual double GetMassDensity() const
Get the mass density [g/cm3].
Definition: Medium.cc:101
std::vector< double > m_bFields
Definition: Medium.hh:537
void ResetIonDiffusion()
Definition: Medium.hh:433
virtual void GetElectronMomentum(const double e, double &px, double &py, double &pz, int &band)
Definition: Medium.cc:473
void SetTemperature(const double t)
Set the temperature [K].
Definition: Medium.cc:71
void ResetElectronLorentzAngle()
Definition: Medium.hh:417
double m_pressure
Definition: Medium.hh:506
virtual bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:534
unsigned int m_intpMob
Definition: Medium.hh:591
virtual double ScaleAttachment(const double eta) const
Definition: Medium.hh:469
bool GetEntry(const size_t i, const size_t j, const size_t k, const std::string &fcn, const std::vector< std::vector< std::vector< double > > > &tab, double &val) const
Definition: Medium.cc:884
void SetExtrapolationMethodIonDissociation(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1095
void SetInterpolationMethodIonDissociation(const unsigned int intrp)
Definition: Medium.cc:1181
double Interpolate1D(const double e, const std::vector< double > &table, const std::vector< double > &fields, const unsigned int intpMeth, const std::pair< unsigned int, unsigned int > &extr) const
Definition: Medium.cc:1221
double GetAngle(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const double e, const double b) const
Definition: Medium.cc:1185
virtual double ScaleDiffusion(const double d) const
Definition: Medium.hh:466
virtual double ElectronMobility()
Low-field mobility [cm2 V-1 ns-1].
Definition: Medium.cc:453
virtual bool HoleVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:513
virtual double ScaleLorentzAngle(const double lor) const
Definition: Medium.hh:470
unsigned int m_intpDis
Definition: Medium.hh:592
virtual bool ElectronLorentzAngle(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
Lorentz angle.
Definition: Medium.cc:429
bool Diffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &difL, const std::vector< std::vector< std::vector< double > > > &difT, double &dl, double &dt) const
Definition: Medium.cc:269
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition: Medium.hh:546
void ResetHoleVelocity()
Definition: Medium.hh:419
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition: Medium.hh:558
void ResetIonDissociation()
Definition: Medium.hh:437
virtual bool GetDielectricFunction(const double e, double &eps1, double &eps2, const unsigned int i=0)
Get the complex dielectric function at a given energy.
Definition: Medium.cc:636
virtual void SetNumberDensity(const double n)
Set the number density [cm-3].
Definition: Medium.cc:134
std::string m_name
Definition: Medium.hh:502
unsigned int m_intpDif
Definition: Medium.hh:587
virtual bool ElectronVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:379
virtual bool IonVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:565
std::pair< unsigned int, unsigned int > m_extrAtt
Definition: Medium.hh:580
void SetExtrapolationMethodTownsend(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1080
virtual void SetAtomicNumber(const double z)
Set the effective atomic number.
Definition: Medium.cc:114
virtual bool GetElectronCollision(const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< std::pair< int, double > > &secondaries, int &ndxc, int &band)
Sample the collision type. Update energy and direction vector.
Definition: Medium.cc:491
std::vector< std::vector< std::vector< double > > > m_iDifT
Definition: Medium.hh:566
void SetExtrapolationMethodIonMobility(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1090
unsigned int m_hThrAtt
Definition: Medium.hh:573
std::pair< unsigned int, unsigned int > m_extrVel
Definition: Medium.hh:577
virtual double ScaleElectricField(const double e) const
Definition: Medium.hh:463
virtual double ScaleDiffusionTensor(const double d) const
Definition: Medium.hh:467
std::vector< std::vector< std::vector< double > > > m_hDifL
Definition: Medium.hh:556
virtual double GetElectronNullCollisionRate(const int band=0)
Null-collision rate [ns-1].
Definition: Medium.cc:480
size_t SetThreshold(const std::vector< std::vector< std::vector< double > > > &tab) const
Definition: Medium.cc:1137
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition: Medium.hh:541
std::vector< std::vector< std::vector< double > > > m_eVelX
Definition: Medium.hh:542
std::vector< std::vector< std::vector< double > > > m_eDifL
Definition: Medium.hh:544
void Init(const size_t nE, const size_t nB, const size_t nA, std::vector< std::vector< std::vector< double > > > &tab, const double val)
Definition: Medium.cc:1295
virtual double ScaleTownsend(const double alpha) const
Definition: Medium.hh:468
unsigned int m_intpVel
Definition: Medium.hh:586
virtual double GetPhotonCollisionRate(const double e)
Definition: Medium.cc:675
void SetDielectricConstant(const double eps)
Set the relative static dielectric constant.
Definition: Medium.cc:91
virtual bool GetDeexcitationProduct(const unsigned int i, double &t, double &s, int &type, double &energy) const
Definition: Medium.cc:504
virtual double GetElectronCollisionRate(const double e, const int band=0)
Collision rate [ns-1] for given electron energy.
Definition: Medium.cc:485
bool GetExtrapolationIndex(std::string str, unsigned int &nb) const
Definition: Medium.cc:1120
void SetInterpolationMethodVelocity(const unsigned int intrp)
Set the degree of polynomial interpolation (usually 2).
Definition: Medium.cc:1161
virtual bool ElectronDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Longitudinal and transverse diffusion coefficients [cm1/2].
Definition: Medium.cc:387
void SetInterpolationMethodDiffusion(const unsigned int intrp)
Definition: Medium.cc:1165
double m_epsilon
Definition: Medium.hh:508
void SetExtrapolationMethodDiffusion(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1075
std::vector< double > m_eFields
Definition: Medium.hh:536
void GetFieldGrid(std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles)
Get the fields and E-B angles used in the transport tables.
Definition: Medium.cc:860
static int m_idCounter
Definition: Medium.hh:495
std::vector< std::vector< std::vector< double > > > m_hDifT
Definition: Medium.hh:557
virtual bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:403
void ResetElectronDiffusion()
Definition: Medium.hh:410
void SetFieldGrid(double emin, double emax, const size_t ne, bool logE, double bmin=0., double bmax=0., const size_t nb=1, double amin=HalfPi, double amax=HalfPi, const size_t na=1)
Set the range of fields to be covered by the transport tables.
Definition: Medium.cc:693
virtual bool GetOpticalDataRange(double &emin, double &emax, const unsigned int i=0)
Get the energy range [eV] of the available optical data.
Definition: Medium.cc:624
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition: Medium.hh:547
virtual void SetMassDensity(const double rho)
Set the mass density [g/cm3].
Definition: Medium.cc:144
virtual double ScaleDissociation(const double diss) const
Definition: Medium.hh:471
bool Alpha(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &tab, unsigned int intp, const unsigned int thr, const std::pair< unsigned int, unsigned int > &extr, double &alpha) const
Definition: Medium.cc:348
std::vector< double > m_bAngles
Definition: Medium.hh:538
Medium()
Constructor.
Definition: Medium.cc:60
std::vector< std::vector< std::vector< double > > > m_iDifL
Definition: Medium.hh:565
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition: Medium.hh:553
std::vector< std::vector< std::vector< double > > > m_eLor
Definition: Medium.hh:548
std::vector< std::vector< std::vector< double > > > m_eDifT
Definition: Medium.hh:545
virtual double GetElectronEnergy(const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
Dispersion relation (energy vs. wave vector)
Definition: Medium.cc:458
void SetInterpolationMethodIonMobility(const unsigned int intrp)
Definition: Medium.cc:1177
unsigned int m_intpAtt
Definition: Medium.hh:589
void SetExtrapolationMethodVelocity(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1070
virtual bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:416
virtual double HoleMobility()
Low-field mobility [cm2 V-1 ns-1].
Definition: Medium.cc:560
std::pair< unsigned int, unsigned int > m_extrDis
Definition: Medium.hh:583
std::vector< std::vector< std::vector< double > > > m_hAtt
Definition: Medium.hh:559
void SetPressure(const double p)
Definition: Medium.cc:81
void SetExtrapolationMethod(const std::string &low, const std::string &high, std::pair< unsigned int, unsigned int > &extr, const std::string &fcn)
Definition: Medium.cc:1100
std::pair< unsigned int, unsigned int > m_extrAlp
Definition: Medium.hh:579
std::vector< std::vector< std::vector< double > > > m_eVelB
Definition: Medium.hh:543
unsigned int m_eThrAtt
Definition: Medium.hh:571
virtual void SetAtomicWeight(const double a)
Set the effective atomic weight.
Definition: Medium.cc:124
void SetExtrapolationMethodAttachment(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1085
std::pair< unsigned int, unsigned int > m_extrLor
Definition: Medium.hh:581
unsigned int m_nComponents
Definition: Medium.hh:500
unsigned int m_intpLor
Definition: Medium.hh:590
std::string m_className
Definition: Medium.hh:493
std::pair< unsigned int, unsigned int > m_extrDif
Definition: Medium.hh:578
void SetInterpolationMethodAttachment(const unsigned int intrp)
Definition: Medium.cc:1173
virtual bool HoleDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Longitudinal and transverse diffusion coefficients [cm1/2].
Definition: Medium.cc:521
unsigned int m_iThrDis
Definition: Medium.hh:574
virtual bool IonDissociation(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &diss)
Dissociation coefficient.
Definition: Medium.cc:607
virtual void ResetTables()
Reset all tables of transport parameters.
Definition: Medium.cc:901
std::vector< std::vector< std::vector< double > > > m_hVelX
Definition: Medium.hh:554
virtual ~Medium()
Destructor.
Definition: Medium.cc:69
std::pair< unsigned int, unsigned int > m_extrMob
Definition: Medium.hh:582
void Clone(std::vector< std::vector< std::vector< double > > > &tab, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int intp, const std::pair< unsigned int, unsigned int > &extr, const double init, const std::string &label)
Definition: Medium.cc:918
void ResetHoleDiffusion()
Definition: Medium.hh:424
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
Definition: Medium.hh:550
std::vector< std::vector< std::vector< double > > > m_hVelB
Definition: Medium.hh:555
bool m_isChanged
Definition: Medium.hh:527
bool Velocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &velE, const std::vector< std::vector< std::vector< double > > > &velB, const std::vector< std::vector< std::vector< double > > > &velX, const double q, double &vx, double &vy, double &vz) const
Definition: Medium.cc:160
std::vector< std::vector< std::vector< std::vector< double > > > > m_hDifM
Definition: Medium.hh:561
unsigned int m_hThrAlp
Definition: Medium.hh:572
std::vector< std::vector< std::vector< double > > > m_iMob
Definition: Medium.hh:564
bool SetEntry(const size_t i, const size_t j, const size_t k, const std::string &fcn, std::vector< std::vector< std::vector< double > > > &tab, const double val)
Definition: Medium.cc:868
void ResetElectronTownsend()
Definition: Medium.hh:415
unsigned int m_eThrAlp
Definition: Medium.hh:570
virtual bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:547
void ResetHoleTownsend()
Definition: Medium.hh:429
void ResetElectronVelocity()
Definition: Medium.hh:405
void SetInterpolationMethodTownsend(const unsigned int intrp)
Definition: Medium.cc:1169
virtual bool GetPhotoAbsorptionCrossSection(const double e, double &sigma, const unsigned int i=0)
Definition: Medium.cc:654
virtual bool IonDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Longitudinal and transverse diffusion coefficients [cm1/2].
Definition: Medium.cc:600
double m_temperature
Definition: Medium.hh:504
std::vector< std::vector< std::vector< double > > > m_iDis
Definition: Medium.hh:567
bool SetIonMobility(const std::vector< double > &fields, const std::vector< double > &mobilities)
Definition: Medium.cc:1039
unsigned int m_intpAlp
Definition: Medium.hh:588
bool Boxin3(const std::vector< std::vector< std::vector< double > > > &value, const std::vector< double > &xAxis, const std::vector< double > &yAxis, const std::vector< double > &zAxis, const int nx, const int ny, const int nz, const double xx, const double yy, const double zz, double &f, const int iOrder)
Definition: Numerics.cc:1496
double Divdif(const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
Definition: Numerics.cc:1206
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:107