Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Sensor.hh
Go to the documentation of this file.
1#ifndef G_SENSOR_H
2#define G_SENSOR_H
3
4#include <fstream>
5#include <mutex>
6#include <tuple>
7#include <utility>
8#include <vector>
9
10#include "Component.hh"
11#include "Shaper.hh"
12
13namespace Garfield {
14
15/// %Sensor
16
17class Sensor {
18 public:
19 /// Constructor
20 Sensor() = default;
21 /// Destructor
23
24 /// Add a component.
25 void AddComponent(Component* comp);
26 /// Get the number of components attached to the sensor.
27 size_t GetNumberOfComponents() const { return m_components.size(); }
28 /// Retrieve the pointer to a given component.
29 Component* GetComponent(const unsigned int i);
30 /// Activate/deactivate a given component.
31 void EnableComponent(const unsigned int i, const bool on);
32 /// Activate/deactivate use of the magnetic field of a given component.
33 void EnableMagneticField(const unsigned int i, const bool on);
34
35 /// Add an electrode.
36 void AddElectrode(Component* comp, const std::string& label);
37 /// Get the number of electrodes attached to the sensor.
38 size_t GetNumberOfElectrodes() const { return m_electrodes.size(); }
39 /// Remove all components, electrodes and reset the sensor.
40 void Clear();
41
42 /// Get the drift field and potential at (x, y, z).
43 void ElectricField(const double x, const double y, const double z, double& ex,
44 double& ey, double& ez, double& v, Medium*& medium,
45 int& status);
46 /// Get the drift field at (x, y, z).
47 void ElectricField(const double x, const double y, const double z, double& ex,
48 double& ey, double& ez, Medium*& medium, int& status);
49
50 /// Get the magnetic field at (x, y, z).
51 void MagneticField(const double x, const double y, const double z, double& bx,
52 double& by, double& bz, int& status);
53
54 /// Get the weighting field at (x, y, z).
55 void WeightingField(const double x, const double y, const double z,
56 double& wx, double& wy, double& wz,
57 const std::string& label);
58 /// Get the weighting potential at (x, y, z).
59 double WeightingPotential(const double x, const double y, const double z,
60 const std::string& label);
61
62 /// Get the medium at (x, y, z).
63 bool GetMedium(const double x, const double y, const double z,
64 Medium*& medium);
65
66 /// Switch debugging messages on/off.
67 void EnableDebugging(const bool on = true) { m_debug = on; }
68
69 /// Set the user area to the default.
70 bool SetArea();
71 /// Set the user area explicitly.
72 bool SetArea(const double xmin, const double ymin, const double zmin,
73 const double xmax, const double ymax, const double zmax);
74 /// Return the current user area.
75 bool GetArea(double& xmin, double& ymin, double& zmin, double& xmax,
76 double& ymax, double& zmax);
77 /// Check if a point is inside the user area.
78 bool IsInArea(const double x, const double y, const double z);
79
80 /// Return the voltage range.
81 bool GetVoltageRange(double& vmin, double& vmax);
82
83 /// Start a new event, when computing the average signal over multiple events.
84 void NewSignal() { ++m_nEvents; }
85 /// Reset signals and induced charges of all electrodes.
86 void ClearSignal();
87
88 /** Set the time window and binning for the signal calculation.
89 * \param tstart start time [ns]
90 * \param tstep bin width [ns]
91 * \param nsteps number of bins
92 */
93 void SetTimeWindow(const double tstart, const double tstep,
94 const unsigned int nsteps);
95 /// Retrieve the time window and binning.
96 void GetTimeWindow(double& tstart, double& tstep,
97 unsigned int& nsteps) const {
98 tstart = m_tStart;
99 tstep = m_tStep;
100 nsteps = m_nTimeBins;
101 }
102
103 /// Compute the component of the signal due to the delayed weighting field.
104 void EnableDelayedSignal(const bool on = true) { m_delayedSignal = on; }
105 /// Set the points in time at which to evaluate the delayed weighting field.
106 void SetDelayedSignalTimes(const std::vector<double>& ts);
107 /// Set the number of points to be used when averaging the delayed
108 /// signal vector over a time bin (default: 0).
109 /// The averaging is done with a \f$2\times navg + 1\f$ point
110 /// Newton-Raphson integration.
111 void SetDelayedSignalAveragingOrder(const unsigned int navg) {
112 m_nAvgDelayedSignal = navg;
113 }
114
115 /// Set/override the signal in a given time bin explicitly.
116 void SetSignal(const std::string& label, const unsigned int bin,
117 const double signal);
118 /// Retrieve the total signal for a given electrode and time bin.
119 double GetSignal(const std::string& label, const unsigned int bin);
120 /// Retrieve the electron signal for a given electrode and time bin.
121 double GetSignal(const std::string& label, const unsigned int bin,
122 const int comp);
123 /// Retrieve the prompt signal for a given electrode and time bin.
124 double GetPromptSignal(const std::string& label, const unsigned int bin);
125 /// Retrieve the delayed signal for a given electrode and time bin.
126 double GetDelayedSignal(const std::string& label, const unsigned int bin);
127 /// Retrieve the electron signal for a given electrode and time bin.
128 double GetElectronSignal(const std::string& label, const unsigned int bin);
129 /// Retrieve the ion or hole signal for a given electrode and time bin.
130 double GetIonSignal(const std::string& label, const unsigned int bin);
131 /// Retrieve the delayed electron signal for a given electrode and time bin.
132 double GetDelayedElectronSignal(const std::string& label,
133 const unsigned int bin);
134 /// Retrieve the delayed ion/hole signal for a given electrode and time bin.
135 double GetDelayedIonSignal(const std::string& label, const unsigned int bin);
136 /// Calculated using the weighting potentials at the start and end points.
137 double GetInducedCharge(const std::string& label);
138
139 /// Set the function to be used for evaluating the transfer function.
140 void SetTransferFunction(double (*f)(double t));
141 /// Set the points to be used for interpolating the transfer function.
142 void SetTransferFunction(const std::vector<double>& times,
143 const std::vector<double>& values);
144 /// Set the transfer function using a Shaper object.
145 void SetTransferFunction(Shaper& shaper);
146 /// Evaluate the transfer function at a given time.
147 double GetTransferFunction(const double t);
148 /// Print some information about the presently set transfer function.
150 /// Cache integral and FFT of the transfer function
151 /// instead of recomputing it at every call (default: on).
152 void EnableTransferFunctionCache(const bool on = true) {
153 m_cacheTransferFunction = on;
154 }
155 /// Convolute the induced current on a given electrode
156 /// with the transfer function.
157 bool ConvoluteSignal(const std::string& label, const bool fft = false);
158 /// Convolute all induced currents with the transfer function.
159 bool ConvoluteSignals(const bool fft = false);
160 /// Replace the signal on a given electrode by its integral.
161 bool IntegrateSignal(const std::string& label);
162 /// Replace the signals on all electrodes curve by their integrals.
163 bool IntegrateSignals();
164 /// Return whether the signal has been integrated/convoluted.
165 bool IsIntegrated(const std::string& label) const;
166
167 /// Delay the signal and subtract an attenuated copy
168 /// (modelling a constant fraction discriminator).
169 /// \f[
170 /// v_{out} = v_{in}\left(t - t_d\right) - f v_{in}.
171 /// \f]
172 bool DelayAndSubtractFraction(const double td, const double f);
173 /// Set the function to be used for evaluating the noise component.
174 void SetNoiseFunction(double (*f)(double t));
175 /// Add noise to the induced signal.
176 void AddNoise(const bool total = true, const bool electron = false,
177 const bool ion = false);
178 /** Add white noise to the induced signal, given a desired output ENC.
179 * \param label name of the electrode
180 * \param enc Equivalent Noise Charge, in electrons.
181 * \param poisson flag to sample the number of noise pulses from a
182 * Poisson distribution. Otherwise the noise charge in each
183 * bin is sampled from a Gaussian distribution.
184 * \param q0 amplitude of the noise delta pulses (in electrons).
185 */
186 void AddWhiteNoise(const std::string& label, const double enc,
187 const bool poisson = true, const double q0 = 1.);
188 /// Add white noise to the induced signals on all electrodes.
189 void AddWhiteNoise(const double enc, const bool poisson = true,
190 const double q0 = 1.);
191
192 /** Determine the threshold crossings of the current signal curve.
193 * \param thr threshold value
194 * \param label electrode for which to compute the threshold crossings
195 * \param n number of threshold crossings
196 */
197 bool ComputeThresholdCrossings(const double thr, const std::string& label,
198 int& n);
199 /// Get the number of threshold crossings
200 /// (after having called ComputeThresholdCrossings).
202 return m_thresholdCrossings.size();
203 }
204 /** Retrieve the time and type of a given threshold crossing (after having
205 * called ComputeThresholdCrossings.
206 * \param i index
207 * \param time threshold crossing time [ns]
208 * \param level threshold (should correspond to the value requested).
209 * \param rise flag whether the crossing is on a rising or falling slope.
210 */
211 bool GetThresholdCrossing(const unsigned int i, double& time, double& level,
212 bool& rise) const;
213
214 /// Add the signal from a charge-carrier step.
215 void AddSignal(const double q, const double t0, const double t1,
216 const double x0, const double y0, const double z0,
217 const double x1, const double y1, const double z1,
218 const bool integrateWeightingField,
219 const bool useWeightingPotential = false);
220
221 /// Add the signal from a drift line.
222 void AddSignal(const double q, const std::vector<double>& ts,
223 const std::vector<std::array<double, 3> >& xs,
224 const std::vector<std::array<double, 3> >& vs,
225 const std::vector<double>& ns, const int navg,
226 const bool useWeightingPotential = false);
227
228 /// Exporting induced signal to a csv file.
229 void ExportSignal(const std::string& label, const std::string& filename);
230
231 /// Add the induced charge from a charge carrier drift between two points.
232 void AddInducedCharge(const double q, const double x0, const double y0,
233 const double z0, const double x1, const double y1,
234 const double z1);
235
236 /// Determine whether a line between two points has crossed a wire,
237 /// calls Component::IsWireCrossed.
238 bool IsWireCrossed(const double x0, const double y0, const double z0,
239 const double x1, const double y1, const double z1,
240 double& xc, double& yc, double& zc, const bool centre,
241 double& rc);
242 bool IsInTrapRadius(const double q0, const double x0, const double y0,
243 const double z0, double& xw, double& yw, double& rw);
244
245 /// Integrate the electric field flux through a line from
246 /// (x0,y0,z0) to (x1,y1,z1) along a direction (xp,yp,zp).
247 double IntegrateFluxLine(const double x0, const double y0, const double z0,
248 const double x1, const double y1, const double z1,
249 const double xp, const double yp, const double zp,
250 const unsigned int nI, const int isign = 0);
251 // TODO!
252 double GetTotalInducedCharge(const std::string label);
253 private:
254 std::string m_className = "Sensor";
255 /// Mutex.
256 std::mutex m_mutex;
257
258 /// Components
259 std::vector<std::tuple<Component*, bool, bool> > m_components;
260 Component* m_lastComponent = nullptr;
261
262 struct Electrode {
263 Component* comp;
264 std::string label;
265 std::vector<double> signal;
266 std::vector<double> delayedSignal;
267 std::vector<double> electronSignal;
268 std::vector<double> ionSignal;
269 std::vector<double> delayedElectronSignal;
270 std::vector<double> delayedIonSignal;
271 double charge;
272 bool integrated;
273 };
274 /// Electrodes
275 std::vector<Electrode> m_electrodes;
276
277 // Time window for signals
278 double m_tStart = 0.;
279 double m_tStep = 10.;
280 unsigned int m_nTimeBins = 200;
281 unsigned int m_nEvents = 0;
282 bool m_delayedSignal = false;
283 std::vector<double> m_delayedSignalTimes;
284 unsigned int m_nAvgDelayedSignal = 0;
285
286 // Transfer function
287 double (*m_fTransfer)(double t) = nullptr;
288 Shaper* m_shaper = nullptr;
289 std::vector<std::pair<double, double> > m_fTransferTab;
290 bool m_cacheTransferFunction = true;
291 // Integral of the transfer function squared.
292 double m_fTransferSq = -1.;
293 // FFT of the transfer function.
294 std::vector<double> m_fTransferFFT;
295
296 // Noise
297 double (*m_fNoise)(double t) = nullptr;
298
299 std::vector<std::pair<double, bool> > m_thresholdCrossings;
300 double m_thresholdLevel = 0.;
301
302 // User bounding box
303 double m_xMinUser = 0., m_yMinUser = 0., m_zMinUser = 0.;
304 double m_xMaxUser = 0., m_yMaxUser = 0., m_zMaxUser = 0.;
305 bool m_hasUserArea = false;
306
307 // Switch on/off debugging messages
308 bool m_debug = false;
309
310 // Return the current sensor size
311 bool GetBoundingBox(double& xmin, double& ymin, double& zmin, double& xmax,
312 double& ymax, double& zmax);
313
314 void FillSignal(Electrode& electrode, const double q,
315 const std::vector<double>& ts, const std::vector<double>& is,
316 const int navg, const bool delayed = false);
317 void FillBin(Electrode& electrode, const unsigned int bin,
318 const double signal, const bool electron, const bool delayed) {
319 std::lock_guard<std::mutex> guard(m_mutex);
320 electrode.signal[bin] += signal;
321 if (delayed) electrode.delayedSignal[bin] += signal;
322 if (electron) {
323 electrode.electronSignal[bin] += signal;
324 if (delayed) electrode.delayedElectronSignal[bin] += signal;
325 } else {
326 electrode.ionSignal[bin] += signal;
327 if (delayed) electrode.delayedIonSignal[bin] += signal;
328 }
329 }
330
331 void IntegrateSignal(Electrode& electrode);
332 void ConvoluteSignal(Electrode& electrode, const std::vector<double>& tab);
333 bool ConvoluteSignalFFT();
334 bool ConvoluteSignalFFT(const std::string& label);
335 void ConvoluteSignalFFT(Electrode& electrode, const std::vector<double>& tab,
336 const unsigned int nn);
337 // Evaluate the integral over the transfer function squared.
338 double TransferFunctionSq();
339 double InterpolateTransferFunctionTable(const double t) const;
340 void MakeTransferFunctionTable(std::vector<double>& tab);
341 void FFT(std::vector<double>& data, const bool inverse, const int nn);
342};
343} // namespace Garfield
344
345#endif
Abstract base class for components.
Definition: Component.hh:13
Abstract base class for media.
Definition: Medium.hh:13
void MagneticField(const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
Get the magnetic field at (x, y, z).
Definition: Sensor.cc:116
void EnableMagneticField(const unsigned int i, const bool on)
Activate/deactivate use of the magnetic field of a given component.
Definition: Sensor.cc:341
double GetElectronSignal(const std::string &label, const unsigned int bin)
Retrieve the electron signal for a given electrode and time bin.
Definition: Sensor.cc:857
bool ConvoluteSignal(const std::string &label, const bool fft=false)
Definition: Sensor.cc:1082
double GetTotalInducedCharge(const std::string label)
Definition: Sensor.cc:1700
bool GetVoltageRange(double &vmin, double &vmax)
Return the voltage range.
Definition: Sensor.cc:396
void SetDelayedSignalTimes(const std::vector< double > &ts)
Set the points in time at which to evaluate the delayed weighting field.
Definition: Sensor.cc:443
void SetTimeWindow(const double tstart, const double tstep, const unsigned int nsteps)
Definition: Sensor.cc:820
void ExportSignal(const std::string &label, const std::string &filename)
Exporting induced signal to a csv file.
Definition: Sensor.cc:1663
void SetSignal(const std::string &label, const unsigned int bin, const double signal)
Set/override the signal in a given time bin explicitly.
Definition: Sensor.cc:900
void EnableDelayedSignal(const bool on=true)
Compute the component of the signal due to the delayed weighting field.
Definition: Sensor.hh:104
double GetDelayedIonSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed ion/hole signal for a given electrode and time bin.
Definition: Sensor.cc:889
bool IsInTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
Definition: Sensor.cc:290
bool DelayAndSubtractFraction(const double td, const double f)
Definition: Sensor.cc:1283
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label)
Get the weighting field at (x, y, z).
Definition: Sensor.cc:131
double GetTransferFunction(const double t)
Evaluate the transfer function at a given time.
Definition: Sensor.cc:1041
void EnableComponent(const unsigned int i, const bool on)
Activate/deactivate a given component.
Definition: Sensor.cc:333
~Sensor()
Destructor.
Definition: Sensor.hh:22
void GetTimeWindow(double &tstart, double &tstep, unsigned int &nsteps) const
Retrieve the time window and binning.
Definition: Sensor.hh:96
void AddElectrode(Component *comp, const std::string &label)
Add an electrode.
Definition: Sensor.cc:349
bool IntegrateSignal(const std::string &label)
Replace the signal on a given electrode by its integral.
Definition: Sensor.cc:1244
void AddComponent(Component *comp)
Add a component.
Definition: Sensor.cc:316
void PrintTransferFunction()
Print some information about the presently set transfer function.
Definition: Sensor.cc:1050
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&medium, int &status)
Get the drift field and potential at (x, y, z).
Definition: Sensor.cc:65
bool SetArea()
Set the user area to the default.
Definition: Sensor.cc:183
double GetSignal(const std::string &label, const unsigned int bin)
Retrieve the total signal for a given electrode and time bin.
Definition: Sensor.cc:912
size_t GetNumberOfComponents() const
Get the number of components attached to the sensor.
Definition: Sensor.hh:27
void SetNoiseFunction(double(*f)(double t))
Set the function to be used for evaluating the noise component.
Definition: Sensor.cc:1301
double GetPromptSignal(const std::string &label, const unsigned int bin)
Retrieve the prompt signal for a given electrode and time bin.
Definition: Sensor.cc:946
bool IntegrateSignals()
Replace the signals on all electrodes curve by their integrals.
Definition: Sensor.cc:1234
bool ConvoluteSignals(const bool fft=false)
Convolute all induced currents with the transfer function.
Definition: Sensor.cc:1105
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
Definition: Sensor.cc:254
double GetInducedCharge(const std::string &label)
Calculated using the weighting potentials at the start and end points.
Definition: Sensor.cc:969
void AddNoise(const bool total=true, const bool electron=false, const bool ion=false)
Add noise to the induced signal.
Definition: Sensor.cc:1309
size_t GetNumberOfElectrodes() const
Get the number of electrodes attached to the sensor.
Definition: Sensor.hh:38
double IntegrateFluxLine(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const double xp, const double yp, const double zp, const unsigned int nI, const int isign=0)
Definition: Sensor.cc:301
double WeightingPotential(const double x, const double y, const double z, const std::string &label)
Get the weighting potential at (x, y, z).
Definition: Sensor.cc:147
Component * GetComponent(const unsigned int i)
Retrieve the pointer to a given component.
Definition: Sensor.cc:325
void AddWhiteNoise(const std::string &label, const double enc, const bool poisson=true, const double q0=1.)
Definition: Sensor.cc:1328
bool ComputeThresholdCrossings(const double thr, const std::string &label, int &n)
Definition: Sensor.cc:1439
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:159
void Clear()
Remove all components, electrodes and reset the sensor.
Definition: Sensor.cc:379
void AddSignal(const double q, const double t0, const double t1, const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const bool integrateWeightingField, const bool useWeightingPotential=false)
Add the signal from a charge-carrier step.
Definition: Sensor.cc:452
void SetTransferFunction(double(*f)(double t))
Set the function to be used for evaluating the transfer function.
Definition: Sensor.cc:979
bool GetThresholdCrossing(const unsigned int i, double &time, double &level, bool &rise) const
Definition: Sensor.cc:1553
size_t GetNumberOfThresholdCrossings() const
Definition: Sensor.hh:201
Sensor()=default
Constructor.
double GetDelayedElectronSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed electron signal for a given electrode and time bin.
Definition: Sensor.cc:878
void EnableTransferFunctionCache(const bool on=true)
Definition: Sensor.hh:152
bool IsWireCrossed(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc, const bool centre, double &rc)
Definition: Sensor.cc:276
void ClearSignal()
Reset signals and induced charges of all electrodes.
Definition: Sensor.cc:429
double GetDelayedSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed signal for a given electrode and time bin.
Definition: Sensor.cc:958
void AddInducedCharge(const double q, const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Add the induced charge from a charge carrier drift between two points.
Definition: Sensor.cc:798
double GetIonSignal(const std::string &label, const unsigned int bin)
Retrieve the ion or hole signal for a given electrode and time bin.
Definition: Sensor.cc:868
void EnableDebugging(const bool on=true)
Switch debugging messages on/off.
Definition: Sensor.hh:67
void SetDelayedSignalAveragingOrder(const unsigned int navg)
Definition: Sensor.hh:111
void NewSignal()
Start a new event, when computing the average signal over multiple events.
Definition: Sensor.hh:84
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
Definition: Sensor.cc:228
bool IsIntegrated(const std::string &label) const
Return whether the signal has been integrated/convoluted.
Definition: Sensor.cc:1276
Class for signal processing.
Definition: Shaper.hh:11