Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation 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
14/// Abstract base class for media.
15
16class Medium {
17 public:
18 /// Constructor
19 Medium();
20 /// Destructor
21 virtual ~Medium();
22
23 /// Return the id number of the class instance.
24 int GetId() const { return m_id; }
25 /// Get the medium name/identifier.
26 const std::string& GetName() const { return m_name; }
27 /// Is this medium a gas?
28 virtual bool IsGas() const { return false; }
29 /// Is this medium a semiconductor?
30 virtual bool IsSemiconductor() const { return false; }
31 /// Is this medium a conductor?
32 virtual bool IsConductor() const { return false; }
33
34 /// Set the temperature [K].
35 void SetTemperature(const double t);
36 /// Get the temperature [K].
37 double GetTemperature() const { return m_temperature; }
38 // Set the pressure [Torr].
39 void SetPressure(const double p);
40 // Get the pressure [Torr].
41 double GetPressure() const { return m_pressure; }
42 /// Set the relative static dielectric constant.
43 void SetDielectricConstant(const double eps);
44 /// Get the relative static dielectric constant.
45 double GetDielectricConstant() const { return m_epsilon; }
46
47 /// Get number of components of the medium.
48 unsigned int GetNumberOfComponents() const { return m_nComponents; }
49 /// Get the name and fraction of a given component.
50 virtual void GetComponent(const unsigned int i, std::string& label,
51 double& f);
52 /// Set the effective atomic number.
53 virtual void SetAtomicNumber(const double z);
54 /// Get the effective atomic number.
55 virtual double GetAtomicNumber() const { return m_z; }
56 /// Set the effective atomic weight.
57 virtual void SetAtomicWeight(const double a);
58 /// Get the effective atomic weight.
59 virtual double GetAtomicWeight() const { return m_a; }
60 /// Set the number density [cm-3].
61 virtual void SetNumberDensity(const double n);
62 /// Get the number density [cm-3].
63 virtual double GetNumberDensity() const { return m_density; }
64 /// Set the mass density [g/cm3].
65 virtual void SetMassDensity(const double rho);
66 /// Get the mass density [g/cm3].
67 virtual double GetMassDensity() const;
68
69 /// Switch electron/ion/hole transport on/off.
70 virtual void EnableDrift(const bool on = true) { m_driftable = on; }
71 /// Make the medium ionisable or non-ionisable.
72 virtual void EnablePrimaryIonisation(const bool on = true) {
73 m_ionisable = on;
74 }
75
76 /// Is charge carrier transport enabled in this medium?
77 bool IsDriftable() const { return m_driftable; }
78 /// Does the medium have electron scattering rates?
79 bool IsMicroscopic() const { return m_microscopic; }
80 /// Is charge deposition by charged particles/photon enabled in this medium?
81 bool IsIonisable() const { return m_ionisable; }
82
83 /// Set the W value (average energy to produce an electron/ion or e/h pair).
84 void SetW(const double w) { m_w = w; }
85 /// Get the W value.
86 double GetW() const { return m_w; }
87 /// Set the Fano factor.
88 void SetFanoFactor(const double f) { m_fano = f; }
89 /// Get the Fano factor.
90 double GetFanoFactor() const { return m_fano; }
91
92 /// Plot the drift velocity as function of the electric field.
93 void PlotVelocity(const std::string& carriers, TPad* pad);
94 /// Plot the diffusion coefficients as function of the electric field.
95 void PlotDiffusion(const std::string& carriers, TPad* pad);
96 /// Plot the Townsend coefficient(s) as function of the electric field.
97 void PlotTownsend(const std::string& carriers, TPad* pad);
98 /// Plot the attachment coefficient(s) as function of the electric field.
99 void PlotAttachment(const std::string& carriers, TPad* pad);
100 /// Plot Townsend and attachment coefficients.
101 void PlotAlphaEta(const std::string& carriers, TPad* pad);
102
103 // Transport parameters for electrons
104 /// Drift velocity [cm / ns]
105 virtual bool ElectronVelocity(const double ex, const double ey,
106 const double ez, const double bx,
107 const double by, const double bz, double& vx,
108 double& vy, double& vz);
109 /// Longitudinal and transverse diffusion coefficients [cm1/2]
110 virtual bool ElectronDiffusion(const double ex, const double ey,
111 const double ez, const double bx,
112 const double by, const double bz, double& dl,
113 double& dt);
114 /// Diffusion tensor: diagonal elements are the diffusion
115 /// coefficients [cm] along e, btrans, e x b,
116 /// off-diagonal elements are the covariances
117 virtual bool ElectronDiffusion(const double ex, const double ey,
118 const double ez, const double bx,
119 const double by, const double bz,
120 double cov[3][3]);
121 /// Ionisation coefficient [cm-1]
122 virtual bool ElectronTownsend(const double ex, const double ey,
123 const double ez, const double bx,
124 const double by, const double bz,
125 double& alpha);
126 /// Attachment coefficient [cm-1]
127 virtual bool ElectronAttachment(const double ex, const double ey,
128 const double ez, const double bx,
129 const double by, const double bz,
130 double& eta);
131 /// Lorentz angle
132 virtual bool ElectronLorentzAngle(const double ex, const double ey,
133 const double ez, const double bx,
134 const double by, const double bz,
135 double& lor);
136 /// Low-field mobility [cm2 V-1 ns-1]
137 virtual double ElectronMobility();
138 // Microscopic electron transport properties
139
140 /// Dispersion relation (energy vs. wave vector)
141 virtual double GetElectronEnergy(const double px, const double py,
142 const double pz, double& vx, double& vy,
143 double& vz, const int band = 0);
144 /// Sample the momentum vector for a given energy
145 /// (only meaningful in semiconductors).
146 virtual void GetElectronMomentum(const double e, double& px, double& py,
147 double& pz, int& band);
148
149 /// Null-collision rate [ns-1]
150 virtual double GetElectronNullCollisionRate(const int band = 0);
151 /// Collision rate [ns-1] for given electron energy
152 virtual double GetElectronCollisionRate(const double e, const int band = 0);
153 /// Sample the collision type. Update energy and direction vector.
154 virtual bool ElectronCollision(
155 const double e, int& type, int& level, double& e1,
156 double& dx, double& dy, double& dz,
157 std::vector<std::pair<Particle, double> >& secondaries, int& ndxc,
158 int& band);
159 virtual unsigned int GetNumberOfDeexcitationProducts() const { return 0; }
160 virtual bool GetDeexcitationProduct(const unsigned int i, double& t,
161 double& s, int& type,
162 double& energy) const;
163
164 // Transport parameters for holes
165 /// Drift velocity [cm / ns]
166 virtual bool HoleVelocity(const double ex, const double ey, const double ez,
167 const double bx, const double by, const double bz,
168 double& vx, double& vy, double& vz);
169 /// Longitudinal and transverse diffusion coefficients [cm1/2]
170 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
171 const double bx, const double by, const double bz,
172 double& dl, double& dt);
173 /// Diffusion tensor
174 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
175 const double bx, const double by, const double bz,
176 double cov[3][3]);
177 /// Ionisation coefficient [cm-1]
178 virtual bool HoleTownsend(const double ex, const double ey, const double ez,
179 const double bx, const double by, const double bz,
180 double& alpha);
181 /// Attachment coefficient [cm-1]
182 virtual bool HoleAttachment(const double ex, const double ey, const double ez,
183 const double bx, const double by, const double bz,
184 double& eta);
185 /// Low-field mobility [cm2 V-1 ns-1]
186 virtual double HoleMobility();
187
188 // Transport parameters for ions
189 /// Ion drift velocity [cm / ns]
190 virtual bool IonVelocity(const double ex, const double ey, const double ez,
191 const double bx, const double by, const double bz,
192 double& vx, double& vy, double& vz);
193 /// Longitudinal and transverse diffusion coefficients [cm1/2]
194 virtual bool IonDiffusion(const double ex, const double ey, const double ez,
195 const double bx, const double by, const double bz,
196 double& dl, double& dt);
197 /// Dissociation coefficient
198 virtual bool IonDissociation(const double ex, const double ey,
199 const double ez, const double bx,
200 const double by, const double bz, double& diss);
201 /// Low-field ion mobility [cm2 V-1 ns-1]
202 virtual double IonMobility();
203
204 /// Negative ion drift velocity [cm / ns]
205 virtual bool NegativeIonVelocity(
206 const double ex, const double ey, const double ez,
207 const double bx, const double by, const double bz,
208 double& vx, double& vy, double& vz);
209 /// Low-field negative ion mobility [cm2 V-1 ns-1]
210 virtual double NegativeIonMobility();
211
212 /// Set the range of fields to be covered by the transport tables.
213 void SetFieldGrid(double emin, double emax, const size_t ne, bool logE,
214 double bmin = 0., double bmax = 0., const size_t nb = 1,
215 double amin = HalfPi, double amax = HalfPi,
216 const size_t na = 1);
217 /// Set the fields and E-B angles to be used in the transport tables.
218 void SetFieldGrid(const std::vector<double>& efields,
219 const std::vector<double>& bfields,
220 const std::vector<double>& angles);
221 /// Get the fields and E-B angles used in the transport tables.
222 void GetFieldGrid(std::vector<double>& efields, std::vector<double>& bfields,
223 std::vector<double>& angles);
224
225 /// Set an entry in the table of drift speeds along E.
226 bool SetElectronVelocityE(const size_t ie, const size_t ib,
227 const size_t ia, const double v) {
228 return SetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
229 }
230 /// Get an entry in the table of drift speeds along E.
231 bool GetElectronVelocityE(const size_t ie, const size_t ib,
232 const size_t ia, double& v) {
233 return GetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
234 }
235 /// Set an entry in the table of drift speeds along ExB.
236 bool SetElectronVelocityExB(const size_t ie, const size_t ib,
237 const size_t ia, const double v) {
238 return SetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
239 }
240 /// Get an entry in the table of drift speeds along ExB.
241 bool GetElectronVelocityExB(const size_t ie, const size_t ib,
242 const size_t ia, double& v) {
243 return GetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
244 }
245 /// Set an entry in the table of drift speeds along Btrans.
246 bool SetElectronVelocityB(const size_t ie, const size_t ib,
247 const size_t ia, const double v) {
248 return SetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
249 }
250 /// Get an entry in the table of drift speeds along Btrans.
251 bool GetElectronVelocityB(const size_t ie, const size_t ib,
252 const size_t ia, double& v) {
253 return GetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
254 }
255 /// Set an entry in the table of longitudinal diffusion coefficients.
256 bool SetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
257 const size_t ia, const double dl) {
258 return SetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
259 }
260 /// Get an entry in the table of longitudinal diffusion coefficients.
261 bool GetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
262 const size_t ia, double& dl) {
263 return GetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
264 }
265 /// Set an entry in the table of transverse diffusion coefficients.
266 bool SetElectronTransverseDiffusion(const size_t ie, const size_t ib,
267 const size_t ia, const double dt) {
268 return SetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
269 }
270 /// Get an entry in the table of transverse diffusion coefficients.
271 bool GetElectronTransverseDiffusion(const size_t ie, const size_t ib,
272 const size_t ia, double& dt) {
273 return GetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
274 }
275 /// Set an entry in the table of Townsend coefficients.
276 bool SetElectronTownsend(const size_t ie, const size_t ib,
277 const size_t ia, const double alpha) {
278 return SetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
279 }
280 /// Get an entry in the table of Townsend coefficients.
281 bool GetElectronTownsend(const size_t ie, const size_t ib,
282 const size_t ia, double& alpha) {
283 return GetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
284 }
285 /// Set an entry in the table of attachment coefficients.
286 bool SetElectronAttachment(const size_t ie, const size_t ib,
287 const size_t ia, const double eta) {
288 return SetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
289 }
290 /// Get an entry in the table of attachment coefficients.
291 bool GetElectronAttachment(const size_t ie, const size_t ib,
292 const size_t ia, double& eta) {
293 return GetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
294 }
295
296 /// Set an entry in the table of Lorentz angles.
297 bool SetElectronLorentzAngle(const size_t ie, const size_t ib,
298 const size_t ia, const double lor) {
299 return SetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
300 }
301 /// Get an entry in the table of Lorentz angles.
302 bool GetElectronLorentzAngle(const size_t ie, const size_t ib,
303 const size_t ia, double& lor) {
304 return GetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
305 }
306
307 /// Set an entry in the table of drift speeds along E.
308 bool SetHoleVelocityE(const size_t ie, const size_t ib,
309 const size_t ia, const double v) {
310 return SetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
311 }
312 /// Get an entry in the table of drift speeds along E.
313 bool GetHoleVelocityE(const size_t ie, const size_t ib,
314 const size_t ia, double& v) {
315 return GetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
316 }
317 /// Set an entry in the table of drift speeds along ExB.
318 bool SetHoleVelocityExB(const size_t ie, const size_t ib,
319 const size_t ia, const double v) {
320 return SetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
321 }
322 /// Get an entry in the table of drift speeds along ExB.
323 bool GetHoleVelocityExB(const size_t ie, const size_t ib,
324 const size_t ia, double& v) {
325 return GetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
326 }
327 /// Set an entry in the table of drift speeds along Btrans.
328 bool SetHoleVelocityB(const size_t ie, const size_t ib,
329 const size_t ia, const double v) {
330 return SetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
331 }
332 /// Get an entry in the table of drift speeds along Btrans.
333 bool GetHoleVelocityB(const size_t ie, const size_t ib,
334 const size_t ia, double& v) {
335 return GetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
336 }
337
338 /// Set an entry in the table of longitudinal diffusion coefficients.
339 bool SetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
340 const size_t ia, const double dl) {
341 return SetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
342 }
343 /// Get an entry in the table of longitudinal diffusion coefficients.
344 bool GetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
345 const size_t ia, double& dl) {
346 return GetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
347 }
348 /// Set an entry in the table of transverse diffusion coefficients.
349 bool SetHoleTransverseDiffusion(const size_t ie, const size_t ib,
350 const size_t ia, const double dt) {
351 return SetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
352 }
353 /// Get an entry in the table of transverse diffusion coefficients.
354 bool GetHoleTransverseDiffusion(const size_t ie, const size_t ib,
355 const size_t ia, double& dt) {
356 return GetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
357 }
358 /// Set an entry in the table of Townsend coefficients.
359 bool SetHoleTownsend(const size_t ie, const size_t ib,
360 const size_t ia, const double alpha) {
361 return SetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
362 }
363 /// Get an entry in the table of Townsend coefficients.
364 bool GetHoleTownsend(const size_t ie, const size_t ib,
365 const size_t ia, double& alpha) {
366 return GetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
367 }
368 /// Set an entry in the table of attachment coefficients.
369 bool SetHoleAttachment(const size_t ie, const size_t ib,
370 const size_t ia, const double eta) {
371 return SetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
372 }
373 /// Get an entry in the table of attachment coefficients.
374 bool GetHoleAttachment(const size_t ie, const size_t ib,
375 const size_t ia, double& eta) {
376 return GetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
377 }
378
379 /// Initialise the table of ion mobilities from a list of
380 /// electric fields and corresponding mobilities.
381 /// The mobilities will be interpolated at the electric fields
382 /// of the currently set grid.
383 bool SetIonMobility(const std::vector<double>& fields,
384 const std::vector<double>& mobilities,
385 const bool negativeIons = false);
386 /// Set an entry in the table of ion mobilities.
387 bool SetIonMobility(const size_t ie, const size_t ib,
388 const size_t ia, const double mu);
389 /// Get an entry in the table of ion mobilities.
390 bool GetIonMobility(const size_t ie, const size_t ib,
391 const size_t ia, double& mu) {
392 return GetEntry(ie, ib, ia, "IonMobility", m_iMob, mu);
393 }
394
395 /// Set an entry in the table of longitudinal diffusion coefficients.
396 bool SetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
397 const size_t ia, const double dl) {
398 return SetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
399 }
400 /// Get an entry in the table of longitudinal diffusion coefficients.
401 bool GetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
402 const size_t ia, double& dl) {
403 return GetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
404 }
405 /// Set an entry in the table of transverse diffusion coefficients.
406 bool SetIonTransverseDiffusion(const size_t ie, const size_t ib,
407 const size_t ia, const double dt) {
408 return SetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
409 }
410 /// Get an entry in the table of transverse diffusion coefficients.
411 bool GetIonTransverseDiffusion(const size_t ie, const size_t ib,
412 const size_t ia, double& dt) {
413 return GetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
414 }
415 /// Set an entry in the table of dissociation coefficients.
416 bool SetIonDissociation(const size_t ie, const size_t ib,
417 const size_t ia, const double diss) {
418 return SetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
419 }
420 /// Get an entry in the table of dissociation coefficients.
421 bool GetIonDissociation(const size_t ie, const size_t ib,
422 const size_t ia, double& diss) {
423 return GetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
424 }
425
426 /// Set an entry in the table of negative ion mobilities.
427 bool SetNegativeIonMobility(const size_t ie, const size_t ib,
428 const size_t ia, const double mu) {
429 return SetEntry(ie, ib, ia, "NegativeIonMobility", m_nMob, mu);
430 }
431 /// Get an entry in the table of negative ion mobilities.
432 bool GetNegativeIonMobility(const size_t ie, const size_t ib,
433 const size_t ia, double& mu) {
434 return GetEntry(ie, ib, ia, "NegativeIonMobility", m_nMob, mu);
435 }
436
437 /// Reset all tables of transport parameters.
438 virtual void ResetTables();
439
441 m_eVelE.clear();
442 m_eVelB.clear();
443 m_eVelX.clear();
444 }
446 m_eDifL.clear();
447 m_eDifT.clear();
448 m_eDifM.clear();
449 }
450 void ResetElectronTownsend() { m_eAlp.clear(); }
451 void ResetElectronAttachment() { m_eAtt.clear(); }
453
455 m_hVelE.clear();
456 m_hVelB.clear();
457 m_hVelX.clear();
458 }
460 m_hDifL.clear();
461 m_hDifT.clear();
462 m_hDifM.clear();
463 }
464 void ResetHoleTownsend() { m_hAlp.clear(); }
465 void ResetHoleAttachment() { m_hAtt.clear(); }
466
467 void ResetIonMobility() { m_iMob.clear(); }
469 m_iDifL.clear();
470 m_iDifT.clear();
471 }
472 void ResetIonDissociation() { m_iDis.clear(); }
474
475 /// Select the extrapolation method for fields below/above the table range.
476 /// Possible options are "constant", "linear", and "exponential".
477 void SetExtrapolationMethodVelocity(const std::string& extrLow,
478 const std::string& extrHigh);
479 void SetExtrapolationMethodDiffusion(const std::string& extrLow,
480 const std::string& extrHigh);
481 void SetExtrapolationMethodTownsend(const std::string& extrLow,
482 const std::string& extrHigh);
483 void SetExtrapolationMethodAttachment(const std::string& extrLow,
484 const std::string& extrHigh);
485 void SetExtrapolationMethodIonMobility(const std::string& extrLow,
486 const std::string& extrHigh);
487 void SetExtrapolationMethodIonDissociation(const std::string& extrLow,
488 const std::string& extrHigh);
489
490 /// Set the degree of polynomial interpolation (usually 2).
491 void SetInterpolationMethodVelocity(const unsigned int intrp);
492 void SetInterpolationMethodDiffusion(const unsigned int intrp);
493 void SetInterpolationMethodTownsend(const unsigned int intrp);
494 void SetInterpolationMethodAttachment(const unsigned int intrp);
495 void SetInterpolationMethodIonMobility(const unsigned int intrp);
496 void SetInterpolationMethodIonDissociation(const unsigned int intrp);
497
498 // Scaling of fields and transport parameters.
499 virtual double ScaleElectricField(const double e) const { return e; }
500 virtual double UnScaleElectricField(const double e) const { return e; }
501 virtual double ScaleVelocity(const double v) const { return v; }
502 virtual double ScaleDiffusion(const double d) const { return d; }
503 virtual double ScaleDiffusionTensor(const double d) const { return d; }
504 virtual double ScaleTownsend(const double alpha) const { return alpha; }
505 virtual double ScaleAttachment(const double eta) const { return eta; }
506 virtual double ScaleLorentzAngle(const double lor) const { return lor; }
507 virtual double ScaleDissociation(const double diss) const { return diss; }
508
509 // Optical properties
510 /// Get the energy range [eV] of the available optical data.
511 virtual bool GetOpticalDataRange(double& emin, double& emax,
512 const unsigned int i = 0);
513 /// Get the complex dielectric function at a given energy.
514 virtual bool GetDielectricFunction(const double e, double& eps1, double& eps2,
515 const unsigned int i = 0);
516 // Get the photoabsorption cross-section [cm2] at a given energy.
517 virtual bool GetPhotoAbsorptionCrossSection(const double e, double& sigma,
518 const unsigned int i = 0);
519 virtual double GetPhotonCollisionRate(const double e);
520 virtual bool GetPhotonCollision(const double e, int& type, int& level,
521 double& e1, double& ctheta, int& nsec,
522 double& esec);
523
524 /// Switch on/off debugging messages
525 void EnableDebugging() { m_debug = true; }
526 void DisableDebugging() { m_debug = false; }
527
528 protected:
529 std::string m_className = "Medium";
530
531 static int m_idCounter;
532
533 // Id number
534 int m_id;
535 // Number of components
536 unsigned int m_nComponents = 1;
537 // Name
538 std::string m_name = "";
539 // Temperature [K]
540 double m_temperature = 293.15;
541 // Pressure [Torr]
542 double m_pressure = 760.;
543 // Static dielectric constant
544 double m_epsilon = 1.;
545 // (Effective) atomic number Z
546 double m_z = 1.;
547 // Atomic weight A
548 double m_a = 0.;
549 // Number density [cm-3]
550 double m_density = 0.;
551
552 // W value
553 double m_w = 0.;
554 // Fano factor
555 double m_fano = 0.;
556
557 // Transport flags
558 bool m_driftable = false;
559 bool m_microscopic = false;
560 bool m_ionisable = false;
561
562 // Update flag
563 bool m_isChanged = true;
564
565 // Switch on/off debugging messages
566 bool m_debug = false;
567
568 // Tables of transport parameters
569 bool m_tab2d = false;
570
571 // Field grids
572 std::vector<double> m_eFields;
573 std::vector<double> m_bFields;
574 std::vector<double> m_bAngles;
575
576 // Electrons
577 std::vector<std::vector<std::vector<double> > > m_eVelE;
578 std::vector<std::vector<std::vector<double> > > m_eVelX;
579 std::vector<std::vector<std::vector<double> > > m_eVelB;
580 std::vector<std::vector<std::vector<double> > > m_eDifL;
581 std::vector<std::vector<std::vector<double> > > m_eDifT;
582 std::vector<std::vector<std::vector<double> > > m_eAlp;
583 std::vector<std::vector<std::vector<double> > > m_eAtt;
584 std::vector<std::vector<std::vector<double> > > m_eLor;
585
586 std::vector<std::vector<std::vector<std::vector<double> > > > m_eDifM;
587
588 // Holes
589 std::vector<std::vector<std::vector<double> > > m_hVelE;
590 std::vector<std::vector<std::vector<double> > > m_hVelX;
591 std::vector<std::vector<std::vector<double> > > m_hVelB;
592 std::vector<std::vector<std::vector<double> > > m_hDifL;
593 std::vector<std::vector<std::vector<double> > > m_hDifT;
594 std::vector<std::vector<std::vector<double> > > m_hAlp;
595 std::vector<std::vector<std::vector<double> > > m_hAtt;
596
597 std::vector<std::vector<std::vector<std::vector<double> > > > m_hDifM;
598
599 // Ions
600 std::vector<std::vector<std::vector<double> > > m_iMob;
601 std::vector<std::vector<std::vector<double> > > m_iDifL;
602 std::vector<std::vector<std::vector<double> > > m_iDifT;
603 std::vector<std::vector<std::vector<double> > > m_iDis;
604 // Negative ions
605 std::vector<std::vector<std::vector<double> > > m_nMob;
606
607 // Thresholds for Townsend, attachment and dissociation coefficients.
608 unsigned int m_eThrAlp = 0;
609 unsigned int m_eThrAtt = 0;
610 unsigned int m_hThrAlp = 0;
611 unsigned int m_hThrAtt = 0;
612 unsigned int m_iThrDis = 0;
613
614 // Extrapolation methods (TODO: enum).
615 std::pair<unsigned int, unsigned int> m_extrVel = {0, 1};
616 std::pair<unsigned int, unsigned int> m_extrDif = {0, 1};
617 std::pair<unsigned int, unsigned int> m_extrAlp = {0, 1};
618 std::pair<unsigned int, unsigned int> m_extrAtt = {0, 1};
619 std::pair<unsigned int, unsigned int> m_extrLor = {0, 1};
620 std::pair<unsigned int, unsigned int> m_extrMob = {0, 1};
621 std::pair<unsigned int, unsigned int> m_extrDis = {0, 1};
622
623 // Interpolation methods
624 unsigned int m_intpVel = 2;
625 unsigned int m_intpDif = 2;
626 unsigned int m_intpAlp = 2;
627 unsigned int m_intpAtt = 2;
628 unsigned int m_intpLor = 2;
629 unsigned int m_intpMob = 2;
630 unsigned int m_intpDis = 2;
631
632 bool Velocity(const double ex, const double ey, const double ez,
633 const double bx, const double by, const double bz,
634 const std::vector<std::vector<std::vector<double> > >& velE,
635 const std::vector<std::vector<std::vector<double> > >& velB,
636 const std::vector<std::vector<std::vector<double> > >& velX,
637 const double q, double& vx, double& vy, double& vz) const;
638 static void Langevin(const double ex, const double ey, const double ez,
639 double bx, double by, double bz, const double mu,
640 double& vx, double& vy, double& vz);
641 static void Langevin(const double ex, const double ey, const double ez,
642 double bx, double by, double bz,
643 const double mu, const double muH,
644 double& vx, double& vy, double& vz);
645 bool Diffusion(const double ex, const double ey, const double ez,
646 const double bx, const double by, const double bz,
647 const std::vector<std::vector<std::vector<double> > >& difL,
648 const std::vector<std::vector<std::vector<double> > >& difT,
649 double& dl, double& dt) const;
650 bool Diffusion(const double ex, const double ey, const double ez,
651 const double bx, const double by, const double bz,
652 const std::vector<std::vector<std::vector<std::vector<double> > > >& diff,
653 double cov[3][3]) const;
654 bool Alpha(const double ex, const double ey, const double ez,
655 const double bx, const double by, const double bz,
656 const std::vector<std::vector<std::vector<double> > >& tab,
657 unsigned int intp, const unsigned int thr,
658 const std::pair<unsigned int, unsigned int>& extr,
659 double& alpha) const;
660 double GetAngle(const double ex, const double ey, const double ez,
661 const double bx, const double by, const double bz,
662 const double e, const double b) const;
663 bool Interpolate(const double e, const double b, const double a,
664 const std::vector<std::vector<std::vector<double> > >& table,
665 double& y, const unsigned int intp,
666 const std::pair<unsigned int, unsigned int>& extr) const;
667
668 double Interpolate1D(const double e, const std::vector<double>& table,
669 const std::vector<double>& fields,
670 const unsigned int intpMeth,
671 const std::pair<unsigned int, unsigned int>& extr) const;
672
673 bool SetEntry(const size_t i, const size_t j, const size_t k,
674 const std::string& fcn,
675 std::vector<std::vector<std::vector<double> > >& tab,
676 const double val);
677 bool GetEntry(const size_t i, const size_t j, const size_t k,
678 const std::string& fcn,
679 const std::vector<std::vector<std::vector<double> > >& tab,
680 double& val) const;
681
682 void SetExtrapolationMethod(const std::string& low, const std::string& high,
683 std::pair<unsigned int, unsigned int>& extr,
684 const std::string& fcn);
685 bool GetExtrapolationIndex(std::string str, unsigned int& nb) const;
686 size_t SetThreshold(
687 const std::vector<std::vector<std::vector<double> > >& tab) const;
688
689 void Clone(std::vector<std::vector<std::vector<double> > >& tab,
690 const std::vector<double>& efields,
691 const std::vector<double>& bfields,
692 const std::vector<double>& angles, const unsigned int intp,
693 const std::pair<unsigned int, unsigned int>& extr,
694 const double init, const std::string& label);
695 void Clone(
696 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
697 const size_t n, const std::vector<double>& efields,
698 const std::vector<double>& bfields, const std::vector<double>& angles,
699 const unsigned int intp,
700 const std::pair<unsigned int, unsigned int>& extr, const double init,
701 const std::string& label);
702
703 void Init(const size_t nE, const size_t nB, const size_t nA,
704 std::vector<std::vector<std::vector<double> > >& tab,
705 const double val);
706 void Init(
707 const size_t nE, const size_t nB, const size_t nA, const size_t nT,
708 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
709 const double val);
710
711};
712}
713
714#endif
virtual double IonMobility()
Low-field ion mobility [cm2 V-1 ns-1].
Definition Medium.cc:680
void ResetHoleAttachment()
Definition Medium.hh:465
virtual void GetComponent(const unsigned int i, std::string &label, double &f)
Get the name and fraction of a given component.
Definition Medium.cc:106
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.
Definition Medium.hh:427
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:1314
double m_density
Definition Medium.hh:550
void ResetIonMobility()
Definition Medium.hh:467
virtual bool GetPhotonCollision(const double e, int &type, int &level, double &e1, double &ctheta, int &nsec, double &esec)
Definition Medium.cc:780
double GetW() const
Get the W value.
Definition Medium.hh:86
virtual double UnScaleElectricField(const double e) const
Definition Medium.hh:500
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:344
virtual double ScaleVelocity(const double v) const
Definition Medium.hh:501
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:364
void ResetElectronAttachment()
Definition Medium.hh:451
virtual double GetMassDensity() const
Get the mass density [g/cm3].
Definition Medium.cc:102
std::vector< double > m_bFields
Definition Medium.hh:573
void ResetIonDiffusion()
Definition Medium.hh:468
virtual void GetElectronMomentum(const double e, double &px, double &py, double &pz, int &band)
Definition Medium.cc:538
bool IsMicroscopic() const
Does the medium have electron scattering rates?
Definition Medium.hh:79
void SetTemperature(const double t)
Set the temperature [K].
Definition Medium.cc:72
void ResetElectronLorentzAngle()
Definition Medium.hh:452
double m_pressure
Definition Medium.hh:542
void SetFanoFactor(const double f)
Set the Fano factor.
Definition Medium.hh:88
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:599
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:302
unsigned int m_intpMob
Definition Medium.hh:629
virtual double ScaleAttachment(const double eta) const
Definition Medium.hh:505
bool GetEntry(const size_t i, const size_t j, const size_t k, const std::string &fcn, const std::vector< std::vector< std::vector< double > > > &tab, double &val) const
Definition Medium.cc:982
void SetExtrapolationMethodIonDissociation(const std::string &extrLow, const std::string &extrHigh)
Definition Medium.cc:1210
void SetInterpolationMethodIonDissociation(const unsigned int intrp)
Definition Medium.cc:1296
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:256
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:1334
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:1300
virtual double ScaleDiffusion(const double d) const
Definition Medium.hh:502
virtual double ElectronMobility()
Low-field mobility [cm2 V-1 ns-1].
Definition Medium.cc:518
virtual double GetNumberDensity() const
Get the number density [cm-3].
Definition Medium.hh:63
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:374
virtual unsigned int GetNumberOfDeexcitationProducts() const
Definition Medium.hh:159
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:578
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:246
virtual double ScaleLorentzAngle(const double lor) const
Definition Medium.hh:506
unsigned int m_intpDis
Definition Medium.hh:630
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:494
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:411
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:334
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition Medium.hh:582
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:271
void ResetHoleVelocity()
Definition Medium.hh:454
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition Medium.hh:594
void ResetIonDissociation()
Definition Medium.hh:472
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:734
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:286
virtual void EnableDrift(const bool on=true)
Switch electron/ion/hole transport on/off.
Definition Medium.hh:70
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:308
double GetDielectricConstant() const
Get the relative static dielectric constant.
Definition Medium.hh:45
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:339
virtual void SetNumberDensity(const double n)
Set the number density [cm-3].
Definition Medium.cc:135
std::string m_name
Definition Medium.hh:538
double GetPressure() const
Definition Medium.hh:41
unsigned int m_intpDif
Definition Medium.hh:625
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:444
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].
Definition Medium.cc:630
virtual void EnablePrimaryIonisation(const bool on=true)
Make the medium ionisable or non-ionisable.
Definition Medium.hh:72
std::pair< unsigned int, unsigned int > m_extrAtt
Definition Medium.hh:618
int GetId() const
Return the id number of the class instance.
Definition Medium.hh:24
void SetExtrapolationMethodTownsend(const std::string &extrLow, const std::string &extrHigh)
Definition Medium.cc:1195
virtual void SetAtomicNumber(const double z)
Set the effective atomic number.
Definition Medium.cc:115
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:369
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:276
virtual double GetAtomicWeight() const
Get the effective atomic weight.
Definition Medium.hh:59
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:354
std::vector< std::vector< std::vector< double > > > m_iDifT
Definition Medium.hh:602
void SetExtrapolationMethodIonMobility(const std::string &extrLow, const std::string &extrHigh)
Definition Medium.cc:1205
void PlotVelocity(const std::string &carriers, TPad *pad)
Plot the drift velocity as function of the electric field.
Definition Medium.cc:161
unsigned int m_hThrAtt
Definition Medium.hh:611
std::pair< unsigned int, unsigned int > m_extrVel
Definition Medium.hh:615
virtual double ScaleElectricField(const double e) const
Definition Medium.hh:499
virtual double ScaleDiffusionTensor(const double d) const
Definition Medium.hh:503
unsigned int GetNumberOfComponents() const
Get number of components of the medium.
Definition Medium.hh:48
std::vector< std::vector< std::vector< double > > > m_hDifL
Definition Medium.hh:592
virtual double GetElectronNullCollisionRate(const int band=0)
Null-collision rate [ns-1].
Definition Medium.cc:545
virtual double GetAtomicNumber() const
Get the effective atomic number.
Definition Medium.hh:55
size_t SetThreshold(const std::vector< std::vector< std::vector< double > > > &tab) const
Definition Medium.cc:1252
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition Medium.hh:577
std::vector< std::vector< std::vector< double > > > m_eVelX
Definition Medium.hh:578
std::vector< std::vector< std::vector< double > > > m_eDifL
Definition Medium.hh:580
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:1408
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:323
virtual bool ElectronCollision(const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< std::pair< Particle, double > > &secondaries, int &ndxc, int &band)
Sample the collision type. Update energy and direction vector.
Definition Medium.cc:556
virtual double ScaleTownsend(const double alpha) const
Definition Medium.hh:504
unsigned int m_intpVel
Definition Medium.hh:624
virtual double GetPhotonCollisionRate(const double e)
Definition Medium.cc:773
void SetDielectricConstant(const double eps)
Set the relative static dielectric constant.
Definition Medium.cc:92
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:291
virtual bool GetDeexcitationProduct(const unsigned int i, double &t, double &s, int &type, double &energy) const
Definition Medium.cc:569
virtual double GetElectronCollisionRate(const double e, const int band=0)
Collision rate [ns-1] for given electron energy.
Definition Medium.cc:550
bool GetExtrapolationIndex(std::string str, unsigned int &nb) const
Definition Medium.cc:1235
void SetInterpolationMethodVelocity(const unsigned int intrp)
Set the degree of polynomial interpolation (usually 2).
Definition Medium.cc:1276
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:452
void SetInterpolationMethodDiffusion(const unsigned int intrp)
Definition Medium.cc:1280
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:328
double m_epsilon
Definition Medium.hh:544
void SetExtrapolationMethodDiffusion(const std::string &extrLow, const std::string &extrHigh)
Definition Medium.cc:1190
std::vector< double > m_eFields
Definition Medium.hh:572
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:231
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:958
static int m_idCounter
Definition Medium.hh:531
std::vector< std::vector< std::vector< double > > > m_hDifT
Definition Medium.hh:593
void DisableDebugging()
Definition Medium.hh:526
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:468
void ResetElectronDiffusion()
Definition Medium.hh:445
virtual bool IsSemiconductor() const
Is this medium a semiconductor?
Definition Medium.hh:30
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:791
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:722
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:297
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition Medium.hh:583
virtual void SetMassDensity(const double rho)
Set the mass density [g/cm3].
Definition Medium.cc:145
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:333
void PlotAlphaEta(const std::string &carriers, TPad *pad)
Plot Townsend and attachment coefficients.
Definition Medium.cc:189
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:281
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:236
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:396
virtual double ScaleDissociation(const double diss) const
Definition Medium.hh:507
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:413
const std::string & GetName() const
Get the medium name/identifier.
Definition Medium.hh:26
std::vector< double > m_bAngles
Definition Medium.hh:574
Medium()
Constructor.
Definition Medium.cc:61
std::vector< std::vector< std::vector< double > > > m_iDifL
Definition Medium.hh:601
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition Medium.hh:589
virtual bool IsGas() const
Is this medium a gas?
Definition Medium.hh:28
std::vector< std::vector< std::vector< double > > > m_eLor
Definition Medium.hh:584
std::vector< std::vector< std::vector< double > > > m_eDifT
Definition Medium.hh:581
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:523
void SetInterpolationMethodIonMobility(const unsigned int intrp)
Definition Medium.cc:1292
unsigned int m_intpAtt
Definition Medium.hh:627
void SetExtrapolationMethodVelocity(const std::string &extrLow, const std::string &extrHigh)
Definition Medium.cc:1185
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:390
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:481
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:241
virtual double HoleMobility()
Low-field mobility [cm2 V-1 ns-1].
Definition Medium.cc:625
virtual double NegativeIonMobility()
Low-field negative ion mobility [cm2 V-1 ns-1].
Definition Medium.cc:718
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:266
std::pair< unsigned int, unsigned int > m_extrDis
Definition Medium.hh:621
void PlotDiffusion(const std::string &carriers, TPad *pad)
Plot the diffusion coefficients as function of the electric field.
Definition Medium.cc:168
std::vector< std::vector< std::vector< double > > > m_hAtt
Definition Medium.hh:595
void SetPressure(const double p)
Definition Medium.cc:82
void SetExtrapolationMethod(const std::string &low, const std::string &high, std::pair< unsigned int, unsigned int > &extr, const std::string &fcn)
Definition Medium.cc:1215
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].
Definition Medium.cc:684
std::pair< unsigned int, unsigned int > m_extrAlp
Definition Medium.hh:617
void EnableDebugging()
Switch on/off debugging messages.
Definition Medium.hh:525
std::vector< std::vector< std::vector< double > > > m_eVelB
Definition Medium.hh:579
unsigned int m_eThrAtt
Definition Medium.hh:609
bool IsIonisable() const
Is charge deposition by charged particles/photon enabled in this medium?
Definition Medium.hh:81
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:421
void PlotAttachment(const std::string &carriers, TPad *pad)
Plot the attachment coefficient(s) as function of the electric field.
Definition Medium.cc:182
virtual void SetAtomicWeight(const double a)
Set the effective atomic weight.
Definition Medium.cc:125
bool SetIonMobility(const std::vector< double > &fields, const std::vector< double > &mobilities, const bool negativeIons=false)
Definition Medium.cc:1139
void SetExtrapolationMethodAttachment(const std::string &extrLow, const std::string &extrHigh)
Definition Medium.cc:1200
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:251
void ResetNegativeIonMobility()
Definition Medium.hh:473
std::pair< unsigned int, unsigned int > m_extrLor
Definition Medium.hh:619
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:401
unsigned int m_nComponents
Definition Medium.hh:536
virtual bool IsConductor() const
Is this medium a conductor?
Definition Medium.hh:32
unsigned int m_intpLor
Definition Medium.hh:628
std::string m_className
Definition Medium.hh:529
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:318
std::pair< unsigned int, unsigned int > m_extrDif
Definition Medium.hh:616
void SetInterpolationMethodAttachment(const unsigned int intrp)
Definition Medium.cc:1288
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:586
double GetTemperature() const
Get the temperature [K].
Definition Medium.hh:37
unsigned int m_iThrDis
Definition Medium.hh:612
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:667
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:432
virtual void ResetTables()
Reset all tables of transport parameters.
Definition Medium.cc:999
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:261
std::vector< std::vector< std::vector< double > > > m_hVelX
Definition Medium.hh:590
virtual ~Medium()
Destructor.
Definition Medium.cc:70
std::pair< unsigned int, unsigned int > m_extrMob
Definition Medium.hh:620
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:1018
void ResetHoleDiffusion()
Definition Medium.hh:459
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
Definition Medium.hh:586
std::vector< std::vector< std::vector< double > > > m_hVelB
Definition Medium.hh:591
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:196
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:359
std::vector< std::vector< std::vector< std::vector< double > > > > m_hDifM
Definition Medium.hh:597
unsigned int m_hThrAlp
Definition Medium.hh:610
void SetW(const double w)
Set the W value (average energy to produce an electron/ion or e/h pair).
Definition Medium.hh:84
std::vector< std::vector< std::vector< double > > > m_iMob
Definition Medium.hh:600
bool SetEntry(const size_t i, const size_t j, const size_t k, const std::string &fcn, std::vector< std::vector< std::vector< double > > > &tab, const double val)
Definition Medium.cc:966
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:416
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:349
void ResetElectronTownsend()
Definition Medium.hh:450
unsigned int m_eThrAlp
Definition Medium.hh:608
double GetFanoFactor() const
Get the Fano factor.
Definition Medium.hh:90
void PlotTownsend(const std::string &carriers, TPad *pad)
Plot the Townsend coefficient(s) as function of the electric field.
Definition Medium.cc:175
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:313
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:612
bool IsDriftable() const
Is charge carrier transport enabled in this medium?
Definition Medium.hh:77
void ResetHoleTownsend()
Definition Medium.hh:464
void ResetElectronVelocity()
Definition Medium.hh:440
void SetInterpolationMethodTownsend(const unsigned int intrp)
Definition Medium.cc:1284
virtual bool GetPhotoAbsorptionCrossSection(const double e, double &sigma, const unsigned int i=0)
Definition Medium.cc:752
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:226
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:660
double m_temperature
Definition Medium.hh:540
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:406
std::vector< std::vector< std::vector< double > > > m_iDis
Definition Medium.hh:603
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)
Definition Medium.cc:301
std::vector< std::vector< std::vector< double > > > m_nMob
Definition Medium.hh:605
unsigned int m_intpAlp
Definition Medium.hh:626