Garfield++ 3.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 <cmath>
2#include <iomanip>
3#include <iostream>
4#include <utility>
5
8#include "Garfield/Medium.hh"
10#include "Garfield/Random.hh"
11
12namespace {
13
14void PrintNotImplemented(const std::string& cls, const std::string& fcn) {
15 std::cerr << cls << "::" << fcn << ": Function is not implemented.\n";
16}
17
18void PrintOutOfRange(const std::string& cls, const std::string& fcn,
19 const unsigned int i, const unsigned int j,
20 const unsigned int 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> eFieldsNew(ne);
770 std::vector<double> bFieldsNew(nb);
771 std::vector<double> bAnglesNew(na);
772 for (size_t i = 0; i < ne; ++i) {
773 eFieldsNew[i] = logE ? emin * pow(estep, i) : emin + i * estep;
774 }
775 for (size_t i = 0; i < nb; ++i) {
776 bFieldsNew[i] = bmin + i * bstep;
777 }
778 if (na == 1 && nb == 1 && fabs(bmin) < 1.e-4) {
779 bAnglesNew[0] = HalfPi;
780 } else {
781 for (size_t i = 0; i < na; ++i) {
782 bAnglesNew[i] = amin + i * astep;
783 }
784 }
785 SetFieldGrid(eFieldsNew, bFieldsNew, bAnglesNew);
786}
787
788void Medium::SetFieldGrid(const std::vector<double>& efields,
789 const std::vector<double>& bfields,
790 const std::vector<double>& angles) {
791 const std::string hdr = m_className + "::SetFieldGrid";
792 if (!CheckFields(efields, hdr, "E-fields")) return;
793 if (!CheckFields(bfields, hdr, "B-fields")) return;
794 if (!CheckFields(angles, hdr, "angles")) return;
795
796 if (m_debug) {
797 std::cout << m_className << "::SetFieldGrid:\n E-fields:\n";
798 for (const auto efield : efields) std::cout << " " << efield << "\n";
799 std::cout << " B-fields:\n";
800 for (const auto bfield : bfields) std::cout << " " << bfield << "\n";
801 std::cout << " Angles:\n";
802 for (const auto angle : angles) std::cout << " " << angle << "\n";
803 }
804
805 // Clone the existing tables.
806 // Electrons
807 Clone(m_eVelE, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
808 "electron velocity along E");
809 Clone(m_eVelB, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
810 "electron velocity along Bt");
811 Clone(m_eVelX, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
812 "electron velocity along ExB");
813 Clone(m_eDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
814 "electron longitudinal diffusion");
815 Clone(m_eDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
816 "electron transverse diffusion");
817 Clone(m_eAlp, efields, bfields, angles, m_intpAlp, m_extrAlp, -30.,
818 "electron Townsend coefficient");
819 Clone(m_eAtt, efields, bfields, angles, m_intpAtt, m_extrAtt, -30.,
820 "electron attachment coefficient");
821 Clone(m_eLor, efields, bfields, angles, m_intpLor, m_extrLor, 0.,
822 "electron Lorentz angle");
823 if (!m_eDifM.empty()) {
824 Clone(m_eDifM, 6, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
825 "electron diffusion tensor");
826 }
827
828 // Holes
829 Clone(m_hVelE, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
830 "hole velocity along E");
831 Clone(m_hVelB, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
832 "hole velocity along Bt");
833 Clone(m_hVelX, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
834 "hole velocity along ExB");
835 Clone(m_hDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
836 "hole longitudinal diffusion");
837 Clone(m_hDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
838 "hole transverse diffusion");
839 Clone(m_hAlp, efields, bfields, angles, m_intpAlp, m_extrAlp, -30.,
840 "hole Townsend coefficient");
841 Clone(m_hAtt, efields, bfields, angles, m_intpAtt, m_extrAtt, -30.,
842 "hole attachment coefficient");
843 if (!m_hDifM.empty()) {
844 Clone(m_hDifM, 6, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
845 "hole diffusion tensor");
846 }
847
848 // Ions
849 Clone(m_iMob, efields, bfields, angles, m_intpMob, m_extrMob, 0.,
850 "ion mobility");
851 Clone(m_iDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
852 "ion longitudinal diffusion");
853 Clone(m_iDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
854 "ion transverse diffusion");
855 Clone(m_iDis, efields, bfields, angles, m_intpDis, m_extrDis, -30.,
856 "ion dissociation");
857
858 if (bfields.size() > 1 || angles.size() > 1) m_tab2d = true;
859 m_eFields = efields;
860 m_bFields = bfields;
861 m_bAngles = angles;
862}
863
864void Medium::GetFieldGrid(std::vector<double>& efields,
865 std::vector<double>& bfields,
866 std::vector<double>& angles) {
867 efields = m_eFields;
868 bfields = m_bFields;
869 angles = m_bAngles;
870}
871
872bool Medium::SetEntry(const unsigned int i, const unsigned int j,
873 const unsigned int k, const std::string& fcn,
874 std::vector<std::vector<std::vector<double> > >& tab,
875 const double val) {
876
877 if (i >= m_eFields.size() || j >= m_bFields.size() || k >= m_bAngles.size()) {
878 PrintOutOfRange(m_className, fcn, i, j, k);
879 return false;
880 }
881 if (tab.empty()) {
882 Init(m_eFields.size(), m_bFields.size(), m_bAngles.size(), tab, val);
883 }
884 tab[k][j][i] = val;
885 return true;
886}
887
888bool Medium::GetEntry(const unsigned int i, const unsigned int j,
889 const unsigned int k, const std::string& fcn,
890 const std::vector<std::vector<std::vector<double> > >& tab,
891 double& val) const {
892 val = 0.;
893 if (i >= m_eFields.size() || j >= m_bFields.size() || k >= m_bAngles.size()) {
894 PrintOutOfRange(m_className, fcn, i, j, k);
895 return false;
896 }
897 if (tab.empty()) {
898 if (m_debug) PrintDataNotAvailable(m_className, fcn);
899 return false;
900 }
901 val = tab[k][j][i];
902 return true;
903}
904
911
916
920}
921
922void Medium::Clone(std::vector<std::vector<std::vector<double> > >& tab,
923 const std::vector<double>& efields,
924 const std::vector<double>& bfields,
925 const std::vector<double>& angles,
926 const unsigned int intp,
927 const std::pair<unsigned int, unsigned int>& extr,
928 const double init, const std::string& lbl) {
929 if (m_debug) {
930 std::cout << m_className << "::Clone: Copying " << lbl << " to new grid.\n";
931 }
932
933 if (tab.empty()) {
934 if (m_debug) std::cout << m_className << "::Clone: Table is empty.\n";
935 return;
936 }
937 // Get the dimensions of the new grid.
938 const auto nE = efields.size();
939 const auto nB = bfields.size();
940 const auto nA = angles.size();
941
942 // Create a temporary table to store the values at the new grid points.
943 std::vector<std::vector<std::vector<double> > > tabClone;
944 Init(nE, nB, nA, tabClone, init);
945
946 // Fill the temporary table.
947 for (size_t i = 0; i < nE; ++i) {
948 const double e = efields[i];
949 for (size_t j = 0; j < nB; ++j) {
950 const double b = bfields[j];
951 for (size_t k = 0; k < nA; ++k) {
952 const double a = angles[k];
953 double val = 0.;
954 if (!Interpolate(e, b, a, tab, val, intp, extr)) {
955 std::cerr << m_className << "::Clone:\n"
956 << " Interpolation of " << lbl << " failed.\n"
957 << " Cannot copy value to new grid at E = " << e
958 << ", B = " << b << ", angle: " << a << "\n";
959 continue;
960 }
961 tabClone[k][j][i] = val;
962 }
963 }
964 }
965 // Copy the values to the original table.
966 tab.swap(tabClone);
967 tabClone.clear();
968}
969
971 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
972 const unsigned int n, const std::vector<double>& efields,
973 const std::vector<double>& bfields, const std::vector<double>& angles,
974 const unsigned int intp, const std::pair<unsigned int, unsigned int>& extr,
975 const double init, const std::string& lbl) {
976 // If the table does not exist, do nothing.
977 if (tab.empty()) return;
978
979 // Get the dimensions of the new grid.
980 const unsigned int nE = efields.size();
981 const unsigned int nB = bfields.size();
982 const unsigned int nA = angles.size();
983
984 // Create a temporary table to store the values at the new grid points.
985 std::vector<std::vector<std::vector<std::vector<double> > > > tabClone;
986 Init(nE, nB, nA, n, tabClone, init);
987
988 // Fill the temporary table.
989 for (unsigned int l = 0; l < n; ++l) {
990 for (unsigned int i = 0; i < nE; ++i) {
991 const double e = efields[i];
992 for (unsigned int j = 0; j < nB; ++j) {
993 const double b = bfields[j];
994 for (unsigned int k = 0; k < nA; ++k) {
995 const double a = angles[k];
996 double val = 0.;
997 if (!Interpolate(e, b, a, tab[l], val, intp, extr)) {
998 std::cerr << m_className << "::Clone:\n"
999 << " Interpolation of " << lbl << " failed.\n"
1000 << " Cannot copy value to new grid at index " << l
1001 << ", E = " << e << ", B = " << b << ", angle: " << a
1002 << "\n";
1003 continue;
1004 }
1005 tabClone[l][k][j][i] = val;
1006 }
1007 }
1008 }
1009 }
1010 // Copy the values to the original table.
1011 tab.swap(tabClone);
1012}
1013
1014bool Medium::SetIonMobility(const unsigned int ie, const unsigned int ib,
1015 const unsigned int ia, const double mu) {
1016 // Check the index.
1017 if (ie >= m_eFields.size() || ib >= m_bFields.size() ||
1018 ia >= m_bAngles.size()) {
1019 PrintOutOfRange(m_className, "SetIonMobility", ie, ib, ia);
1020 return false;
1021 }
1022
1023 if (m_iMob.empty()) {
1024 std::cerr << m_className << "::SetIonMobility:\n"
1025 << " Ion mobility table not initialised.\n";
1026 return false;
1027 }
1028
1029 if (mu == 0.) {
1030 std::cerr << m_className << "::SetIonMobility: Zero value not allowed.\n";
1031 return false;
1032 }
1033
1034 m_iMob[ia][ib][ie] = mu;
1035 if (m_debug) {
1036 std::cout << m_className << "::SetIonMobility:\n Ion mobility at E = "
1037 << m_eFields[ie] << " V/cm, B = "
1038 << m_bFields[ib] << " T, angle "
1039 << m_bAngles[ia] << " set to " << mu << " cm2/(V ns).\n";
1040 }
1041 return true;
1042}
1043
1044bool Medium::SetIonMobility(const std::vector<double>& efields,
1045 const std::vector<double>& mobs) {
1046 const int ne = efields.size();
1047 const int nm = mobs.size();
1048 if (ne != nm) {
1049 std::cerr << m_className << "::SetIonMobility:\n"
1050 << " E-field and mobility arrays have different sizes.\n";
1051 return false;
1052 }
1053
1055 const unsigned int nEfields = m_eFields.size();
1056 const unsigned int nBfields = m_bFields.size();
1057 const unsigned int nAngles = m_bAngles.size();
1058 Init(nEfields, nBfields, nAngles, m_iMob, 0.);
1059 for (unsigned int i = 0; i < nEfields; ++i) {
1060 const double e = m_eFields[i];
1061 const double mu = Interpolate1D(e, mobs, efields, m_intpMob, m_extrMob);
1062 m_iMob[0][0][i] = mu;
1063 }
1064
1065 if (m_tab2d) {
1066 for (unsigned int i = 0; i < nAngles; ++i) {
1067 for (unsigned int j = 0; j < nBfields; ++j) {
1068 for (unsigned int k = 0; k < nEfields; ++k) {
1069 m_iMob[i][j][k] = m_iMob[0][0][k];
1070 }
1071 }
1072 }
1073 }
1074 return true;
1075}
1076
1077void Medium::SetExtrapolationMethodVelocity(const std::string& low,
1078 const std::string& high) {
1079 SetExtrapolationMethod(low, high, m_extrVel, "Velocity");
1080}
1081
1083 const std::string& high) {
1084 SetExtrapolationMethod(low, high, m_extrDif, "Diffusion");
1085}
1086
1087void Medium::SetExtrapolationMethodTownsend(const std::string& low,
1088 const std::string& high) {
1089 SetExtrapolationMethod(low, high, m_extrAlp, "Townsend");
1090}
1091
1093 const std::string& high) {
1094 SetExtrapolationMethod(low, high, m_extrAtt, "Attachment");
1095}
1096
1098 const std::string& high) {
1099 SetExtrapolationMethod(low, high, m_extrMob, "IonMobility");
1100}
1101
1103 const std::string& high) {
1104 SetExtrapolationMethod(low, high, m_extrDis, "IonDissociation");
1105}
1106
1107void Medium::SetExtrapolationMethod(const std::string& low,
1108 const std::string& high,
1109 std::pair<unsigned int, unsigned int>& extr,
1110 const std::string& fcn) {
1111 unsigned int i = 0;
1112 if (GetExtrapolationIndex(low, i)) {
1113 extr.first = i;
1114 } else {
1115 std::cerr << m_className << "::SetExtrapolationMethod" << fcn << ":\n"
1116 << " Unknown extrapolation method (" << low << ")\n";
1117 }
1118 unsigned int j = 0;
1119 if (GetExtrapolationIndex(high, j)) {
1120 extr.second = j;
1121 } else {
1122 std::cerr << m_className << "::SetExtrapolationMethod" << fcn << ":\n"
1123 << " Unknown extrapolation method (" << high << ")\n";
1124 }
1125}
1126
1127bool Medium::GetExtrapolationIndex(std::string str, unsigned int& nb) const {
1128 // Convert to upper-case
1129 for (unsigned int i = 0; i < str.length(); ++i) {
1130 str[i] = toupper(str[i]);
1131 }
1132
1133 if (str == "CONST" || str == "CONSTANT") {
1134 nb = 0;
1135 } else if (str == "LIN" || str == "LINEAR") {
1136 nb = 1;
1137 } else if (str == "EXP" || str == "EXPONENTIAL") {
1138 nb = 2;
1139 } else {
1140 return false;
1141 }
1142
1143 return true;
1144}
1145
1147 const std::vector<std::vector<std::vector<double> > >& tab) const {
1148
1149 if (tab.empty()) return 0;
1150 const unsigned int nE = m_eFields.size();
1151 const unsigned int nB = m_bFields.size();
1152 const unsigned int nA = m_bAngles.size();
1153 for (unsigned int i = 0; i < nE; ++i) {
1154 bool below = false;
1155 for (unsigned int k = 0; k < nA; ++k) {
1156 for (unsigned int j = 0; j < nB; ++j) {
1157 if (tab[k][j][i] < -20.) {
1158 below = true;
1159 break;
1160 }
1161 }
1162 if (below) break;
1163 }
1164 if (below) continue;
1165 return i;
1166 }
1167 return nE - 1;
1168}
1169
1170void Medium::SetInterpolationMethodVelocity(const unsigned int intrp) {
1171 if (intrp > 0) m_intpVel = intrp;
1172}
1173
1174void Medium::SetInterpolationMethodDiffusion(const unsigned int intrp) {
1175 if (intrp > 0) m_intpDif = intrp;
1176}
1177
1178void Medium::SetInterpolationMethodTownsend(const unsigned int intrp) {
1179 if (intrp > 0) m_intpAlp = intrp;
1180}
1181
1182void Medium::SetInterpolationMethodAttachment(const unsigned int intrp) {
1183 if (intrp > 0) m_intpAtt = intrp;
1184}
1185
1186void Medium::SetInterpolationMethodIonMobility(const unsigned int intrp) {
1187 if (intrp > 0) m_intpMob = intrp;
1188}
1189
1191 if (intrp > 0) m_intpDis = intrp;
1192}
1193
1194double Medium::GetAngle(const double ex, const double ey, const double ez,
1195 const double bx, const double by, const double bz,
1196 const double e, const double b) const {
1197 const double eb = e * b;
1198 if (eb <= 0.) return m_bAngles[0];
1199 const double einb = fabs(ex * bx + ey * by + ez * bz);
1200 if (einb > 0.2 * eb) {
1201 const double ebxy = ex * by - ey * bx;
1202 const double ebxz = ex * bz - ez * bx;
1203 const double ebzy = ez * by - ey * bz;
1204 return asin(
1205 std::min(1., sqrt(ebxy * ebxy + ebxz * ebxz + ebzy * ebzy) / eb));
1206 }
1207 return acos(std::min(1., einb / eb));
1208}
1209
1211 const double e, const double b, const double a,
1212 const std::vector<std::vector<std::vector<double> > >& table, double& y,
1213 const unsigned int intp,
1214 const std::pair<unsigned int, unsigned int>& extr) const {
1215 if (table.empty()) {
1216 y = 0.;
1217 return false; // TODO: true!
1218 }
1219
1220 if (m_tab2d) {
1222 m_bAngles.size(), m_bFields.size(), m_eFields.size(),
1223 a, b, e, y, intp);
1224 } else {
1225 y = Interpolate1D(e, table[0][0], m_eFields, intp, extr);
1226 }
1227 return true;
1228}
1229
1231 const double e, const std::vector<double>& table,
1232 const std::vector<double>& fields, const unsigned int intpMeth,
1233 const std::pair<unsigned int, unsigned int>& extr) const {
1234 // This function is a generalized version of the Fortran functions
1235 // GASVEL, GASVT1, GASVT2, GASLOR, GASMOB, GASDFT, and GASDFL
1236 // for the case of a 1D table. All variables are generic.
1237
1238 const int nSizeTable = fields.size();
1239
1240 if (e < 0. || nSizeTable < 1) return 0.;
1241
1242 double result = 0.;
1243
1244 if (nSizeTable == 1) {
1245 // Only one point
1246 result = table[0];
1247 } else if (e < fields[0]) {
1248 // Extrapolation towards small fields
1249 if (fields[0] >= fields[1]) {
1250 if (m_debug) {
1251 std::cerr << m_className << "::Interpolate1D:\n";
1252 std::cerr << " First two field values coincide.\n";
1253 std::cerr << " No extrapolation to lower fields.\n";
1254 }
1255 result = table[0];
1256 } else if (extr.first == 1) {
1257 // Linear extrapolation
1258 const double extr4 = (table[1] - table[0]) / (fields[1] - fields[0]);
1259 const double extr3 = table[0] - extr4 * fields[0];
1260 result = extr3 + extr4 * e;
1261 } else if (extr.first == 2) {
1262 // Logarithmic extrapolation
1263 const double extr4 = log(table[1] / table[0]) / (fields[1] - fields[0]);
1264 const double extr3 = log(table[0] - extr4 * fields[0]);
1265 result = std::exp(std::min(50., extr3 + extr4 * e));
1266 } else {
1267 result = table[0];
1268 }
1269 } else if (e > fields[nSizeTable - 1]) {
1270 // Extrapolation towards large fields
1271 if (fields[nSizeTable - 1] <= fields[nSizeTable - 2]) {
1272 if (m_debug) {
1273 std::cerr << m_className << "::Interpolate1D:\n";
1274 std::cerr << " Last two field values coincide.\n";
1275 std::cerr << " No extrapolation to higher fields.\n";
1276 }
1277 result = table[nSizeTable - 1];
1278 } else if (extr.second == 1) {
1279 // Linear extrapolation
1280 const double extr2 = (table[nSizeTable - 1] - table[nSizeTable - 2]) /
1281 (fields[nSizeTable - 1] - fields[nSizeTable - 2]);
1282 const double extr1 =
1283 table[nSizeTable - 1] - extr2 * fields[nSizeTable - 1];
1284 result = extr1 + extr2 * e;
1285 } else if (extr.second == 2) {
1286 // Logarithmic extrapolation
1287 const double extr2 = log(table[nSizeTable - 1] / table[nSizeTable - 2]) /
1288 (fields[nSizeTable - 1] - fields[nSizeTable - 2]);
1289 const double extr1 =
1290 log(table[nSizeTable - 1]) - extr2 * fields[nSizeTable - 1];
1291 result = exp(std::min(50., extr1 + extr2 * e));
1292 } else {
1293 result = table[nSizeTable - 1];
1294 }
1295 } else {
1296 // Intermediate points, spline interpolation (not implemented).
1297 // Intermediate points, Newtonian interpolation
1298 result = Numerics::Divdif(table, fields, nSizeTable, e, intpMeth);
1299 }
1300
1301 return result;
1302}
1303
1304void Medium::Init(const size_t nE, const size_t nB, const size_t nA,
1305 std::vector<std::vector<std::vector<double> > >& tab,
1306 const double val) {
1307 if (nE == 0 || nB == 0 || nA == 0) {
1308 std::cerr << m_className << "::Init: Invalid grid.\n";
1309 return;
1310 }
1311 tab.assign(
1312 nA, std::vector<std::vector<double> >(nB, std::vector<double>(nE, val)));
1313}
1314
1316 const size_t nE, const size_t nB, const size_t nA, const size_t nT,
1317 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
1318 const double val) {
1319 if (nE == 0 || nB == 0 || nA == 0 || nT == 0) {
1320 std::cerr << m_className << "::Init: Invalid grid.\n";
1321 return;
1322 }
1323
1324 tab.assign(nT, std::vector<std::vector<std::vector<double> > >(
1325 nA, std::vector<std::vector<double> >(
1326 nB, std::vector<double>(nE, val))));
1327}
1328}
virtual double IonMobility()
Low-field mobility [cm2 V-1 ns-1].
Definition: Medium.cc:620
void ResetHoleAttachment()
Definition: Medium.hh:443
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:1210
double m_density
Definition: Medium.hh:527
void ResetIonMobility()
Definition: Medium.hh:445
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:477
void ResetElectronAttachment()
Definition: Medium.hh:429
virtual double GetMassDensity() const
Get the mass density [g/cm3].
Definition: Medium.cc:101
std::vector< double > m_bFields
Definition: Medium.hh:547
void ResetIonDiffusion()
Definition: Medium.hh:446
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:430
double m_pressure
Definition: Medium.hh:517
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:603
virtual double ScaleAttachment(const double eta) const
Definition: Medium.hh:482
void SetExtrapolationMethodIonDissociation(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1102
void SetInterpolationMethodIonDissociation(const unsigned int intrp)
Definition: Medium.cc:1190
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:1230
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:1194
virtual double ScaleDiffusion(const double d) const
Definition: Medium.hh:479
virtual double ElectronMobility()
Low-field mobility [cm2 V-1 ns-1].
Definition: Medium.cc:453
unsigned int SetThreshold(const std::vector< std::vector< std::vector< double > > > &tab) const
Definition: Medium.cc:1146
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:483
unsigned int m_intpDis
Definition: Medium.hh:604
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:558
void ResetHoleVelocity()
Definition: Medium.hh:432
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition: Medium.hh:570
void ResetIonDissociation()
Definition: Medium.hh:450
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:513
unsigned int m_intpDif
Definition: Medium.hh:599
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:592
void SetExtrapolationMethodTownsend(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1087
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:578
bool GetEntry(const unsigned int i, const unsigned int j, const unsigned int k, const std::string &fcn, const std::vector< std::vector< std::vector< double > > > &tab, double &val) const
Definition: Medium.cc:888
void SetExtrapolationMethodIonMobility(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1097
unsigned int m_hThrAtt
Definition: Medium.hh:585
std::pair< unsigned int, unsigned int > m_extrVel
Definition: Medium.hh:589
virtual double ScaleElectricField(const double e) const
Definition: Medium.hh:476
bool SetEntry(const unsigned int i, const unsigned int j, const unsigned int k, const std::string &fcn, std::vector< std::vector< std::vector< double > > > &tab, const double val)
Definition: Medium.cc:872
virtual double ScaleDiffusionTensor(const double d) const
Definition: Medium.hh:480
std::vector< std::vector< std::vector< double > > > m_hDifL
Definition: Medium.hh:568
virtual double GetElectronNullCollisionRate(const int band=0)
Null-collision rate [ns-1].
Definition: Medium.cc:480
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition: Medium.hh:553
std::vector< std::vector< std::vector< double > > > m_eVelX
Definition: Medium.hh:554
std::vector< std::vector< std::vector< double > > > m_eDifL
Definition: Medium.hh:556
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:1304
virtual double ScaleTownsend(const double alpha) const
Definition: Medium.hh:481
unsigned int m_intpVel
Definition: Medium.hh:598
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:1127
void SetInterpolationMethodVelocity(const unsigned int intrp)
Set the degree of polynomial interpolation (usually 2).
Definition: Medium.cc:1170
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:1174
double m_epsilon
Definition: Medium.hh:519
void SetExtrapolationMethodDiffusion(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1082
std::vector< double > m_eFields
Definition: Medium.hh:546
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:864
static int m_idCounter
Definition: Medium.hh:508
std::vector< std::vector< std::vector< double > > > m_hDifT
Definition: Medium.hh:569
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:423
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:559
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:484
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:548
Medium()
Constructor.
Definition: Medium.cc:60
std::vector< std::vector< std::vector< double > > > m_iDifL
Definition: Medium.hh:577
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition: Medium.hh:565
std::vector< std::vector< std::vector< double > > > m_eLor
Definition: Medium.hh:560
std::vector< std::vector< std::vector< double > > > m_eDifT
Definition: Medium.hh:557
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:1186
unsigned int m_intpAtt
Definition: Medium.hh:601
void SetExtrapolationMethodVelocity(const std::string &extrLow, const std::string &extrHigh)
Definition: Medium.cc:1077
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:595
std::vector< std::vector< std::vector< double > > > m_hAtt
Definition: Medium.hh:571
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:1107
std::pair< unsigned int, unsigned int > m_extrAlp
Definition: Medium.hh:591
std::vector< std::vector< std::vector< double > > > m_eVelB
Definition: Medium.hh:555
unsigned int m_eThrAtt
Definition: Medium.hh:583
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:1092
std::pair< unsigned int, unsigned int > m_extrLor
Definition: Medium.hh:593
unsigned int m_nComponents
Definition: Medium.hh:521
unsigned int m_intpLor
Definition: Medium.hh:602
std::string m_className
Definition: Medium.hh:506
std::pair< unsigned int, unsigned int > m_extrDif
Definition: Medium.hh:590
void SetInterpolationMethodAttachment(const unsigned int intrp)
Definition: Medium.cc:1182
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:586
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:905
std::vector< std::vector< std::vector< double > > > m_hVelX
Definition: Medium.hh:566
virtual ~Medium()
Destructor.
Definition: Medium.cc:69
std::pair< unsigned int, unsigned int > m_extrMob
Definition: Medium.hh:594
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:922
void ResetHoleDiffusion()
Definition: Medium.hh:437
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
Definition: Medium.hh:562
std::vector< std::vector< std::vector< double > > > m_hVelB
Definition: Medium.hh:567
bool m_isChanged
Definition: Medium.hh:540
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:573
unsigned int m_hThrAlp
Definition: Medium.hh:584
std::vector< std::vector< std::vector< double > > > m_iMob
Definition: Medium.hh:576
void ResetElectronTownsend()
Definition: Medium.hh:428
unsigned int m_eThrAlp
Definition: Medium.hh:582
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:442
void ResetElectronVelocity()
Definition: Medium.hh:418
void SetInterpolationMethodTownsend(const unsigned int intrp)
Definition: Medium.cc:1178
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:515
std::vector< std::vector< std::vector< double > > > m_iDis
Definition: Medium.hh:579
bool SetIonMobility(const std::vector< double > &fields, const std::vector< double > &mobilities)
Definition: Medium.cc:1044
unsigned int m_intpAlp
Definition: Medium.hh:600
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:1495
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