Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Medium.hh
Go to the documentation of this file.
1#ifndef G_MEDIUM_H
2#define G_MEDIUM_H
3
4#include <string>
5#include <vector>
6
9
10class TPad;
11
12namespace Garfield {
13
14class MediumGPU;
15
17class Medium {
18 public:
22 virtual ~Medium();
23
25 int GetId() const { return m_id; }
26
28 const std::string& GetName() const { return m_name; }
30 virtual bool IsGas() const { return false; }
32 virtual bool IsSemiconductor() const { return false; }
34 virtual bool IsConductor() const { return false; }
35
37 void SetTemperature(const double t);
39 double GetTemperature() const { return m_temperature; }
40 // Set the pressure [Torr].
41 void SetPressure(const double p);
42 // Get the pressure [Torr].
43 double GetPressure() const { return m_pressure; }
45 void SetDielectricConstant(const double eps);
47 double GetDielectricConstant() const { return m_epsilon; }
48
50 unsigned int GetNumberOfComponents() const { return m_nComponents; }
52 virtual void GetComponent(const unsigned int i, std::string& label,
53 double& f);
55 virtual void SetAtomicNumber(const double z);
57 virtual double GetAtomicNumber() const { return m_z; }
59 virtual void SetAtomicWeight(const double a);
61 virtual double GetAtomicWeight() const { return m_a; }
63 virtual void SetNumberDensity(const double n);
65 virtual double GetNumberDensity() const { return m_density; }
67 virtual void SetMassDensity(const double rho);
69 virtual double GetMassDensity() const;
70
72 virtual void EnableDrift(const bool on = true) { m_driftable = on; }
74 virtual void EnablePrimaryIonisation(const bool on = true) {
75 m_ionisable = on;
76 }
77
79 bool IsDriftable() const { return m_driftable; }
81 bool IsMicroscopic() const { return m_microscopic; }
82
84 bool IsIonisable() const { return m_ionisable; }
85
87 void SetW(const double w) { m_w = w; }
89 double GetW() const { return m_w; }
91 void SetFanoFactor(const double f) { m_fano = f; }
93 double GetFanoFactor() const { return m_fano; }
94
96 void PlotVelocity(const std::string& carriers, TPad* pad);
98 void PlotDiffusion(const std::string& carriers, TPad* pad);
100 void PlotTownsend(const std::string& carriers, TPad* pad);
102 void PlotAttachment(const std::string& carriers, TPad* pad);
104 void PlotAlphaEta(const std::string& carriers, TPad* pad);
105
106 // Transport parameters for electrons
108 virtual bool ElectronVelocity(const double ex, const double ey,
109 const double ez, const double bx,
110 const double by, const double bz, double& vx,
111 double& vy, double& vz);
114 virtual bool ElectronVelocityFluxBulk(const double ex, const double ey,
115 const double ez, const double bx,
116 const double by, const double bz,
117 double& wv, double& wr);
119 virtual bool ElectronDiffusion(const double ex, const double ey,
120 const double ez, const double bx,
121 const double by, const double bz, double& dl,
122 double& dt);
126 virtual bool ElectronDiffusion(const double ex, const double ey,
127 const double ez, const double bx,
128 const double by, const double bz,
129 double cov[3][3]);
131 virtual bool ElectronTownsend(const double ex, const double ey,
132 const double ez, const double bx,
133 const double by, const double bz,
134 double& alpha);
136 virtual bool ElectronAttachment(const double ex, const double ey,
137 const double ez, const double bx,
138 const double by, const double bz,
139 double& eta);
141 virtual bool ElectronTOFIonisation(const double ex, const double ey,
142 const double ez, const double bx,
143 const double by, const double bz,
144 double& riontof);
146 virtual bool ElectronTOFAttachment(const double ex, const double ey,
147 const double ez, const double bx,
148 const double by, const double bz,
149 double& ratttof);
151 virtual bool ElectronLorentzAngle(const double ex, const double ey,
152 const double ez, const double bx,
153 const double by, const double bz,
154 double& lor);
156 virtual double ElectronMobility();
157 // Microscopic electron transport properties
158
160 virtual double GetElectronEnergy(const double px, const double py,
161 const double pz, double& vx, double& vy,
162 double& vz, const int band = 0);
165 virtual void GetElectronMomentum(const double e, double& px, double& py,
166 double& pz, int& band);
167
169 virtual double GetElectronNullCollisionRate(const int band = 0);
170
172 virtual double GetElectronCollisionRate(const double e, const int band = 0);
173 struct Secondary {
175 double energy = 0.;
176 double time = 0.;
177 double distance = 0.;
178 };
179
180 virtual bool ElectronCollision(const double e, int& type, int& level,
181 double& e1, double& dx, double& dy, double& dz,
182 std::vector<Secondary>& secondaries,
183 int& band);
184
185 // Transport parameters for holes
187 virtual bool HoleVelocity(const double ex, const double ey, const double ez,
188 const double bx, const double by, const double bz,
189 double& vx, double& vy, double& vz);
191 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
192 const double bx, const double by, const double bz,
193 double& dl, double& dt);
195 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
196 const double bx, const double by, const double bz,
197 double cov[3][3]);
199 virtual bool HoleTownsend(const double ex, const double ey, const double ez,
200 const double bx, const double by, const double bz,
201 double& alpha);
203 virtual bool HoleAttachment(const double ex, const double ey, const double ez,
204 const double bx, const double by, const double bz,
205 double& eta);
207 virtual double HoleMobility();
208
209 // Transport parameters for ions
211 virtual bool IonVelocity(const double ex, const double ey, const double ez,
212 const double bx, const double by, const double bz,
213 double& vx, double& vy, double& vz);
214 bool HasIonVelocity() const { return !(m_iVel.empty() && m_iMob.empty()); }
216 virtual bool IonDiffusion(const double ex, const double ey, const double ez,
217 const double bx, const double by, const double bz,
218 double& dl, double& dt);
220 virtual bool IonDissociation(const double ex, const double ey,
221 const double ez, const double bx,
222 const double by, const double bz, double& diss);
224 virtual double IonMobility();
225
227 virtual bool NegativeIonVelocity(const double ex, const double ey,
228 const double ez, const double bx,
229 const double by, const double bz, double& vx,
230 double& vy, double& vz);
232 virtual double NegativeIonMobility();
233
235 void SetFieldGrid(double emin, double emax, const size_t ne, bool logE,
236 double bmin = 0., double bmax = 0., const size_t nb = 1,
237 double amin = HalfPi, double amax = HalfPi,
238 const size_t na = 1);
240 void SetFieldGrid(const std::vector<double>& efields,
241 const std::vector<double>& bfields,
242 const std::vector<double>& angles);
244 void GetFieldGrid(std::vector<double>& efields, std::vector<double>& bfields,
245 std::vector<double>& angles);
246
248 bool SetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia,
249 const double v) {
250 return SetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
251 }
252
253 bool GetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia,
254 double& v) {
255 return GetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
256 }
257
258 bool SetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia,
259 const double v) {
260 return SetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
261 }
262
263 bool GetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia,
264 double& v) {
265 return GetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
266 }
267
268 bool SetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia,
269 const double v) {
270 return SetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
271 }
272
273 bool GetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia,
274 double& v) {
275 return GetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
276 }
277
278 bool SetElectronFluxVelocity(const size_t ie, const size_t ib,
279 const size_t ia, const double v) {
280 return SetEntry(ie, ib, ia, "ElectronFluxVelocity", m_eVelWv, v);
281 }
282
283 bool GetElectronFluxVelocity(const size_t ie, const size_t ib,
284 const size_t ia, double& v) {
285 return GetEntry(ie, ib, ia, "ElectronFluxVelocity", m_eVelWv, v);
286 }
287
288 bool SetElectronBulkVelocity(const size_t ie, const size_t ib,
289 const size_t ia, const double v) {
290 return SetEntry(ie, ib, ia, "ElectronBulkVelocity", m_eVelWr, v);
291 }
292
293 bool GetElectronBulkVelocity(const size_t ie, const size_t ib,
294 const size_t ia, double& v) {
295 return GetEntry(ie, ib, ia, "ElectronBulkVelocity", m_eVelWr, v);
296 }
297
298 bool SetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
299 const size_t ia, const double dl) {
300 return SetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
301 }
302
303 bool GetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
304 const size_t ia, double& dl) {
305 return GetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
306 }
307
308 bool SetElectronTransverseDiffusion(const size_t ie, const size_t ib,
309 const size_t ia, const double dt) {
310 return SetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
311 }
312
313 bool GetElectronTransverseDiffusion(const size_t ie, const size_t ib,
314 const size_t ia, double& dt) {
315 return GetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
316 }
317
318 bool SetElectronTownsend(const size_t ie, const size_t ib, const size_t ia,
319 const double alpha) {
320 return SetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
321 }
322
323 bool GetElectronTownsend(const size_t ie, const size_t ib, const size_t ia,
324 double& alpha) {
325 return GetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
326 }
327
328 bool SetElectronAttachment(const size_t ie, const size_t ib, const size_t ia,
329 const double eta) {
330 return SetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
331 }
332
333 bool GetElectronAttachment(const size_t ie, const size_t ib, const size_t ia,
334 double& eta) {
335 return GetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
336 }
337
338 bool SetElectronTOFIonisation(const size_t ie, const size_t ib,
339 const size_t ia, const double v) {
340 return SetEntry(ie, ib, ia, "ElectronTOFIonisation", m_eRIon, v);
341 }
342
343 bool GetElectronTOFIonisation(const size_t ie, const size_t ib,
344 const size_t ia, double& v) {
345 return GetEntry(ie, ib, ia, "ElectronTOFIonisation", m_eRIon, v);
346 }
347
348 bool SetElectronTOFAttachment(const size_t ie, const size_t ib,
349 const size_t ia, const double v) {
350 return SetEntry(ie, ib, ia, "ElectronTOFAttachment", m_eRAtt, v);
351 }
352
353 bool GetElectronTOFAttachment(const size_t ie, const size_t ib,
354 const size_t ia, double& v) {
355 return GetEntry(ie, ib, ia, "ElectronTOFAttachment", m_eRAtt, v);
356 }
357
359 bool SetElectronLorentzAngle(const size_t ie, const size_t ib,
360 const size_t ia, const double lor) {
361 return SetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
362 }
363
364 bool GetElectronLorentzAngle(const size_t ie, const size_t ib,
365 const size_t ia, double& lor) {
366 return GetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
367 }
368
370 bool SetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia,
371 const double v) {
372 return SetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
373 }
374
375 bool GetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia,
376 double& v) {
377 return GetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
378 }
379
380 bool SetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia,
381 const double v) {
382 return SetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
383 }
384
385 bool GetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia,
386 double& v) {
387 return GetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
388 }
389
390 bool SetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia,
391 const double v) {
392 return SetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
393 }
394
395 bool GetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia,
396 double& v) {
397 return GetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
398 }
399
401 bool SetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
402 const size_t ia, const double dl) {
403 return SetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
404 }
405
406 bool GetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
407 const size_t ia, double& dl) {
408 return GetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
409 }
410
411 bool SetHoleTransverseDiffusion(const size_t ie, const size_t ib,
412 const size_t ia, const double dt) {
413 return SetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
414 }
415
416 bool GetHoleTransverseDiffusion(const size_t ie, const size_t ib,
417 const size_t ia, double& dt) {
418 return GetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
419 }
420
421 bool SetHoleTownsend(const size_t ie, const size_t ib, const size_t ia,
422 const double alpha) {
423 return SetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
424 }
425
426 bool GetHoleTownsend(const size_t ie, const size_t ib, const size_t ia,
427 double& alpha) {
428 return GetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
429 }
430
431 bool SetHoleAttachment(const size_t ie, const size_t ib, const size_t ia,
432 const double eta) {
433 return SetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
434 }
435
436 bool GetHoleAttachment(const size_t ie, const size_t ib, const size_t ia,
437 double& eta) {
438 return GetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
439 }
440
445 bool SetIonMobility(const std::vector<double>& fields,
446 const std::vector<double>& mobilities,
447 const bool negativeIons = false);
449 bool SetIonMobility(const size_t ie, const size_t ib, const size_t ia,
450 const double mu);
452 bool GetIonMobility(const size_t ie, const size_t ib, const size_t ia,
453 double& mu) {
454 return GetEntry(ie, ib, ia, "IonMobility", m_iMob, mu);
455 }
456
458 bool SetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
459 const size_t ia, const double dl) {
460 return SetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
461 }
462
463 bool GetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
464 const size_t ia, double& dl) {
465 return GetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
466 }
467
468 bool SetIonTransverseDiffusion(const size_t ie, const size_t ib,
469 const size_t ia, const double dt) {
470 return SetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
471 }
472
473 bool GetIonTransverseDiffusion(const size_t ie, const size_t ib,
474 const size_t ia, double& dt) {
475 return GetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
476 }
477
478 bool SetIonDissociation(const size_t ie, const size_t ib, const size_t ia,
479 const double diss) {
480 return SetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
481 }
482
483 bool GetIonDissociation(const size_t ie, const size_t ib, const size_t ia,
484 double& diss) {
485 return GetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
486 }
487
489 bool SetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia,
490 const double mu);
492 bool GetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia,
493 double& mu) {
494 return GetEntry(ie, ib, ia, "NegativeIonMobility", m_nMob, mu);
495 }
496
498 virtual void ResetTables();
499
501 m_eVelE.clear();
502 m_eVelB.clear();
503 m_eVelX.clear();
504 m_eVelWv.clear();
505 m_eVelWr.clear();
506 }
508 m_eDifL.clear();
509 m_eDifT.clear();
510 m_eDifM.clear();
511 }
512 void ResetElectronTownsend() { m_eAlp.clear(); }
513 void ResetElectronAttachment() { m_eAtt.clear(); }
515 m_eRIon.clear();
516 m_eRAtt.clear();
517 }
519
521 m_hVelE.clear();
522 m_hVelB.clear();
523 m_hVelX.clear();
524 }
526 m_hDifL.clear();
527 m_hDifT.clear();
528 m_hDifM.clear();
529 }
530 void ResetHoleTownsend() { m_hAlp.clear(); }
531 void ResetHoleAttachment() { m_hAtt.clear(); }
532
534 m_iMob.clear();
535 m_iVel.clear();
536 }
538 m_iDifL.clear();
539 m_iDifT.clear();
540 }
541 void ResetIonDissociation() { m_iDis.clear(); }
543 m_nMob.clear();
544 m_nVel.clear();
545 }
546
548 const std::vector<std::vector<std::vector<double> > >& mob,
549 std::vector<std::vector<std::vector<double> > >& vel);
550
553 void SetExtrapolationMethodVelocity(const std::string& extrLow,
554 const std::string& extrHigh);
555 void SetExtrapolationMethodDiffusion(const std::string& extrLow,
556 const std::string& extrHigh);
557 void SetExtrapolationMethodTownsend(const std::string& extrLow,
558 const std::string& extrHigh);
559 void SetExtrapolationMethodAttachment(const std::string& extrLow,
560 const std::string& extrHigh);
561 void SetExtrapolationMethodIonMobility(const std::string& extrLow,
562 const std::string& extrHigh);
563 void SetExtrapolationMethodIonDissociation(const std::string& extrLow,
564 const std::string& extrHigh);
565
567 void SetInterpolationMethodVelocity(const unsigned int intrp);
568 void SetInterpolationMethodDiffusion(const unsigned int intrp);
569 void SetInterpolationMethodTownsend(const unsigned int intrp);
570 void SetInterpolationMethodAttachment(const unsigned int intrp);
571 void SetInterpolationMethodIonMobility(const unsigned int intrp);
572 void SetInterpolationMethodIonDissociation(const unsigned int intrp);
573
574 // Scaling of fields and transport parameters.
575 virtual double ScaleElectricField(const double e) const { return e; }
576 virtual double UnScaleElectricField(const double e) const { return e; }
577 virtual double ScaleVelocity(const double v) const { return v; }
578 virtual double ScaleDiffusion(const double d) const { return d; }
579 virtual double ScaleDiffusionTensor(const double d) const { return d; }
580 virtual double ScaleTownsend(const double alpha) const { return alpha; }
581 virtual double ScaleAttachment(const double eta) const { return eta; }
582 virtual double ScaleLorentzAngle(const double lor) const { return lor; }
583 virtual double ScaleDissociation(const double diss) const { return diss; }
584
585 // Optical properties
587 virtual bool GetOpticalDataRange(double& emin, double& emax,
588 const unsigned int i = 0);
590 virtual bool GetDielectricFunction(const double e, double& eps1, double& eps2,
591 const unsigned int i = 0);
592 // Get the photoabsorption cross-section [cm2] at a given energy.
593 virtual bool GetPhotoAbsorptionCrossSection(const double e, double& sigma,
594 const unsigned int i = 0);
595 virtual double GetPhotonCollisionRate(const double e);
596 virtual bool PhotonCollision(const double e, int& type, int& level,
597 double& e1, double& ctheta,
598 std::vector<Secondary>& secondaries);
599
601 void EnableDebugging() { m_debug = true; }
602 void DisableDebugging() { m_debug = false; }
603
605 virtual double CreateGPUTransferObject(MediumGPU*& med_gpu);
606
607 protected:
608 std::string m_className = "Medium";
609
610 static int m_idCounter;
611
612 // Number of components
613 unsigned int m_nComponents = 1;
614 // Name
615 std::string m_name = "";
616 // Temperature [K]
617 double m_temperature = 293.15;
618 // Pressure [Torr]
619 double m_pressure = 760.;
620 // Static dielectric constant
621 double m_epsilon = 1.;
622 // (Effective) atomic number Z
623 double m_z = 1.;
624 // Atomic weight A
625 double m_a = 0.;
626 // Number density [cm-3]
627 double m_density = 0.;
628
629 // Id number
630 int m_id;
631
632 // Transport flags
633 bool m_driftable = false;
634 bool m_microscopic = false;
635 bool m_ionisable = false;
636
637 // W value
638 double m_w = 0.;
639 // Fano factor
640 double m_fano = 0.;
641
642 // Update flag
643 bool m_isChanged = true;
644
645 // Switch on/off debugging messages
646 bool m_debug = false;
647
648 // Tables of transport parameters
649 bool m_tab2d = false;
650
651 // Field grids
652 std::vector<double> m_eFields;
653 std::vector<double> m_bFields;
654 std::vector<double> m_bAngles;
655
656 // Electrons
657 std::vector<std::vector<std::vector<double> > > m_eVelE;
658 std::vector<std::vector<std::vector<double> > > m_eVelX;
659 std::vector<std::vector<std::vector<double> > > m_eVelB;
660 std::vector<std::vector<std::vector<double> > > m_eDifL;
661 std::vector<std::vector<std::vector<double> > > m_eDifT;
662 std::vector<std::vector<std::vector<double> > > m_eAlp;
663 std::vector<std::vector<std::vector<double> > > m_eAtt;
664 std::vector<std::vector<std::vector<double> > > m_eLor;
665 std::vector<std::vector<std::vector<double> > > m_eVelWv;
666 std::vector<std::vector<std::vector<double> > > m_eVelWr;
667 std::vector<std::vector<std::vector<double> > > m_eRIon;
668 std::vector<std::vector<std::vector<double> > > m_eRAtt;
669
670 std::vector<std::vector<std::vector<std::vector<double> > > > m_eDifM;
671
672 // Holes
673 std::vector<std::vector<std::vector<double> > > m_hVelE;
674 std::vector<std::vector<std::vector<double> > > m_hVelX;
675 std::vector<std::vector<std::vector<double> > > m_hVelB;
676 std::vector<std::vector<std::vector<double> > > m_hDifL;
677 std::vector<std::vector<std::vector<double> > > m_hDifT;
678 std::vector<std::vector<std::vector<double> > > m_hAlp;
679 std::vector<std::vector<std::vector<double> > > m_hAtt;
680
681 std::vector<std::vector<std::vector<std::vector<double> > > > m_hDifM;
682
683 // Ions
684 std::vector<std::vector<std::vector<double> > > m_iMob;
685 std::vector<std::vector<std::vector<double> > > m_iVel;
686 std::vector<std::vector<std::vector<double> > > m_iDifL;
687 std::vector<std::vector<std::vector<double> > > m_iDifT;
688 std::vector<std::vector<std::vector<double> > > m_iDis;
689 // Negative ions
690 std::vector<std::vector<std::vector<double> > > m_nMob;
691 std::vector<std::vector<std::vector<double> > > m_nVel;
692
693 // Thresholds for Townsend, attachment and dissociation coefficients.
694 unsigned int m_eThrAlp = 0;
695 unsigned int m_eThrAtt = 0;
696 unsigned int m_hThrAlp = 0;
697 unsigned int m_hThrAtt = 0;
698 unsigned int m_iThrDis = 0;
699
700 // Extrapolation methods (TODO: enum).
701 std::pair<unsigned int, unsigned int> m_extrVel = {0, 1};
702 std::pair<unsigned int, unsigned int> m_extrDif = {0, 1};
703 std::pair<unsigned int, unsigned int> m_extrAlp = {0, 1};
704 std::pair<unsigned int, unsigned int> m_extrAtt = {0, 1};
705 std::pair<unsigned int, unsigned int> m_extrLor = {0, 1};
706 std::pair<unsigned int, unsigned int> m_extrMob = {0, 1};
707 std::pair<unsigned int, unsigned int> m_extrDis = {0, 1};
708
709 // Interpolation methods
710 unsigned int m_intpVel = 2;
711 unsigned int m_intpDif = 2;
712 unsigned int m_intpAlp = 2;
713 unsigned int m_intpAtt = 2;
714 unsigned int m_intpLor = 2;
715 unsigned int m_intpMob = 2;
716 unsigned int m_intpDis = 2;
717
718 bool Velocity(const double ex, const double ey, const double ez,
719 const double bx, const double by, const double bz,
720 const std::vector<std::vector<std::vector<double> > >& velE,
721 const std::vector<std::vector<std::vector<double> > >& velB,
722 const std::vector<std::vector<std::vector<double> > >& velX,
723 const double q, double& vx, double& vy, double& vz) const;
725 const double ex, const double ey, const double ez, const double bx,
726 const double by, const double bz,
727 const std::vector<std::vector<std::vector<double> > >& velWv,
728 const std::vector<std::vector<std::vector<double> > >& velWr, double& wv,
729 double& wr) const;
730 static void Langevin(const double ex, const double ey, const double ez,
731 double bx, double by, double bz, const double mu,
732 double& vx, double& vy, double& vz);
733 static void Langevin(const double ex, const double ey, const double ez,
734 double bx, double by, double bz, const double mu,
735 const double muH, double& vx, double& vy, double& vz);
736 bool Diffusion(const double ex, const double ey, const double ez,
737 const double bx, const double by, const double bz,
738 const std::vector<std::vector<std::vector<double> > >& difL,
739 const std::vector<std::vector<std::vector<double> > >& difT,
740 double& dl, double& dt) const;
742 const double ex, const double ey, const double ez, const double bx,
743 const double by, const double bz,
744 const std::vector<std::vector<std::vector<std::vector<double> > > >& diff,
745 double cov[3][3]) const;
746 bool Alpha(const double ex, const double ey, const double ez, const double bx,
747 const double by, const double bz,
748 const std::vector<std::vector<std::vector<double> > >& tab,
749 unsigned int intp, const unsigned int thr,
750 const std::pair<unsigned int, unsigned int>& extr,
751 double& alpha) const;
752 double GetAngle(const double ex, const double ey, const double ez,
753 const double bx, const double by, const double bz,
754 const double e, const double b) const;
755 bool Interpolate(const double e, const double b, const double a,
756 const std::vector<std::vector<std::vector<double> > >& table,
757 double& y, const unsigned int intp,
758 const std::pair<unsigned int, unsigned int>& extr,
759 const bool logval = false) const;
760
761 double Interpolate1D(const double e, const std::vector<double>& table,
762 const std::vector<double>& fields,
763 const unsigned int intpMeth,
764 const std::pair<unsigned int, unsigned int>& extr,
765 const bool logval = false) const;
766
767 bool SetEntry(const size_t i, const size_t j, const size_t k,
768 const std::string& fcn,
769 std::vector<std::vector<std::vector<double> > >& tab,
770 const double val);
771 bool GetEntry(const size_t i, const size_t j, const size_t k,
772 const std::string& fcn,
773 const std::vector<std::vector<std::vector<double> > >& tab,
774 double& val) const;
775
776 void SetExtrapolationMethod(const std::string& low, const std::string& high,
777 std::pair<unsigned int, unsigned int>& extr,
778 const std::string& fcn);
779 bool GetExtrapolationIndex(std::string str, unsigned int& nb) const;
781 const std::vector<std::vector<std::vector<double> > >& tab) const;
782
783 void Clone(std::vector<std::vector<std::vector<double> > >& tab,
784 const std::vector<double>& efields,
785 const std::vector<double>& bfields,
786 const std::vector<double>& angles, const unsigned int intp,
787 const std::pair<unsigned int, unsigned int>& extr,
788 const double init, const std::string& label);
789 void Clone(std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
790 const size_t n, const std::vector<double>& efields,
791 const std::vector<double>& bfields,
792 const std::vector<double>& angles, const unsigned int intp,
793 const std::pair<unsigned int, unsigned int>& extr,
794 const double init, const std::string& label);
795
796 void Init(const size_t nE, const size_t nB, const size_t nA,
797 std::vector<std::vector<std::vector<double> > >& tab,
798 const double val);
799 void Init(const size_t nE, const size_t nB, const size_t nA, const size_t nT,
800 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
801 const double val);
802};
803
804} // namespace Garfield
805#endif
bool GetElectronTOFAttachment(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of attachment rate of TOF.
Definition Medium.hh:353
void ResetHoleAttachment()
Definition Medium.hh:531
bool SetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia, const double mu)
Set an entry in the table of negative ion mobilities.
double m_density
Definition Medium.hh:627
void ResetIonMobility()
Definition Medium.hh:533
double GetW() const
Get the W value.
Definition Medium.hh:89
std::vector< std::vector< std::vector< double > > > m_eRIon
Definition Medium.hh:667
bool HasIonVelocity() const
Definition Medium.hh:214
virtual void SetNumberDensity(const double n)
Set the number density [cm-3].
virtual double UnScaleElectricField(const double e) const
Definition Medium.hh:576
bool GetHoleLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dl)
Get an entry in the table of longitudinal diffusion coefficients.
Definition Medium.hh:406
virtual double ScaleVelocity(const double v) const
Definition Medium.hh:577
bool GetHoleTownsend(const size_t ie, const size_t ib, const size_t ia, double &alpha)
Get an entry in the table of Townsend coefficients.
Definition Medium.hh:426
void ResetElectronAttachment()
Definition Medium.hh:513
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].
std::vector< double > m_bFields
Definition Medium.hh:653
void ResetIonDiffusion()
Definition Medium.hh:537
bool IsMicroscopic() const
Does the medium have electron scattering rates?
Definition Medium.hh:81
virtual void GetElectronMomentum(const double e, double &px, double &py, double &pz, int &band)
Sample the momentum vector for a given energy (only meaningful in semiconductors).
void SetTemperature(const double t)
Set the temperature [K].
virtual bool PhotonCollision(const double e, int &type, int &level, double &e1, double &ctheta, std::vector< Secondary > &secondaries)
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)
Ion drift velocity [cm / ns].
void ResetElectronLorentzAngle()
Definition Medium.hh:518
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].
double m_pressure
Definition Medium.hh:619
void SetFanoFactor(const double f)
Set the Fano factor.
Definition Medium.hh:91
virtual double GetElectronNullCollisionRate(const int band=0)
Null-collision rate [ns-1].
bool GetElectronLorentzAngle(const size_t ie, const size_t ib, const size_t ia, double &lor)
Get an entry in the table of Lorentz angles.
Definition Medium.hh:364
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].
unsigned int m_intpMob
Definition Medium.hh:715
virtual double ScaleAttachment(const double eta) const
Definition Medium.hh:581
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
void SetExtrapolationMethodIonDissociation(const std::string &extrLow, const std::string &extrHigh)
void SetInterpolationMethodIonDissociation(const unsigned int intrp)
bool SetElectronLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, const double dl)
Set an entry in the table of longitudinal diffusion coefficients.
Definition Medium.hh:298
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
virtual double ScaleDiffusion(const double d) const
Definition Medium.hh:578
virtual double GetNumberDensity() const
Get the number density [cm-3].
Definition Medium.hh:65
bool GetHoleAttachment(const size_t ie, const size_t ib, const size_t ia, double &eta)
Get an entry in the table of attachment coefficients.
Definition Medium.hh:436
bool SetIonMobility(const size_t ie, const size_t ib, const size_t ia, const double mu)
Set an entry in the table of ion mobilities.
bool SetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of drift speeds along Btrans.
Definition Medium.hh:268
virtual double ScaleLorentzAngle(const double lor) const
Definition Medium.hh:582
unsigned int m_intpDis
Definition Medium.hh:716
bool GetIonTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dt)
Get an entry in the table of transverse diffusion coefficients.
Definition Medium.hh:473
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].
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
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition Medium.hh:662
bool GetElectronTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dt)
Get an entry in the table of transverse diffusion coefficients.
Definition Medium.hh:313
void ResetHoleVelocity()
Definition Medium.hh:520
virtual void ResetTables()
Reset all tables of transport parameters.
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.
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition Medium.hh:678
void ResetIonDissociation()
Definition Medium.hh:541
virtual void SetAtomicNumber(const double z)
Set the effective atomic number.
bool SetElectronAttachment(const size_t ie, const size_t ib, const size_t ia, const double eta)
Set an entry in the table of attachment coefficients.
Definition Medium.hh:328
virtual void EnableDrift(const bool on=true)
Switch electron/ion/hole transport on/off.
Definition Medium.hh:72
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].
bool SetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of drift speeds along E.
Definition Medium.hh:370
virtual bool ElectronTOFAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &ratttof)
TOF Attachment Rate [ns-1].
bool GetElectronBulkVelocity(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of bulk drift speeds.
Definition Medium.hh:293
double GetDielectricConstant() const
Get the relative static dielectric constant.
Definition Medium.hh:47
bool SetHoleLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, const double dl)
Set an entry in the table of longitudinal diffusion coefficients.
Definition Medium.hh:401
std::string m_name
Definition Medium.hh:615
virtual double HoleMobility()
Low-field mobility [cm2 V-1 ns-1].
static void Langevin(const double ex, const double ey, const double ez, double bx, double by, double bz, const double mu, double &vx, double &vy, double &vz)
std::vector< std::vector< std::vector< double > > > m_eVelWv
Definition Medium.hh:665
std::vector< std::vector< std::vector< double > > > m_eRAtt
Definition Medium.hh:668
double GetPressure() const
Definition Medium.hh:43
unsigned int m_intpDif
Definition Medium.hh:711
virtual void EnablePrimaryIonisation(const bool on=true)
Make the medium ionisable or non-ionisable.
Definition Medium.hh:74
std::pair< unsigned int, unsigned int > m_extrAtt
Definition Medium.hh:704
int GetId() const
Return the id number of the class instance.
Definition Medium.hh:25
void SetExtrapolationMethodTownsend(const std::string &extrLow, const std::string &extrHigh)
virtual double GetElectronCollisionRate(const double e, const int band=0)
Collision rate [ns-1] for given electron energy.
bool SetHoleAttachment(const size_t ie, const size_t ib, const size_t ia, const double eta)
Set an entry in the table of attachment coefficients.
Definition Medium.hh:431
virtual double GetPhotonCollisionRate(const double e)
bool SetElectronTownsend(const size_t ie, const size_t ib, const size_t ia, const double alpha)
Set an entry in the table of Townsend coefficients.
Definition Medium.hh:318
virtual double GetAtomicWeight() const
Get the effective atomic weight.
Definition Medium.hh:61
bool GetHoleTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dt)
Get an entry in the table of transverse diffusion coefficients.
Definition Medium.hh:416
virtual ~Medium()
Destructor.
virtual void SetAtomicWeight(const double a)
Set the effective atomic weight.
std::vector< std::vector< std::vector< double > > > m_iDifT
Definition Medium.hh:687
void SetExtrapolationMethodIonMobility(const std::string &extrLow, const std::string &extrHigh)
void PlotVelocity(const std::string &carriers, TPad *pad)
Plot the drift velocity as function of the electric field.
unsigned int m_hThrAtt
Definition Medium.hh:697
virtual bool GetPhotoAbsorptionCrossSection(const double e, double &sigma, const unsigned int i=0)
std::pair< unsigned int, unsigned int > m_extrVel
Definition Medium.hh:701
virtual double ScaleElectricField(const double e) const
Definition Medium.hh:575
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 bool logval=false) const
virtual double ScaleDiffusionTensor(const double d) const
Definition Medium.hh:579
unsigned int GetNumberOfComponents() const
Get number of components of the medium.
Definition Medium.hh:50
std::vector< std::vector< std::vector< double > > > m_hDifL
Definition Medium.hh:676
virtual double GetAtomicNumber() const
Get the effective atomic number.
Definition Medium.hh:57
virtual bool ElectronVelocityFluxBulk(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &wv, double &wr)
Flux (mean velocity; shorthand: wv) and bulk (center of mass velocity; shorthand: wr) drift velocity ...
size_t SetThreshold(const std::vector< std::vector< std::vector< double > > > &tab) const
virtual bool NegativeIonVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Negative ion drift velocity [cm / ns].
std::vector< std::vector< std::vector< double > > > m_nVel
Definition Medium.hh:691
virtual bool ElectronCollision(const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< Secondary > &secondaries, int &band)
Sample the collision type. Update energy and direction vector.
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition Medium.hh:657
std::vector< std::vector< std::vector< double > > > m_eVelX
Definition Medium.hh:658
std::vector< std::vector< std::vector< double > > > m_eDifL
Definition Medium.hh:660
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)
bool GetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along ExB.
Definition Medium.hh:385
virtual double ScaleTownsend(const double alpha) const
Definition Medium.hh:580
virtual bool HoleDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double cov[3][3])
Diffusion tensor.
unsigned int m_intpVel
Definition Medium.hh:710
void SetDielectricConstant(const double eps)
Set the relative static dielectric constant.
bool GetElectronAttachment(const size_t ie, const size_t ib, const size_t ia, double &eta)
Get an entry in the table of attachment coefficients.
Definition Medium.hh:333
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< std::vector< double > > > > &diff, double cov[3][3]) const
virtual double GetMassDensity() const
Get the mass density [g/cm3].
bool GetExtrapolationIndex(std::string str, unsigned int &nb) const
void SetInterpolationMethodVelocity(const unsigned int intrp)
Set the degree of polynomial interpolation (usually 2).
void VelocityFromMobility(const std::vector< std::vector< std::vector< double > > > &mob, std::vector< std::vector< std::vector< double > > > &vel)
void SetInterpolationMethodDiffusion(const unsigned int intrp)
void ResetElectronTOFRates()
Definition Medium.hh:514
virtual bool GetDielectricFunction(const double e, double &eps1, double &eps2, const unsigned int i=0)
Get the complex dielectric function at a given energy.
bool SetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of drift speeds along Btrans.
Definition Medium.hh:390
double m_epsilon
Definition Medium.hh:621
void SetExtrapolationMethodDiffusion(const std::string &extrLow, const std::string &extrHigh)
std::vector< double > m_eFields
Definition Medium.hh:652
bool GetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along E.
Definition Medium.hh:253
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.
bool SetElectronTOFAttachment(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of attachment rate of TOF.
Definition Medium.hh:348
virtual bool ElectronDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double cov[3][3])
Diffusion tensor: diagonal elements are the diffusion coefficients [cm] along e, btrans,...
static int m_idCounter
Definition Medium.hh:610
std::vector< std::vector< std::vector< double > > > m_hDifT
Definition Medium.hh:677
std::vector< std::vector< std::vector< double > > > m_eVelWr
Definition Medium.hh:666
void DisableDebugging()
Definition Medium.hh:602
void ResetElectronDiffusion()
Definition Medium.hh:507
virtual bool IsSemiconductor() const
Is this medium a semiconductor?
Definition Medium.hh:32
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.
bool SetElectronLorentzAngle(const size_t ie, const size_t ib, const size_t ia, const double lor)
Set an entry in the table of Lorentz angles.
Definition Medium.hh:359
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition Medium.hh:663
bool GetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along Btrans.
Definition Medium.hh:395
void PlotAlphaEta(const std::string &carriers, TPad *pad)
Plot Townsend and attachment coefficients.
bool GetElectronTownsend(const size_t ie, const size_t ib, const size_t ia, double &alpha)
Get an entry in the table of Townsend coefficients.
Definition Medium.hh:323
bool SetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of drift speeds along ExB.
Definition Medium.hh:258
bool SetIonLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, const double dl)
Set an entry in the table of longitudinal diffusion coefficients.
Definition Medium.hh:458
virtual double ScaleDissociation(const double diss) const
Definition Medium.hh:583
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
const std::string & GetName() const
Get the medium name/identifier.
Definition Medium.hh:28
std::vector< double > m_bAngles
Definition Medium.hh:654
Medium()
Constructor.
std::vector< std::vector< std::vector< double > > > m_iDifL
Definition Medium.hh:686
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition Medium.hh:673
virtual double IonMobility()
Low-field ion mobility [cm2 V-1 ns-1].
virtual bool IsGas() const
Is this medium a gas?
Definition Medium.hh:30
std::vector< std::vector< std::vector< double > > > m_eLor
Definition Medium.hh:664
void Clone(std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const size_t n, 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)
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].
std::vector< std::vector< std::vector< double > > > m_eDifT
Definition Medium.hh:661
void SetInterpolationMethodIonMobility(const unsigned int intrp)
unsigned int m_intpAtt
Definition Medium.hh:713
void SetExtrapolationMethodVelocity(const std::string &extrLow, const std::string &extrHigh)
Select the extrapolation method for fields below/above the table range.
bool GetIonMobility(const size_t ie, const size_t ib, const size_t ia, double &mu)
Get an entry in the table of ion mobilities.
Definition Medium.hh:452
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.
bool SetElectronBulkVelocity(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of bulk drift speeds.
Definition Medium.hh:288
bool GetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along ExB.
Definition Medium.hh:263
bool SetElectronTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, const double dt)
Set an entry in the table of transverse diffusion coefficients.
Definition Medium.hh:308
std::pair< unsigned int, unsigned int > m_extrDis
Definition Medium.hh:707
void SetFieldGrid(const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles)
Set the fields and E-B angles to be used in the transport tables.
void PlotDiffusion(const std::string &carriers, TPad *pad)
Plot the diffusion coefficients as function of the electric field.
std::vector< std::vector< std::vector< double > > > m_hAtt
Definition Medium.hh:679
void SetPressure(const double p)
void SetExtrapolationMethod(const std::string &low, const std::string &high, std::pair< unsigned int, unsigned int > &extr, const std::string &fcn)
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].
std::vector< std::vector< std::vector< double > > > m_iVel
Definition Medium.hh:685
std::pair< unsigned int, unsigned int > m_extrAlp
Definition Medium.hh:703
void EnableDebugging()
Switch on/off debugging messages.
Definition Medium.hh:601
std::vector< std::vector< std::vector< double > > > m_eVelB
Definition Medium.hh:659
bool VelocityFluxBulk(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 > > > &velWv, const std::vector< std::vector< std::vector< double > > > &velWr, double &wv, double &wr) const
unsigned int m_eThrAtt
Definition Medium.hh:695
bool IsIonisable() const
Is charge deposition by charged particles/photon enabled in this medium?
Definition Medium.hh:84
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].
bool GetIonDissociation(const size_t ie, const size_t ib, const size_t ia, double &diss)
Get an entry in the table of dissociation coefficients.
Definition Medium.hh:483
void PlotAttachment(const std::string &carriers, TPad *pad)
Plot the attachment coefficient(s) as function of the electric field.
bool SetIonMobility(const std::vector< double > &fields, const std::vector< double > &mobilities, const bool negativeIons=false)
Initialise the table of ion mobilities from a list of electric fields and corresponding mobilities.
void SetExtrapolationMethodAttachment(const std::string &extrLow, const std::string &extrHigh)
bool GetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along Btrans.
Definition Medium.hh:273
virtual void GetComponent(const unsigned int i, std::string &label, double &f)
Get the name and fraction of a given component.
void ResetNegativeIonMobility()
Definition Medium.hh:542
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)
std::pair< unsigned int, unsigned int > m_extrLor
Definition Medium.hh:705
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].
virtual bool GetOpticalDataRange(double &emin, double &emax, const unsigned int i=0)
Get the energy range [eV] of the available optical data.
bool GetIonLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dl)
Get an entry in the table of longitudinal diffusion coefficients.
Definition Medium.hh:463
unsigned int m_nComponents
Definition Medium.hh:613
virtual bool IsConductor() const
Is this medium a conductor?
Definition Medium.hh:34
unsigned int m_intpLor
Definition Medium.hh:714
std::string m_className
Definition Medium.hh:608
bool SetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of drift speeds along ExB.
Definition Medium.hh:380
std::pair< unsigned int, unsigned int > m_extrDif
Definition Medium.hh:702
void SetInterpolationMethodAttachment(const unsigned int intrp)
virtual void SetMassDensity(const double rho)
Set the mass density [g/cm3].
virtual double ElectronMobility()
Low-field mobility [cm2 V-1 ns-1].
double GetTemperature() const
Get the temperature [K].
Definition Medium.hh:39
bool SetElectronFluxVelocity(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of flux drift speeds.
Definition Medium.hh:278
unsigned int m_iThrDis
Definition Medium.hh:698
bool GetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia, double &mu)
Get an entry in the table of negative ion mobilities.
Definition Medium.hh:492
static void Langevin(const double ex, const double ey, const double ez, double bx, double by, double bz, const double mu, const double muH, double &vx, double &vy, double &vz)
bool GetElectronLongitudinalDiffusion(const size_t ie, const size_t ib, const size_t ia, double &dl)
Get an entry in the table of longitudinal diffusion coefficients.
Definition Medium.hh:303
std::vector< std::vector< std::vector< double > > > m_hVelX
Definition Medium.hh:674
std::pair< unsigned int, unsigned int > m_extrMob
Definition Medium.hh:706
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)
void ResetHoleDiffusion()
Definition Medium.hh:525
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
Definition Medium.hh:670
std::vector< std::vector< std::vector< double > > > m_hVelB
Definition Medium.hh:675
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
bool SetHoleTownsend(const size_t ie, const size_t ib, const size_t ia, const double alpha)
Set an entry in the table of Townsend coefficients.
Definition Medium.hh:421
std::vector< std::vector< std::vector< std::vector< double > > > > m_hDifM
Definition Medium.hh:681
virtual double NegativeIonMobility()
Low-field negative ion mobility [cm2 V-1 ns-1].
unsigned int m_hThrAlp
Definition Medium.hh:696
void SetW(const double w)
Set the W value (average energy to produce an electron/ion or e/h pair).
Definition Medium.hh:87
std::vector< std::vector< std::vector< double > > > m_iMob
Definition Medium.hh:684
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)
bool SetIonDissociation(const size_t ie, const size_t ib, const size_t ia, const double diss)
Set an entry in the table of dissociation coefficients.
Definition Medium.hh:478
bool SetHoleTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, const double dt)
Set an entry in the table of transverse diffusion coefficients.
Definition Medium.hh:411
void ResetElectronTownsend()
Definition Medium.hh:512
unsigned int m_eThrAlp
Definition Medium.hh:694
double GetFanoFactor() const
Get the Fano factor.
Definition Medium.hh:93
void PlotTownsend(const std::string &carriers, TPad *pad)
Plot the Townsend coefficient(s) as function of the electric field.
bool GetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of drift speeds along E.
Definition Medium.hh:375
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 bool logval=false) const
bool GetElectronFluxVelocity(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of flux drift speeds.
Definition Medium.hh:283
bool IsDriftable() const
Is charge carrier transport enabled in this medium?
Definition Medium.hh:79
void ResetHoleTownsend()
Definition Medium.hh:530
void ResetElectronVelocity()
Definition Medium.hh:500
void SetInterpolationMethodTownsend(const unsigned int intrp)
virtual double CreateGPUTransferObject(MediumGPU *&med_gpu)
Create and initialise GPU Transfer class.
bool SetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of drift speeds along E.
Definition Medium.hh:248
virtual bool ElectronTOFIonisation(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &riontof)
TOF Ionisation rate [ns-1].
void Init(const size_t nE, const size_t nB, const size_t nA, const size_t nT, std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const double val)
double m_temperature
Definition Medium.hh:617
bool SetIonTransverseDiffusion(const size_t ie, const size_t ib, const size_t ia, const double dt)
Set an entry in the table of transverse diffusion coefficients.
Definition Medium.hh:468
bool SetElectronTOFIonisation(const size_t ie, const size_t ib, const size_t ia, const double v)
Set an entry in the table of ionization rate of TOF.
Definition Medium.hh:338
std::vector< std::vector< std::vector< double > > > m_iDis
Definition Medium.hh:688
bool GetElectronTOFIonisation(const size_t ie, const size_t ib, const size_t ia, double &v)
Get an entry in the table of ionization rate of TOF.
Definition Medium.hh:343
std::vector< std::vector< std::vector< double > > > m_nMob
Definition Medium.hh:690
unsigned int m_intpAlp
Definition Medium.hh:712