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// Include this header if we're compiling with the GPU or this is the first time
2// without
3#if defined(__GPUCOMPILE__) || !defined(G_MEDIUM_H)
4
5#if !defined(__GPUCOMPILE__) && !defined(G_MEDIUM_H)
6#define G_MEDIUM_H
7#endif
8
9#ifdef __GPUCOMPILE__
10#include "GPUInterface.hh"
11#else
12
13#include <string>
14#include <vector>
15
17#endif
18
21
22class TPad;
23
24namespace Garfield {
25
26// setup class names depending on if this is compiling the GPU static version or
27// not
28#ifdef __GPUCOMPILE__
29#else
30class MediumGPU;
31#endif
32
34
35class GARFIELD_CLASS_NAME(Medium) {
36 public:
37#ifdef __GPUCOMPILE__
39 GARFIELD_CLASS_NAME(Medium)() = default;
41 ~GARFIELD_CLASS_NAME(Medium)() {};
42#else
44 GARFIELD_CLASS_NAME(Medium)();
46 virtual ~GARFIELD_CLASS_NAME(Medium)();
47#endif
48
50 __DEVICE__ int GetId() const { return m_id; }
51
52#ifndef __GPUCOMPILE__
54 const std::string& GetName() const { return m_name; }
56 virtual bool IsGas() const { return false; }
58 virtual bool IsSemiconductor() const { return false; }
60 virtual bool IsConductor() const { return false; }
61
63 void SetTemperature(const double t);
65 double GetTemperature() const { return m_temperature; }
66 // Set the pressure [Torr].
67 void SetPressure(const double p);
68 // Get the pressure [Torr].
69 double GetPressure() const { return m_pressure; }
71 void SetDielectricConstant(const double eps);
73 double GetDielectricConstant() const { return m_epsilon; }
74
76 unsigned int GetNumberOfComponents() const { return m_nComponents; }
78 virtual void GetComponent(const unsigned int i, std::string& label,
79 double& f);
81 virtual void SetAtomicNumber(const double z);
83 virtual double GetAtomicNumber() const { return m_z; }
85 virtual void SetAtomicWeight(const double a);
87 virtual double GetAtomicWeight() const { return m_a; }
89 virtual void SetNumberDensity(const double n);
91 virtual double GetNumberDensity() const { return m_density; }
93 virtual void SetMassDensity(const double rho);
95 virtual double GetMassDensity() const;
96
98 virtual void EnableDrift(const bool on = true) { m_driftable = on; }
100 virtual void EnablePrimaryIonisation(const bool on = true) {
101 m_ionisable = on;
102 }
103#endif
104
106 __DEVICE__ bool IsDriftable() const { return m_driftable; }
108 __DEVICE__ bool IsMicroscopic() const { return m_microscopic; }
109
110#ifndef __GPUCOMPILE__
112 bool IsIonisable() const { return m_ionisable; }
113
115 void SetW(const double w) { m_w = w; }
117 double GetW() const { return m_w; }
119 void SetFanoFactor(const double f) { m_fano = f; }
121 double GetFanoFactor() const { return m_fano; }
122
124 void PlotVelocity(const std::string& carriers, TPad* pad);
126 void PlotDiffusion(const std::string& carriers, TPad* pad);
128 void PlotTownsend(const std::string& carriers, TPad* pad);
130 void PlotAttachment(const std::string& carriers, TPad* pad);
132 void PlotAlphaEta(const std::string& carriers, TPad* pad);
133
134 // Transport parameters for electrons
136 virtual bool ElectronVelocity(const double ex, const double ey,
137 const double ez, const double bx,
138 const double by, const double bz, double& vx,
139 double& vy, double& vz);
142 virtual bool ElectronVelocityFluxBulk(const double ex, const double ey,
143 const double ez, const double bx,
144 const double by, const double bz,
145 double& wv, double& wr);
147 virtual bool ElectronDiffusion(const double ex, const double ey,
148 const double ez, const double bx,
149 const double by, const double bz, double& dl,
150 double& dt);
154 virtual bool ElectronDiffusion(const double ex, const double ey,
155 const double ez, const double bx,
156 const double by, const double bz,
157 double cov[3][3]);
159 virtual bool ElectronTownsend(const double ex, const double ey,
160 const double ez, const double bx,
161 const double by, const double bz,
162 double& alpha);
164 virtual bool ElectronAttachment(const double ex, const double ey,
165 const double ez, const double bx,
166 const double by, const double bz,
167 double& eta);
169 virtual bool ElectronTOFIonisation(const double ex, const double ey,
170 const double ez, const double bx,
171 const double by, const double bz,
172 double& riontof);
174 virtual bool ElectronTOFAttachment(const double ex, const double ey,
175 const double ez, const double bx,
176 const double by, const double bz,
177 double& ratttof);
179 virtual bool ElectronLorentzAngle(const double ex, const double ey,
180 const double ez, const double bx,
181 const double by, const double bz,
182 double& lor);
184 virtual double ElectronMobility();
185 // Microscopic electron transport properties
186
188 virtual double GetElectronEnergy(const double px, const double py,
189 const double pz, double& vx, double& vy,
190 double& vz, const int band = 0);
193 virtual void GetElectronMomentum(const double e, double& px, double& py,
194 double& pz, int& band);
195
197 virtual double GetElectronNullCollisionRate(const int band = 0);
198#endif
199
200#ifdef __GPUCOMPILE__
201
202 __device__ cuda_t GetElectronCollisionRate(const cuda_t e, const int band);
203
204 __device__ bool ElectronCollision(const cuda_t e, int& type, int& level,
205 cuda_t& e1, cuda_t& dx, cuda_t& dy,
206 cuda_t& dz, Particle* secondaries_type,
207 cuda_t* secondaries_energy,
208 int& num_secondaries, int& ndxc, int& band);
209#else
211 virtual double GetElectronCollisionRate(const double e, const int band = 0);
212 struct Secondary {
213 Particle type = Particle::Electron;
214 double energy = 0.;
215 double time = 0.;
216 double distance = 0.;
217 };
219 virtual bool ElectronCollision(const double e, int& type, int& level,
220 double& e1, double& dx, double& dy, double& dz,
221 std::vector<Secondary>& secondaries,
222 int& band);
223#endif
224
225#ifndef __GPUCOMPILE__
226 // Transport parameters for holes
228 virtual bool HoleVelocity(const double ex, const double ey, const double ez,
229 const double bx, const double by, const double bz,
230 double& vx, double& vy, double& vz);
232 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
233 const double bx, const double by, const double bz,
234 double& dl, double& dt);
236 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
237 const double bx, const double by, const double bz,
238 double cov[3][3]);
240 virtual bool HoleTownsend(const double ex, const double ey, const double ez,
241 const double bx, const double by, const double bz,
242 double& alpha);
244 virtual bool HoleAttachment(const double ex, const double ey, const double ez,
245 const double bx, const double by, const double bz,
246 double& eta);
248 virtual double HoleMobility();
249
250 // Transport parameters for ions
252 virtual bool IonVelocity(const double ex, const double ey, const double ez,
253 const double bx, const double by, const double bz,
254 double& vx, double& vy, double& vz);
255 bool HasIonVelocity() const { return !(m_iVel.empty() && m_iMob.empty()); }
257 virtual bool IonDiffusion(const double ex, const double ey, const double ez,
258 const double bx, const double by, const double bz,
259 double& dl, double& dt);
261 virtual bool IonDissociation(const double ex, const double ey,
262 const double ez, const double bx,
263 const double by, const double bz, double& diss);
265 virtual double IonMobility();
266
268 virtual bool NegativeIonVelocity(const double ex, const double ey,
269 const double ez, const double bx,
270 const double by, const double bz, double& vx,
271 double& vy, double& vz);
273 virtual double NegativeIonMobility();
274
276 void SetFieldGrid(double emin, double emax, const size_t ne, bool logE,
277 double bmin = 0., double bmax = 0., const size_t nb = 1,
278 double amin = HalfPi, double amax = HalfPi,
279 const size_t na = 1);
281 void SetFieldGrid(const std::vector<double>& efields,
282 const std::vector<double>& bfields,
283 const std::vector<double>& angles);
285 void GetFieldGrid(std::vector<double>& efields, std::vector<double>& bfields,
286 std::vector<double>& angles);
287
289 bool SetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia,
290 const double v) {
291 return SetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
292 }
294 bool GetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia,
295 double& v) {
296 return GetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
297 }
299 bool SetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia,
300 const double v) {
301 return SetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
302 }
304 bool GetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia,
305 double& v) {
306 return GetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
307 }
309 bool SetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia,
310 const double v) {
311 return SetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
312 }
314 bool GetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia,
315 double& v) {
316 return GetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
317 }
319 bool SetElectronFluxVelocity(const size_t ie, const size_t ib,
320 const size_t ia, const double v) {
321 return SetEntry(ie, ib, ia, "ElectronFluxVelocity", m_eVelWv, v);
322 }
324 bool GetElectronFluxVelocity(const size_t ie, const size_t ib,
325 const size_t ia, double& v) {
326 return GetEntry(ie, ib, ia, "ElectronFluxVelocity", m_eVelWv, v);
327 }
329 bool SetElectronBulkVelocity(const size_t ie, const size_t ib,
330 const size_t ia, const double v) {
331 return SetEntry(ie, ib, ia, "ElectronBulkVelocity", m_eVelWr, v);
332 }
334 bool GetElectronBulkVelocity(const size_t ie, const size_t ib,
335 const size_t ia, double& v) {
336 return GetEntry(ie, ib, ia, "ElectronBulkVelocity", m_eVelWr, v);
337 }
339 bool SetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
340 const size_t ia, const double dl) {
341 return SetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
342 }
344 bool GetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
345 const size_t ia, double& dl) {
346 return GetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
347 }
349 bool SetElectronTransverseDiffusion(const size_t ie, const size_t ib,
350 const size_t ia, const double dt) {
351 return SetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
352 }
354 bool GetElectronTransverseDiffusion(const size_t ie, const size_t ib,
355 const size_t ia, double& dt) {
356 return GetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
357 }
359 bool SetElectronTownsend(const size_t ie, const size_t ib, const size_t ia,
360 const double alpha) {
361 return SetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
362 }
364 bool GetElectronTownsend(const size_t ie, const size_t ib, const size_t ia,
365 double& alpha) {
366 return GetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
367 }
369 bool SetElectronAttachment(const size_t ie, const size_t ib, const size_t ia,
370 const double eta) {
371 return SetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
372 }
374 bool GetElectronAttachment(const size_t ie, const size_t ib, const size_t ia,
375 double& eta) {
376 return GetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
377 }
379 bool SetElectronTOFIonisation(const size_t ie, const size_t ib,
380 const size_t ia, const double v) {
381 return SetEntry(ie, ib, ia, "ElectronTOFIonisation", m_eRIon, v);
382 }
384 bool GetElectronTOFIonisation(const size_t ie, const size_t ib,
385 const size_t ia, double& v) {
386 return GetEntry(ie, ib, ia, "ElectronTOFIonisation", m_eRIon, v);
387 }
389 bool SetElectronTOFAttachment(const size_t ie, const size_t ib,
390 const size_t ia, const double v) {
391 return SetEntry(ie, ib, ia, "ElectronTOFAttachment", m_eRAtt, v);
392 }
394 bool GetElectronTOFAttachment(const size_t ie, const size_t ib,
395 const size_t ia, double& v) {
396 return GetEntry(ie, ib, ia, "ElectronTOFAttachment", m_eRAtt, v);
397 }
398
400 bool SetElectronLorentzAngle(const size_t ie, const size_t ib,
401 const size_t ia, const double lor) {
402 return SetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
403 }
405 bool GetElectronLorentzAngle(const size_t ie, const size_t ib,
406 const size_t ia, double& lor) {
407 return GetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
408 }
409
411 bool SetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia,
412 const double v) {
413 return SetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
414 }
416 bool GetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia,
417 double& v) {
418 return GetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
419 }
421 bool SetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia,
422 const double v) {
423 return SetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
424 }
426 bool GetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia,
427 double& v) {
428 return GetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
429 }
431 bool SetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia,
432 const double v) {
433 return SetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
434 }
436 bool GetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia,
437 double& v) {
438 return GetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
439 }
440
442 bool SetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
443 const size_t ia, const double dl) {
444 return SetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
445 }
447 bool GetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
448 const size_t ia, double& dl) {
449 return GetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
450 }
452 bool SetHoleTransverseDiffusion(const size_t ie, const size_t ib,
453 const size_t ia, const double dt) {
454 return SetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
455 }
457 bool GetHoleTransverseDiffusion(const size_t ie, const size_t ib,
458 const size_t ia, double& dt) {
459 return GetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
460 }
462 bool SetHoleTownsend(const size_t ie, const size_t ib, const size_t ia,
463 const double alpha) {
464 return SetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
465 }
467 bool GetHoleTownsend(const size_t ie, const size_t ib, const size_t ia,
468 double& alpha) {
469 return GetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
470 }
472 bool SetHoleAttachment(const size_t ie, const size_t ib, const size_t ia,
473 const double eta) {
474 return SetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
475 }
477 bool GetHoleAttachment(const size_t ie, const size_t ib, const size_t ia,
478 double& eta) {
479 return GetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
480 }
481
486 bool SetIonMobility(const std::vector<double>& fields,
487 const std::vector<double>& mobilities,
488 const bool negativeIons = false);
490 bool SetIonMobility(const size_t ie, const size_t ib, const size_t ia,
491 const double mu);
493 bool GetIonMobility(const size_t ie, const size_t ib, const size_t ia,
494 double& mu) {
495 return GetEntry(ie, ib, ia, "IonMobility", m_iMob, mu);
496 }
497
499 bool SetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
500 const size_t ia, const double dl) {
501 return SetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
502 }
504 bool GetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
505 const size_t ia, double& dl) {
506 return GetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
507 }
509 bool SetIonTransverseDiffusion(const size_t ie, const size_t ib,
510 const size_t ia, const double dt) {
511 return SetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
512 }
514 bool GetIonTransverseDiffusion(const size_t ie, const size_t ib,
515 const size_t ia, double& dt) {
516 return GetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
517 }
519 bool SetIonDissociation(const size_t ie, const size_t ib, const size_t ia,
520 const double diss) {
521 return SetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
522 }
524 bool GetIonDissociation(const size_t ie, const size_t ib, const size_t ia,
525 double& diss) {
526 return GetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
527 }
528
530 bool SetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia,
531 const double mu);
533 bool GetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia,
534 double& mu) {
535 return GetEntry(ie, ib, ia, "NegativeIonMobility", m_nMob, mu);
536 }
537
539 virtual void ResetTables();
540
541 void ResetElectronVelocity() {
542 m_eVelE.clear();
543 m_eVelB.clear();
544 m_eVelX.clear();
545 m_eVelWv.clear();
546 m_eVelWr.clear();
547 }
548 void ResetElectronDiffusion() {
549 m_eDifL.clear();
550 m_eDifT.clear();
551 m_eDifM.clear();
552 }
553 void ResetElectronTownsend() { m_eAlp.clear(); }
554 void ResetElectronAttachment() { m_eAtt.clear(); }
555 void ResetElectronTOFRates() {
556 m_eRIon.clear();
557 m_eRAtt.clear();
558 }
559 void ResetElectronLorentzAngle() { m_eLor.clear(); }
560
561 void ResetHoleVelocity() {
562 m_hVelE.clear();
563 m_hVelB.clear();
564 m_hVelX.clear();
565 }
566 void ResetHoleDiffusion() {
567 m_hDifL.clear();
568 m_hDifT.clear();
569 m_hDifM.clear();
570 }
571 void ResetHoleTownsend() { m_hAlp.clear(); }
572 void ResetHoleAttachment() { m_hAtt.clear(); }
573
574 void ResetIonMobility() {
575 m_iMob.clear();
576 m_iVel.clear();
577 }
578 void ResetIonDiffusion() {
579 m_iDifL.clear();
580 m_iDifT.clear();
581 }
582 void ResetIonDissociation() { m_iDis.clear(); }
583 void ResetNegativeIonMobility() {
584 m_nMob.clear();
585 m_nVel.clear();
586 }
587
588 void VelocityFromMobility(
589 const std::vector<std::vector<std::vector<double> > >& mob,
590 std::vector<std::vector<std::vector<double> > >& vel);
591
594 void SetExtrapolationMethodVelocity(const std::string& extrLow,
595 const std::string& extrHigh);
596 void SetExtrapolationMethodDiffusion(const std::string& extrLow,
597 const std::string& extrHigh);
598 void SetExtrapolationMethodTownsend(const std::string& extrLow,
599 const std::string& extrHigh);
600 void SetExtrapolationMethodAttachment(const std::string& extrLow,
601 const std::string& extrHigh);
602 void SetExtrapolationMethodIonMobility(const std::string& extrLow,
603 const std::string& extrHigh);
604 void SetExtrapolationMethodIonDissociation(const std::string& extrLow,
605 const std::string& extrHigh);
606
608 void SetInterpolationMethodVelocity(const unsigned int intrp);
609 void SetInterpolationMethodDiffusion(const unsigned int intrp);
610 void SetInterpolationMethodTownsend(const unsigned int intrp);
611 void SetInterpolationMethodAttachment(const unsigned int intrp);
612 void SetInterpolationMethodIonMobility(const unsigned int intrp);
613 void SetInterpolationMethodIonDissociation(const unsigned int intrp);
614
615 // Scaling of fields and transport parameters.
616 virtual double ScaleElectricField(const double e) const { return e; }
617 virtual double UnScaleElectricField(const double e) const { return e; }
618 virtual double ScaleVelocity(const double v) const { return v; }
619 virtual double ScaleDiffusion(const double d) const { return d; }
620 virtual double ScaleDiffusionTensor(const double d) const { return d; }
621 virtual double ScaleTownsend(const double alpha) const { return alpha; }
622 virtual double ScaleAttachment(const double eta) const { return eta; }
623 virtual double ScaleLorentzAngle(const double lor) const { return lor; }
624 virtual double ScaleDissociation(const double diss) const { return diss; }
625
626 // Optical properties
628 virtual bool GetOpticalDataRange(double& emin, double& emax,
629 const unsigned int i = 0);
631 virtual bool GetDielectricFunction(const double e, double& eps1, double& eps2,
632 const unsigned int i = 0);
633 // Get the photoabsorption cross-section [cm2] at a given energy.
634 virtual bool GetPhotoAbsorptionCrossSection(const double e, double& sigma,
635 const unsigned int i = 0);
636 virtual double GetPhotonCollisionRate(const double e);
637 virtual bool PhotonCollision(const double e, int& type, int& level,
638 double& e1, double& ctheta,
639 std::vector<Secondary>& secondaries);
640
642 void EnableDebugging() { m_debug = true; }
643 void DisableDebugging() { m_debug = false; }
644
646 virtual double CreateGPUTransferObject(MediumGPU*& med_gpu);
647
648 protected:
649 std::string m_className = "Medium";
650
651 static int m_idCounter;
652
653 // Number of components
654 unsigned int m_nComponents = 1;
655 // Name
656 std::string m_name = "";
657 // Temperature [K]
658 double m_temperature = 293.15;
659 // Pressure [Torr]
660 double m_pressure = 760.;
661 // Static dielectric constant
662 double m_epsilon = 1.;
663 // (Effective) atomic number Z
664 double m_z = 1.;
665 // Atomic weight A
666 double m_a = 0.;
667 // Number density [cm-3]
668 double m_density = 0.;
669#endif
670
671 // Id number
672 int m_id;
673
674 // Transport flags
675 bool m_driftable = false;
676 bool m_microscopic = false;
677 bool m_ionisable = false;
678
679#ifndef __GPUCOMPILE__
680
681 // W value
682 double m_w = 0.;
683 // Fano factor
684 double m_fano = 0.;
685
686 // Update flag
687 bool m_isChanged = true;
688
689 // Switch on/off debugging messages
690 bool m_debug = false;
691
692 // Tables of transport parameters
693 bool m_tab2d = false;
694
695 // Field grids
696 std::vector<double> m_eFields;
697 std::vector<double> m_bFields;
698 std::vector<double> m_bAngles;
699
700 // Electrons
701 std::vector<std::vector<std::vector<double> > > m_eVelE;
702 std::vector<std::vector<std::vector<double> > > m_eVelX;
703 std::vector<std::vector<std::vector<double> > > m_eVelB;
704 std::vector<std::vector<std::vector<double> > > m_eDifL;
705 std::vector<std::vector<std::vector<double> > > m_eDifT;
706 std::vector<std::vector<std::vector<double> > > m_eAlp;
707 std::vector<std::vector<std::vector<double> > > m_eAtt;
708 std::vector<std::vector<std::vector<double> > > m_eLor;
709 std::vector<std::vector<std::vector<double> > > m_eVelWv;
710 std::vector<std::vector<std::vector<double> > > m_eVelWr;
711 std::vector<std::vector<std::vector<double> > > m_eRIon;
712 std::vector<std::vector<std::vector<double> > > m_eRAtt;
713
714 std::vector<std::vector<std::vector<std::vector<double> > > > m_eDifM;
715
716 // Holes
717 std::vector<std::vector<std::vector<double> > > m_hVelE;
718 std::vector<std::vector<std::vector<double> > > m_hVelX;
719 std::vector<std::vector<std::vector<double> > > m_hVelB;
720 std::vector<std::vector<std::vector<double> > > m_hDifL;
721 std::vector<std::vector<std::vector<double> > > m_hDifT;
722 std::vector<std::vector<std::vector<double> > > m_hAlp;
723 std::vector<std::vector<std::vector<double> > > m_hAtt;
724
725 std::vector<std::vector<std::vector<std::vector<double> > > > m_hDifM;
726
727 // Ions
728 std::vector<std::vector<std::vector<double> > > m_iMob;
729 std::vector<std::vector<std::vector<double> > > m_iVel;
730 std::vector<std::vector<std::vector<double> > > m_iDifL;
731 std::vector<std::vector<std::vector<double> > > m_iDifT;
732 std::vector<std::vector<std::vector<double> > > m_iDis;
733 // Negative ions
734 std::vector<std::vector<std::vector<double> > > m_nMob;
735 std::vector<std::vector<std::vector<double> > > m_nVel;
736
737 // Thresholds for Townsend, attachment and dissociation coefficients.
738 unsigned int m_eThrAlp = 0;
739 unsigned int m_eThrAtt = 0;
740 unsigned int m_hThrAlp = 0;
741 unsigned int m_hThrAtt = 0;
742 unsigned int m_iThrDis = 0;
743
744 // Extrapolation methods (TODO: enum).
745 std::pair<unsigned int, unsigned int> m_extrVel = {0, 1};
746 std::pair<unsigned int, unsigned int> m_extrDif = {0, 1};
747 std::pair<unsigned int, unsigned int> m_extrAlp = {0, 1};
748 std::pair<unsigned int, unsigned int> m_extrAtt = {0, 1};
749 std::pair<unsigned int, unsigned int> m_extrLor = {0, 1};
750 std::pair<unsigned int, unsigned int> m_extrMob = {0, 1};
751 std::pair<unsigned int, unsigned int> m_extrDis = {0, 1};
752
753 // Interpolation methods
754 unsigned int m_intpVel = 2;
755 unsigned int m_intpDif = 2;
756 unsigned int m_intpAlp = 2;
757 unsigned int m_intpAtt = 2;
758 unsigned int m_intpLor = 2;
759 unsigned int m_intpMob = 2;
760 unsigned int m_intpDis = 2;
761
762 bool Velocity(const double ex, const double ey, const double ez,
763 const double bx, const double by, const double bz,
764 const std::vector<std::vector<std::vector<double> > >& velE,
765 const std::vector<std::vector<std::vector<double> > >& velB,
766 const std::vector<std::vector<std::vector<double> > >& velX,
767 const double q, double& vx, double& vy, double& vz) const;
768 bool VelocityFluxBulk(
769 const double ex, const double ey, const double ez, const double bx,
770 const double by, const double bz,
771 const std::vector<std::vector<std::vector<double> > >& velWv,
772 const std::vector<std::vector<std::vector<double> > >& velWr,
773 double& wv, double& wr) const;
774 static void Langevin(const double ex, const double ey, const double ez,
775 double bx, double by, double bz, const double mu,
776 double& vx, double& vy, double& vz);
777 static void Langevin(const double ex, const double ey, const double ez,
778 double bx, double by, double bz, const double mu,
779 const double muH, double& vx, double& vy, double& vz);
780 bool Diffusion(const double ex, const double ey, const double ez,
781 const double bx, const double by, const double bz,
782 const std::vector<std::vector<std::vector<double> > >& difL,
783 const std::vector<std::vector<std::vector<double> > >& difT,
784 double& dl, double& dt) const;
785 bool Diffusion(
786 const double ex, const double ey, const double ez, const double bx,
787 const double by, const double bz,
788 const std::vector<std::vector<std::vector<std::vector<double> > > >& diff,
789 double cov[3][3]) const;
790 bool Alpha(const double ex, const double ey, const double ez, const double bx,
791 const double by, const double bz,
792 const std::vector<std::vector<std::vector<double> > >& tab,
793 unsigned int intp, const unsigned int thr,
794 const std::pair<unsigned int, unsigned int>& extr,
795 double& alpha) const;
796 double GetAngle(const double ex, const double ey, const double ez,
797 const double bx, const double by, const double bz,
798 const double e, const double b) const;
799 bool Interpolate(const double e, const double b, const double a,
800 const std::vector<std::vector<std::vector<double> > >& table,
801 double& y, const unsigned int intp,
802 const std::pair<unsigned int, unsigned int>& extr,
803 const bool logval = false) const;
804
805 double Interpolate1D(const double e, const std::vector<double>& table,
806 const std::vector<double>& fields,
807 const unsigned int intpMeth,
808 const std::pair<unsigned int, unsigned int>& extr,
809 const bool logval = false) const;
810
811 bool SetEntry(const size_t i, const size_t j, const size_t k,
812 const std::string& fcn,
813 std::vector<std::vector<std::vector<double> > >& tab,
814 const double val);
815 bool GetEntry(const size_t i, const size_t j, const size_t k,
816 const std::string& fcn,
817 const std::vector<std::vector<std::vector<double> > >& tab,
818 double& val) const;
819
820 void SetExtrapolationMethod(const std::string& low, const std::string& high,
821 std::pair<unsigned int, unsigned int>& extr,
822 const std::string& fcn);
823 bool GetExtrapolationIndex(std::string str, unsigned int& nb) const;
824 size_t SetThreshold(
825 const std::vector<std::vector<std::vector<double> > >& tab) const;
826
827 void Clone(std::vector<std::vector<std::vector<double> > >& tab,
828 const std::vector<double>& efields,
829 const std::vector<double>& bfields,
830 const std::vector<double>& angles, const unsigned int intp,
831 const std::pair<unsigned int, unsigned int>& extr,
832 const double init, const std::string& label);
833 void Clone(std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
834 const size_t n, const std::vector<double>& efields,
835 const std::vector<double>& bfields,
836 const std::vector<double>& angles, const unsigned int intp,
837 const std::pair<unsigned int, unsigned int>& extr,
838 const double init, const std::string& label);
839
840 void Init(const size_t nE, const size_t nB, const size_t nA,
841 std::vector<std::vector<std::vector<double> > >& tab,
842 const double val);
843 void Init(const size_t nE, const size_t nB, const size_t nA, const size_t nT,
844 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
845 const double val);
846
847#else
848
849#include "MediumMagboltz.hh"
850
851 friend class MediumGas;
852 friend class MediumMagboltz;
853
854 // enum to mimic polymorphism
855 enum class MediumType { Medium = 0, MediumGas, MediumMagboltz };
856
857 MediumType m_MediumType{MediumType::Medium};
858
859#endif
860};
861} // namespace Garfield
862
863#endif
#define GARFIELD_CLASS_NAME(name)
#define __DEVICE__
#define __device__
Definition Vector.hh:12
class GARFIELD_CLASS_NAME(Component)
Abstract base class for components.
Definition Component.hh:31