Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
DriftLineRKF.hh
Go to the documentation of this file.
1#ifndef G_DRIFTLINE_RKF_H
2#define G_DRIFTLINE_RKF_H
3
4#include <vector>
5#include <string>
6
7#include "Sensor.hh"
8#include "ViewDrift.hh"
9#include "Medium.hh"
10#include "GeometryBase.hh"
11
12namespace Garfield {
13
14/// Calculation of drift lines based on macroscopic transport coefficients
15/// using Runge-Kutta-Fehlberg integration.
16
18
19 public:
22
23 void SetSensor(Sensor* s);
24
25 void EnablePlotting(ViewDrift* view);
26 void DisablePlotting();
27
28 void SetIntegrationAccuracy(const double a);
29 void SetMaximumStepSize(const double ms);
30
31 void EnableStepSizeLimit() { m_useStepSizeLimit = true; }
32 void DisableStepSizeLimit() { m_useStepSizeLimit = false; }
33
34 void SetMaxSteps(const unsigned int m) { m_maxSteps = m; }
35 void EnableRejectKinks() { m_rejectKinks = true; }
36 void DisableRejectKinks() { m_rejectKinks = false; }
37
38 void SetElectronSignalScalingFactor(const double scale) { m_scaleElectronSignal = scale; }
39 void SetHoleSignalScalingFactor(const double scale) { m_scaleHoleSignal = scale; }
40 void SetIonSignalScalingFactor(const double scale) { m_scaleIonSignal = scale; }
41
42 bool DriftElectron(const double x0, const double y0, const double z0,
43 const double t0);
44 bool DriftHole(const double x0, const double y0, const double z0,
45 const double t0);
46 bool DriftIon(const double x0, const double y0, const double z0,
47 const double t0);
48
49 void GetEndPoint(double& x, double& y, double& z, double& t, int& st) const;
50 unsigned int GetNumberOfDriftLinePoints() const { return m_nPoints; }
51 void GetDriftLinePoint(const unsigned int i, double& x, double& y, double& z, double& t) const;
52
53 double GetArrivalTimeSpread();
54 double GetGain();
55 double GetDriftTime() const {
56 return m_nPoints > 0 ? m_path[m_nPoints - 1].t : 0.;
57 }
58
59 void EnableDebugging() { m_debug = true; }
60 void DisableDebugging() { m_debug = false; }
61
62 void EnableVerbose() { m_verbose = true; }
63 void DisableVerbose() { m_verbose = false; }
64
65 private:
66 static const unsigned int ParticleTypeElectron = 0;
67 static const unsigned int ParticleTypeIon = 1;
68 static const unsigned int ParticleTypeHole = 2;
69
70 std::string m_className;
71
72 Sensor* m_sensor;
73 Medium* m_medium;
74
75 unsigned int m_particleType;
76 double m_maxStepSize;
77 double m_accuracy;
78 unsigned int m_maxSteps;
79 unsigned int m_maxStepsToWire;
80 bool m_rejectKinks;
81 bool m_useStepSizeLimit;
82
83 bool m_usePlotting;
84 ViewDrift* m_view;
85
86 struct step {
87 // Position
88 double x;
89 double y;
90 double z;
91 // Time
92 double t;
93 // Integrated Townsend coefficient
94 double alphaint;
95 };
96 std::vector<step> m_path;
97 int m_status;
98 unsigned int m_nPoints;
99
100 double m_scaleElectronSignal;
101 double m_scaleHoleSignal;
102 double m_scaleIonSignal;
103
104 bool m_debug;
105 bool m_verbose;
106
107 // Calculate a drift line starting at a given position.
108 bool DriftLine(const double x0, const double y0, const double z0,
109 const double t0);
110 // Calculate transport parameters for the respective particle type.
111 bool GetVelocity(const double x, const double y, const double z,
112 double& vx, double& vy, double& vz, int& status);
113 bool GetVelocity(const double ex, const double ey, const double ez,
114 const double bx, const double by, const double bz,
115 double& vx, double& vy, double& vz) const;
116 bool GetDiffusion(const double ex, const double ey, const double ez,
117 const double bx, const double by, const double bz,
118 double& dl, double& dt) const;
119 bool GetTownsend(const double ex, const double ey, const double ez,
120 const double bx, const double by, const double bz,
121 double& alpha) const;
122 // Terminate a drift line at the edge of a boundary.
123 bool EndDriftLine();
124 // Drift a particle to a wire
125 bool DriftToWire(const double xw, const double yw, const double rw);
126 // Determine the longitudinal diffusion over the drift line.
127 double IntegrateDiffusion(const double x, const double y, const double z,
128 const double xe, const double ye,
129 const double ze);
130 // Determine the effective gain over the drift line.
131 double IntegrateTownsend(const double x, const double y, const double z,
132 const double xe, const double ye, const double ze,
133 const double tol);
134 // Calculate the signal for the current drift line.
135 void ComputeSignal(const double q, const double scale) const;
136};
137}
138
139#endif
void SetIonSignalScalingFactor(const double scale)
Definition: DriftLineRKF.hh:40
void SetIntegrationAccuracy(const double a)
Definition: DriftLineRKF.cc:44
bool DriftHole(const double x0, const double y0, const double z0, const double t0)
Definition: DriftLineRKF.cc:94
bool DriftIon(const double x0, const double y0, const double z0, const double t0)
void SetSensor(Sensor *s)
Definition: DriftLineRKF.cc:35
unsigned int GetNumberOfDriftLinePoints() const
Definition: DriftLineRKF.hh:50
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Definition: DriftLineRKF.cc:84
void GetDriftLinePoint(const unsigned int i, double &x, double &y, double &z, double &t) const
void EnablePlotting(ViewDrift *view)
Definition: DriftLineRKF.cc:68
double GetDriftTime() const
Definition: DriftLineRKF.hh:55
void SetMaxSteps(const unsigned int m)
Definition: DriftLineRKF.hh:34
void GetEndPoint(double &x, double &y, double &z, double &t, int &st) const
void SetHoleSignalScalingFactor(const double scale)
Definition: DriftLineRKF.hh:39
void SetMaximumStepSize(const double ms)
Definition: DriftLineRKF.cc:54
void SetElectronSignalScalingFactor(const double scale)
Definition: DriftLineRKF.hh:38
Abstract base class for media.
Definition: Medium.hh:11
Visualize drift lines and tracks.
Definition: ViewDrift.hh:15