Garfield++ 3.0
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 <string>
5#include <vector>
6#include <array>
7
8#include "GeometryBase.hh"
9#include "Medium.hh"
10#include "Sensor.hh"
11#include "ViewDrift.hh"
12
13namespace Garfield {
14
15/// Calculation of drift lines based on macroscopic transport coefficients
16/// using Runge-Kutta-Fehlberg integration.
17
19 public:
20 /// Constructor
22 /// Destructor
24
25 /// Set the sensor.
26 void SetSensor(Sensor* s);
27
28 /// Switch on drift line plotting.
29 void EnablePlotting(ViewDrift* view);
30 /// Switch off drift line plotting.
31 void DisablePlotting();
32
33 /// Switch on/off calculation of induced currents (default: enabled).
34 void EnableSignalCalculation(const bool on = true) { m_doSignal = on; }
35 /// Set the number of points to be used when averaging the delayed
36 /// signal vector over a time bin in the Sensor class.
37 /// The averaging is done with a \f$2\times navg + 1\f$ point
38 /// Newton-Raphson integration. Default: 2.
39 void SetSignalAveragingOrder(const unsigned int navg) { m_navg = navg; }
40
41 /// Set the accuracy of the Runge Kutta Fehlberg drift line integration.
42 void SetIntegrationAccuracy(const double eps);
43 /// Set the maximum step size that is allowed. By default, there is no limit.
44 void SetMaximumStepSize(const double ms);
45 /// Do not apply an upper limit to the step size that is allowed.
46 void UnsetMaximumStepSize() { m_useStepSizeLimit = false; }
47 /// Request (or not) the drift line calculation to be aborted if the
48 /// drift line makes a bend sharper than 90 degrees.
49 void RejectKinks(const bool on = true) { m_rejectKinks = on; }
50
51 /// Set multiplication factor for the signal induced by electrons.
52 void SetElectronSignalScalingFactor(const double scale) { m_scaleE = scale; }
53 /// Set multiplication factor for the signal induced by holes.
54 void SetHoleSignalScalingFactor(const double scale) { m_scaleH = scale; }
55 /// Set multiplication factor for the signal induced by ions.
56 void SetIonSignalScalingFactor(const double scale) { m_scaleI = scale; }
57
58 /// Enable/disable simulation electron multiplication (default: on).
59 void EnableAvalanche(const bool on = true) { m_doAvalanche = on; }
60 /// Enable/disable simulation of the ion tail (default: off).
61 void EnableIonTail(const bool on = true) { m_doIonTail = on; }
62 /// Do not randomize the avalanche size.
63 void SetGainFluctuationsFixed(const double gain = -1.);
64 /// Sample the avalanche size from a Polya distribution with
65 /// shape parameter theta.
66 void SetGainFluctuationsPolya(const double theta, const double mean = -1.);
67
68 /// Simulate the drift line of an electron with a given starting point.
69 bool DriftElectron(const double x0, const double y0, const double z0,
70 const double t0);
71 /// Simulate the drift line of a hole with a given starting point.
72 bool DriftHole(const double x0, const double y0, const double z0,
73 const double t0);
74 /// Simulate the drift line of an ion with a given starting point.
75 bool DriftIon(const double x0, const double y0, const double z0,
76 const double t0);
77 /// Simulate the drift line of an electron with a given starting point,
78 /// assuming that it has positive charge.
79 bool DriftPositron(const double x0, const double y0, const double z0,
80 const double t0);
81 /// Simulate the drift line of an ion with a given starting point,
82 /// assuming that it has negative charge.
83 bool DriftNegativeIon(const double x0, const double y0, const double z0,
84 const double t0);
85
86 /// Print the trajectory of the most recent drift line.
87 void PrintDriftLine() const;
88 /// Get the end point and status flag of the most recent drift line.
89 void GetEndPoint(double& x, double& y, double& z, double& t, int& st) const;
90 /// Get the number of points of the most recent drift line.
91 unsigned int GetNumberOfDriftLinePoints() const { return m_x.size(); }
92 /// Get the coordinates and time of a point along the most recent drift line.
93 void GetDriftLinePoint(const unsigned int i, double& x, double& y, double& z,
94 double& t) const;
95
96 /// Compute the sigma of the arrival time distribution for the current
97 /// drift line by integrating the longitudinal diffusion coefficient.
98 double GetArrivalTimeSpread(const double eps = 1.e-4);
99 /// Compute the multiplication factor for the current drift line.
100 double GetGain(const double eps = 1.e-4);
101 /// Compute the attachment loss factor for the current drift line.
102 double GetLoss(const double eps = 1.e-4);
103 double GetDriftTime() const { return m_t.empty() ? 0. : m_t.back(); }
104
105 void GetAvalancheSize(double& ne, double& ni) const { ne = m_nE; ni = m_nI; }
106
107 void EnableDebugging() { m_debug = true; }
108 void DisableDebugging() { m_debug = false; }
109
110 private:
111 std::string m_className = "DriftLineRKF";
112
113 // Pointer to sensor.
114 Sensor* m_sensor = nullptr;
115
116 enum class Particle { Electron = 0, Ion, Hole, Positron, NegativeIon };
117 // Type of particle (of the most current drift line).
118 Particle m_particle = Particle::Electron;
119
120 // Maximum allowed step size.
121 double m_maxStepSize = 0.;
122 // Precision of the stepping algorithm.
123 double m_accuracy = 1.e-8;
124 // Flag to reject bends > 90 degrees or not.
125 bool m_rejectKinks = true;
126 // Flag to apply a cut on the maximum allowed step size or not.
127 bool m_useStepSizeLimit = false;
128
129 // Pointer to the drift viewer.
130 ViewDrift* m_view = nullptr;
131
132 // Points along the current drift line.
133 std::vector<std::array<double, 3> > m_x;
134 // Times corresponding to the points along the current drift line.
135 std::vector<double> m_t;
136 // Status flag of the current drift line.
137 int m_status = 0;
138
139 // Flag whether to calculate induced signals or not.
140 bool m_doSignal = true;
141 // Averaging order used when projecting the signal on the time bins.
142 unsigned int m_navg = 2;
143
144 // Scaling factor for electron signals.
145 double m_scaleE = 1.;
146 // Scaling factor for hole signals.
147 double m_scaleH = 1.;
148 // Scaling factor for ion signals.
149 double m_scaleI = 1.;
150
151 // Flag wether to simulate electron multiplication or not.
152 bool m_doAvalanche = true;
153 enum class GainFluctuations { None = 0, Polya };
154 // Model to be used for randomizing the avalanche size.
155 GainFluctuations m_gainFluctuations = GainFluctuations::None;
156 // Polya shape parameter.
157 double m_theta = 0.;
158 // Mean avalanche size (only used if < 1).
159 double m_gain = -1.;
160
161 // Flag whether to simulate the ion tail or not.
162 bool m_doIonTail = false;
163 // Avalanche size.
164 double m_nE = 0., m_nI = 0.;
165
166 bool m_debug = false;
167
168 // Calculate a drift line starting at a given position.
169 bool DriftLine(const double x0, const double y0, const double z0,
170 const double t0, const Particle particle,
171 std::vector<double>& ts,
172 std::vector<std::array<double, 3> >& xs, int& status);
173 bool Avalanche(const Particle particle,
174 const std::vector<std::array<double, 3> >& xs,
175 std::vector<double>& ne, std::vector<double>& ni,
176 double& scale);
177 bool AddIonTail(const std::vector<double>& te,
178 const std::vector<std::array<double, 3> >& xe,
179 const std::vector<double>& ni, const double scale);
180 // Compute electric and magnetic field at a given position.
181 int GetField(const std::array<double, 3>& x,
182 double& ex, double& ey, double& ez,
183 double& bx, double& by, double& bz,
184 Medium*& medium) const;
185 // Calculate transport parameters for the respective particle type.
186 bool GetVelocity(const std::array<double, 3>& x, const Particle particle,
187 std::array<double, 3>& v, int& status) const;
188 bool GetVelocity(const double ex, const double ey, const double ez,
189 const double bx, const double by, const double bz,
190 Medium* medium, const Particle particle,
191 std::array<double, 3>& v) const;
192 bool GetDiffusion(const double ex, const double ey, const double ez,
193 const double bx, const double by, const double bz,
194 Medium* medium, const Particle particle,
195 double& dl, double& dt) const;
196 bool GetAlpha(const double ex, const double ey, const double ez,
197 const double bx, const double by, const double bz,
198 Medium* medium, const Particle particle, double& alpha) const;
199 bool GetEta(const double ex, const double ey, const double ez,
200 const double bx, const double by, const double bz,
201 Medium* medium, const Particle particle, double& eta) const;
202
203 // Terminate a drift line at the edge of a boundary.
204 bool Terminate(const std::array<double, 3>& xx0,
205 const std::array<double, 3>& xx1,
206 const Particle particle,
207 std::vector<double>& ts,
208 std::vector<std::array<double, 3> >& xs);
209
210 // Drift a particle to a wire
211 bool DriftToWire(const double xw, const double yw, const double rw,
212 const Particle particle,
213 std::vector<double>& ts,
214 std::vector<std::array<double, 3> >& xs, int& stat);
215
216 // Integrate the longitudinal diffusion over a step.
217 double IntegrateDiffusion(const std::array<double, 3>& xi,
218 const std::array<double, 3>& xe,
219 const Particle particle, const double tol);
220 // Integrate the Townsend coefficient over a step.
221 double IntegrateAlpha(const std::array<double, 3>& xi,
222 const std::array<double, 3>& xe,
223 const Particle particle, const double tol);
224 // Integrate the attachment coefficient over a step.
225 double IntegrateEta(const std::array<double, 3>& xi,
226 const std::array<double, 3>& xe,
227 const Particle particle, const double tol);
228
229 // Calculate the signal for the current drift line.
230 void ComputeSignal(const Particle particle, const double scale,
231 const std::vector<double>& ts,
232 const std::vector<std::array<double, 3> >& xs,
233 const std::vector<double>& ne) const;
234};
235}
236
237#endif
void SetIonSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by ions.
Definition: DriftLineRKF.hh:56
void SetSignalAveragingOrder(const unsigned int navg)
Definition: DriftLineRKF.hh:39
void EnableSignalCalculation(const bool on=true)
Switch on/off calculation of induced currents (default: enabled).
Definition: DriftLineRKF.hh:34
void DisablePlotting()
Switch off drift line plotting.
Definition: DriftLineRKF.cc:85
bool DriftNegativeIon(const double x0, const double y0, const double z0, const double t0)
double GetLoss(const double eps=1.e-4)
Compute the attachment loss factor for the current drift line.
bool DriftHole(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of a hole with a given starting point.
bool DriftIon(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an ion with a given starting point.
void SetSensor(Sensor *s)
Set the sensor.
Definition: DriftLineRKF.cc:50
void UnsetMaximumStepSize()
Do not apply an upper limit to the step size that is allowed.
Definition: DriftLineRKF.hh:46
void SetIntegrationAccuracy(const double eps)
Set the accuracy of the Runge Kutta Fehlberg drift line integration.
Definition: DriftLineRKF.cc:58
unsigned int GetNumberOfDriftLinePoints() const
Get the number of points of the most recent drift line.
Definition: DriftLineRKF.hh:91
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an electron with a given starting point.
void SetGainFluctuationsPolya(const double theta, const double mean=-1.)
double GetGain(const double eps=1.e-4)
Compute the multiplication factor for the current drift line.
void SetGainFluctuationsFixed(const double gain=-1.)
Do not randomize the avalanche size.
Definition: DriftLineRKF.cc:87
void RejectKinks(const bool on=true)
Definition: DriftLineRKF.hh:49
bool DriftPositron(const double x0, const double y0, const double z0, const double t0)
void GetDriftLinePoint(const unsigned int i, double &x, double &y, double &z, double &t) const
Get the coordinates and time of a point along the most recent drift line.
~DriftLineRKF()
Destructor.
Definition: DriftLineRKF.hh:23
void PrintDriftLine() const
Print the trajectory of the most recent drift line.
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
Definition: DriftLineRKF.cc:77
DriftLineRKF()
Constructor.
Definition: DriftLineRKF.cc:45
double GetDriftTime() const
double GetArrivalTimeSpread(const double eps=1.e-4)
void GetAvalancheSize(double &ne, double &ni) const
void GetEndPoint(double &x, double &y, double &z, double &t, int &st) const
Get the end point and status flag of the most recent drift line.
void SetHoleSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by holes.
Definition: DriftLineRKF.hh:54
void EnableAvalanche(const bool on=true)
Enable/disable simulation electron multiplication (default: on).
Definition: DriftLineRKF.hh:59
void SetMaximumStepSize(const double ms)
Set the maximum step size that is allowed. By default, there is no limit.
Definition: DriftLineRKF.cc:67
void EnableIonTail(const bool on=true)
Enable/disable simulation of the ion tail (default: off).
Definition: DriftLineRKF.hh:61
void SetElectronSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by electrons.
Definition: DriftLineRKF.hh:52
Visualize drift lines and tracks.
Definition: ViewDrift.hh:19