1#ifndef G_DRIFTLINE_RKF_H
2#define G_DRIFTLINE_RKF_H
49 void RejectKinks(
const bool on =
true) { m_rejectKinks = on; }
69 bool DriftElectron(
const double x0,
const double y0,
const double z0,
72 bool DriftHole(
const double x0,
const double y0,
const double z0,
75 bool DriftIon(
const double x0,
const double y0,
const double z0,
79 bool DriftPositron(
const double x0,
const double y0,
const double z0,
89 void GetEndPoint(
double& x,
double& y,
double& z,
double& t,
int& st)
const;
100 double GetGain(
const double eps = 1.e-4);
102 double GetLoss(
const double eps = 1.e-4);
111 std::string m_className =
"DriftLineRKF";
114 Sensor* m_sensor =
nullptr;
116 enum class Particle { Electron = 0, Ion, Hole, Positron, NegativeIon };
118 Particle m_particle = Particle::Electron;
121 double m_maxStepSize = 0.;
123 double m_accuracy = 1.e-8;
125 bool m_rejectKinks =
true;
127 bool m_useStepSizeLimit =
false;
130 ViewDrift* m_view =
nullptr;
133 std::vector<std::array<double, 3> > m_x;
135 std::vector<double> m_t;
140 bool m_doSignal =
true;
142 unsigned int m_navg = 2;
145 double m_scaleE = 1.;
147 double m_scaleH = 1.;
149 double m_scaleI = 1.;
152 bool m_doAvalanche =
true;
153 enum class GainFluctuations { None = 0, Polya };
155 GainFluctuations m_gainFluctuations = GainFluctuations::None;
162 bool m_doIonTail =
false;
164 double m_nE = 0., m_nI = 0.;
166 bool m_debug =
false;
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,
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);
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;
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;
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);
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);
217 double IntegrateDiffusion(
const std::array<double, 3>& xi,
218 const std::array<double, 3>& xe,
219 const Particle particle,
const double tol);
221 double IntegrateAlpha(
const std::array<double, 3>& xi,
222 const std::array<double, 3>& xe,
223 const Particle particle,
const double tol);
225 double IntegrateEta(
const std::array<double, 3>& xi,
226 const std::array<double, 3>& xe,
227 const Particle particle,
const double tol);
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;
void SetIonSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by ions.
void SetSignalAveragingOrder(const unsigned int navg)
void EnableSignalCalculation(const bool on=true)
Switch on/off calculation of induced currents (default: enabled).
void DisablePlotting()
Switch off drift line plotting.
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.
void UnsetMaximumStepSize()
Do not apply an upper limit to the step size that is allowed.
void SetIntegrationAccuracy(const double eps)
Set the accuracy of the Runge Kutta Fehlberg drift line integration.
unsigned int GetNumberOfDriftLinePoints() const
Get the number of points of the most recent drift line.
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.
void RejectKinks(const bool on=true)
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.
void PrintDriftLine() const
Print the trajectory of the most recent drift line.
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
DriftLineRKF()
Constructor.
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.
void EnableAvalanche(const bool on=true)
Enable/disable simulation electron multiplication (default: on).
void SetMaximumStepSize(const double ms)
Set the maximum step size that is allowed. By default, there is no limit.
void EnableIonTail(const bool on=true)
Enable/disable simulation of the ion tail (default: off).
void SetElectronSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by electrons.
Visualize drift lines and tracks.