Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization 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 <array>
5#include <functional>
6#include <mutex>
7#include <string>
8#include <tuple>
9#include <utility>
10#include <vector>
11
12class TPad;
13
14namespace Garfield {
15
16class SensorGPU;
17class Shaper;
18class Component;
19class Medium;
20
22class Sensor {
23 public:
25 Sensor() = default;
27 ~Sensor() = default;
28
34 size_t GetNumberOfComponents() const { return m_components.size(); }
36 Component* GetComponent(const unsigned int i);
38 void EnableComponent(const unsigned int i, const bool on);
40 void EnableMagneticField(const unsigned int i, const bool on);
42 bool HasMagneticField() const;
43
45 void AddElectrode(Component* comp, const std::string& label);
47 size_t GetNumberOfElectrodes() const { return m_electrodes.size(); }
51 void Clear();
52
54 void ElectricField(const double x, const double y, const double z, double& ex,
55 double& ey, double& ez, double& v, Medium*& medium,
56 int& status);
58 void ElectricField(const double x, const double y, const double z, double& ex,
59 double& ey, double& ez, Medium*& medium, int& status);
60
62 void MagneticField(const double x, const double y, const double z, double& bx,
63 double& by, double& bz, int& status);
64
66 void WeightingField(const double x, const double y, const double z,
67 double& wx, double& wy, double& wz,
68 const std::string& label);
70 double WeightingPotential(const double x, const double y, const double z,
71 const std::string& label);
72
74 double DelayedWeightingPotential(const double x, const double y,
75 const double z, const double t,
76 const std::string& label);
77
79 Medium* GetMedium(const double x, const double y, const double z);
80
82 void EnableDebugging(const bool on = true) { m_debug = on; }
83
85 bool SetArea(const bool verbose = false);
87 bool SetArea(const double xmin, const double ymin, const double zmin,
88 const double xmax, const double ymax, const double zmax);
90 bool GetArea(double& xmin, double& ymin, double& zmin, double& xmax,
91 double& ymax, double& zmax);
93 bool IsInArea(const double x, const double y, const double z);
94
96 bool IsInside(const double x, const double y, const double z);
97
99 bool GetVoltageRange(double& vmin, double& vmax);
100
102 void NewSignal() { ++m_nEvents; }
105
111 void SetTimeWindow(const double tstart, const double tstep,
112 const unsigned int nsteps);
114 void GetTimeWindow(double& tstart, double& tstep,
115 unsigned int& nsteps) const {
116 tstart = m_tStart;
117 tstep = m_tStep;
118 nsteps = m_nTimeBins;
119 }
120
122 void EnableDelayedSignal(const bool on = true) { m_delayedSignal = on; }
124 void SetDelayedSignalTimes(const std::vector<double>& ts);
129 void SetDelayedSignalAveragingOrder(const unsigned int navg) {
130 m_nAvgDelayedSignal = navg;
131 }
132
134 void SetSignal(const std::string& label, const unsigned int bin,
135 const double signal);
137 void SetSignal(const std::string& label, const std::vector<double>& ts,
138 const std::vector<double>& is);
140 double GetSignal(const std::string& label, const unsigned int bin);
142 double GetSignal(const std::string& label, const unsigned int bin,
143 const int comp);
145 double GetPromptSignal(const std::string& label, const unsigned int bin);
147 double GetDelayedSignal(const std::string& label, const unsigned int bin);
149 double GetElectronSignal(const std::string& label, const unsigned int bin);
151 double GetIonSignal(const std::string& label, const unsigned int bin);
153 double GetDelayedElectronSignal(const std::string& label,
154 const unsigned int bin);
156 double GetDelayedIonSignal(const std::string& label, const unsigned int bin);
158 double GetInducedCharge(const std::string& label);
159
161 void SetTransferFunction(std::function<double(double)>);
163 void SetTransferFunction(const std::vector<double>& times,
164 const std::vector<double>& values);
168 double GetTransferFunction(const double t);
175 void EnableTransferFunctionCache(const bool on = true) {
177 }
178
180 bool ConvoluteSignal(const std::string& label, const bool fft = false);
182 bool ConvoluteSignals(const bool fft = false);
184 bool IntegrateSignal(const std::string& label);
188 bool IsIntegrated(const std::string& label) const;
189
195 bool DelayAndSubtractFraction(const double td, const double f);
197 void SetNoiseFunction(double (*f)(double t));
199 void AddNoise(const bool total = true, const bool electron = false,
200 const bool ion = false);
209 void AddWhiteNoise(const std::string& label, const double enc,
210 const bool poisson = true, const double q0 = 1.);
212 void AddWhiteNoise(const double enc, const bool poisson = true,
213 const double q0 = 1.);
214
220 bool ComputeThresholdCrossings(const double thr, const std::string& label,
221 int& n);
225 return m_thresholdCrossings.size();
226 }
227
234 bool GetThresholdCrossing(const unsigned int i, double& time, double& level,
235 bool& rise) const;
236
240 const double q, const std::vector<double>& ts,
241 const std::vector<std::array<double, 3> >& xs);
245 const double q, const std::vector<double>& ts,
246 const std::vector<std::array<double, 3> >& xs,
247 const std::vector<double>& qs);
250 void AddSignalWeightingField(const double q, const std::vector<double>& ts,
251 const std::vector<std::array<double, 3> >& xs,
252 const bool integrateWeightingField);
256 void AddSignalWeightingField(const double q, const std::vector<double>& ts,
257 const std::vector<std::array<double, 3> >& xs,
258 const std::vector<std::array<double, 3> >& vs,
259 const std::vector<double>& ns, const int navg);
260
262 void PlotSignal(const std::string& label, TPad* pad,
263 const std::string optTotal = "t",
264 const std::string optPrompt = "",
265 const std::string optDelayed = "");
267 void ExportSignal(const std::string& label, const std::string& filename,
268 const bool chargeCariers = false) const;
269
271 void AddInducedCharge(const double q, const double x0, const double y0,
272 const double z0, const double x1, const double y1,
273 const double z1);
274
277 bool CrossedWire(const double x0, const double y0, const double z0,
278 const double x1, const double y1, const double z1,
279 double& xc, double& yc, double& zc, const bool centre,
280 double& rc);
282 bool InTrapRadius(const double q0, const double x0, const double y0,
283 const double z0, double& xw, double& yw, double& rw);
286 bool CrossedPlane(const double x0, const double y0, const double z0,
287 const double x1, const double y1, const double z1,
288 double& xc, double& yc, double& zc);
289
292 double IntegrateFluxLine(const double x0, const double y0, const double z0,
293 const double x1, const double y1, const double z1,
294 const double xp, const double yp, const double zp,
295 const unsigned int nI, const int isign = 0);
296 // TODO!
297 double GetTotalInducedCharge(const std::string& label);
298
299 double StepSizeHint();
300
302 double CreateGPUTransferObject(SensorGPU*& sensor_gpu);
303#if USEGPU
305 void TransferGPUElectrodeSignals(SensorGPU*& sensor_gpu);
306#endif
307
308 private:
309 std::string m_className = "Sensor";
311 std::mutex m_mutex;
312
314 std::vector<std::tuple<Component*, bool, bool> > m_components;
315
316 struct Electrode {
318 std::string label;
319 std::vector<double> signal;
320 std::vector<double> delayedSignal;
321 std::vector<double> electronSignal;
322 std::vector<double> ionSignal;
323 std::vector<double> delayedElectronSignal;
324 std::vector<double> delayedIonSignal;
325 double charge;
327 };
328
329 std::vector<Electrode> m_electrodes;
330
331 // Time window for signals
332 double m_tStart = 0.;
333 double m_tStep = 10.;
334 unsigned int m_nTimeBins = 200;
335 unsigned int m_nEvents = 0;
336
337 bool m_delayedSignal = false;
338 std::vector<double> m_delayedSignalTimes;
339 unsigned int m_nAvgDelayedSignal = 0;
340
341 // Transfer function
342 std::function<double(double)> m_fTransfer;
343 Shaper* m_shaper = nullptr;
344 std::vector<std::pair<double, double> > m_fTransferTab;
346 // Integral of the transfer function squared.
347 double m_fTransferSq = -1.;
348 // FFT of the transfer function.
349 std::vector<double> m_fTransferFFT;
350
351 // Noise
352 double (*m_fNoise)(double t) = nullptr;
353
354 std::vector<std::pair<double, bool> > m_thresholdCrossings;
355 double m_thresholdLevel = 0.;
356
357 // User bounding box
358 double m_xMinUser = 0., m_yMinUser = 0., m_zMinUser = 0.;
359 double m_xMaxUser = 0., m_yMaxUser = 0., m_zMaxUser = 0.;
360 bool m_hasUserArea = false;
361
362 // Switch on/off debugging messages
363 bool m_debug = false;
364
365 // Return the current sensor size
366 bool GetBoundingBox(double& xmin, double& ymin, double& zmin, double& xmax,
367 double& ymax, double& zmax);
368
369 void FillSignal(Electrode& electrode, const double q,
370 const std::vector<double>& ts, const std::vector<double>& is,
371 const int navg, const bool delayed = false);
372 void FillBin(Electrode& electrode, const unsigned int bin,
373 const double signal, const bool electron, const bool delayed) {
374 std::lock_guard<std::mutex> guard(m_mutex);
375 electrode.signal[bin] += signal;
376 if (delayed) electrode.delayedSignal[bin] += signal;
377 if (electron) {
378 electrode.electronSignal[bin] += signal;
379 if (delayed) electrode.delayedElectronSignal[bin] += signal;
380 } else {
381 electrode.ionSignal[bin] += signal;
382 if (delayed) electrode.delayedIonSignal[bin] += signal;
383 }
384 }
385
386 void IntegrateSignal(Electrode& electrode);
387 void ConvoluteSignal(Electrode& electrode, const std::vector<double>& tab);
389 bool ConvoluteSignalFFT(const std::string& label);
390 void ConvoluteSignalFFT(Electrode& electrode, const std::vector<double>& tab,
391 const unsigned int nn);
392 // Evaluate the integral over the transfer function squared.
394 double InterpolateTransferFunctionTable(const double t) const;
395 void MakeTransferFunctionTable(std::vector<double>& tab);
396 void FFT(std::vector<double>& data, const bool inverse, const int nn);
397};
398
399} // namespace Garfield
400#endif
Abstract base class for components.
Definition Component.hh:16
Abstract base class for components.
Definition Medium.hh:17
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).
double TransferFunctionSq()
void EnableMagneticField(const unsigned int i, const bool on)
Activate/deactivate use of the magnetic field of a given component.
double GetElectronSignal(const std::string &label, const unsigned int bin)
Retrieve the electron signal for a given electrode and time bin.
void PlotSignal(const std::string &label, TPad *pad, const std::string optTotal="t", const std::string optPrompt="", const std::string optDelayed="")
Plot the induced signal.
bool ConvoluteSignal(const std::string &label, const bool fft=false)
Convolute the induced current on a given electrode with the transfer function.
std::vector< std::tuple< Component *, bool, bool > > m_components
Components.
Definition Sensor.hh:314
void FillBin(Electrode &electrode, const unsigned int bin, const double signal, const bool electron, const bool delayed)
Definition Sensor.hh:372
Medium * GetMedium(const double x, const double y, const double z)
Get the medium at (x, y, z).
bool GetVoltageRange(double &vmin, double &vmax)
Return the voltage range.
bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
double GetTotalInducedCharge(const std::string &label)
double StepSizeHint()
void AddWhiteNoise(const double enc, const bool poisson=true, const double q0=1.)
Add white noise to the induced signals on all electrodes.
void SetDelayedSignalTimes(const std::vector< double > &ts)
Set the points in time at which to evaluate the delayed weighting field.
void SetTimeWindow(const double tstart, const double tstep, const unsigned int nsteps)
Set the time window and binning for the signal calculation.
void SetSignal(const std::string &label, const unsigned int bin, const double signal)
Set/override the signal in a given time bin explicitly.
unsigned int m_nEvents
Definition Sensor.hh:335
void EnableDelayedSignal(const bool on=true)
Compute the component of the signal due to the delayed weighting field.
Definition Sensor.hh:122
double GetDelayedIonSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed ion/hole signal for a given electrode and time bin.
std::vector< double > m_delayedSignalTimes
Definition Sensor.hh:338
bool DelayAndSubtractFraction(const double td, const double f)
Delay the signal and subtract an attenuated copy (modelling a constant fraction discriminator).
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).
Sensor(Component *comp)
Constructor from a single component.
void AddSignalWeightingField(const double q, const std::vector< double > &ts, const std::vector< std::array< double, 3 > > &xs, const std::vector< std::array< double, 3 > > &vs, const std::vector< double > &ns, const int navg)
Calculate the signal from a drift line using the weighting field method, given the drift velocities a...
double GetTransferFunction(const double t)
Evaluate the transfer function at a given time.
double m_xMaxUser
Definition Sensor.hh:359
bool ConvoluteSignalFFT()
void EnableComponent(const unsigned int i, const bool on)
Activate/deactivate a given component.
double InterpolateTransferFunctionTable(const double t) const
bool CrossedWire(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)
Determine whether a line between two points crosses a wire, calls Component::CrossedWire.
bool SetArea(const double xmin, const double ymin, const double zmin, const double xmax, const double ymax, const double zmax)
Set the user area explicitly.
void SetTransferFunction(const std::vector< double > &times, const std::vector< double > &values)
Set the points to be used for interpolating the transfer function.
void FillSignal(Electrode &electrode, const double q, const std::vector< double > &ts, const std::vector< double > &is, const int navg, const bool delayed=false)
double(* m_fNoise)(double t)
Definition Sensor.hh:352
void GetTimeWindow(double &tstart, double &tstep, unsigned int &nsteps) const
Retrieve the time window and binning.
Definition Sensor.hh:114
double m_tStart
Definition Sensor.hh:332
void AddElectrode(Component *comp, const std::string &label)
Add an electrode.
bool IntegrateSignal(const std::string &label)
Replace the signal on a given electrode by its integral.
void AddComponent(Component *comp)
Add a component.
void SetTransferFunction(std::function< double(double)>)
Set the function to be used for evaluating the transfer function.
void PrintTransferFunction()
Print some information about the presently set transfer function.
double DelayedWeightingPotential(const double x, const double y, const double z, const double t, const std::string &label)
Get the delayed weighting potential at (x, y, z).
void ConvoluteSignalFFT(Electrode &electrode, const std::vector< double > &tab, const unsigned int nn)
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).
void AddSignalWeightingPotential(const double q, const std::vector< double > &ts, const std::vector< std::array< double, 3 > > &xs, const std::vector< double > &qs)
Calculate the signal from an avalanche using the weighting potential method.
Shaper * m_shaper
Definition Sensor.hh:343
double GetSignal(const std::string &label, const unsigned int bin)
Retrieve the total signal for a given electrode and time bin.
size_t GetNumberOfComponents() const
Get the number of components attached to the sensor.
Definition Sensor.hh:34
bool m_cacheTransferFunction
Definition Sensor.hh:345
std::function< double(double)> m_fTransfer
Definition Sensor.hh:342
void ConvoluteSignal(Electrode &electrode, const std::vector< double > &tab)
void SetNoiseFunction(double(*f)(double t))
Set the function to be used for evaluating the noise component.
void AddSignalWeightingField(const double q, const std::vector< double > &ts, const std::vector< std::array< double, 3 > > &xs, const bool integrateWeightingField)
Calculate the signal from a drift line using the weighting field method.
std::mutex m_mutex
Mutex.
Definition Sensor.hh:311
void IntegrateSignal(Electrode &electrode)
bool m_delayedSignal
Definition Sensor.hh:337
std::string m_className
Definition Sensor.hh:309
double GetPromptSignal(const std::string &label, const unsigned int bin)
Retrieve the prompt signal for a given electrode and time bin.
bool IntegrateSignals()
Replace the signals on all electrodes curve by their integrals.
bool HasMagneticField() const
Does the sensor have a non-zero magnetic field?
bool ConvoluteSignals(const bool fft=false)
Convolute all induced currents with the transfer function.
std::vector< Electrode > m_electrodes
Electrodes.
Definition Sensor.hh:329
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
bool CrossedPlane(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc)
Determine whether a line between two points crosses a plane, calls Component::CrossedPlane.
std::vector< std::pair< double, double > > m_fTransferTab
Definition Sensor.hh:344
double GetInducedCharge(const std::string &label)
Calculated using the weighting potentials at the start and end points.
void AddNoise(const bool total=true, const bool electron=false, const bool ion=false)
Add noise to the induced signal.
double GetSignal(const std::string &label, const unsigned int bin, const int comp)
Retrieve the electron signal for a given electrode and time bin.
void AddSignalWeightingPotential(const double q, const std::vector< double > &ts, const std::vector< std::array< double, 3 > > &xs)
Calculate the signal from a drift line using the weighting potential method.
void MakeTransferFunctionTable(std::vector< double > &tab)
size_t GetNumberOfElectrodes() const
Get the number of electrodes attached to the sensor.
Definition Sensor.hh:47
unsigned int m_nAvgDelayedSignal
Definition Sensor.hh:339
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)
Integrate the electric field flux through a line from (x0,y0,z0) to (x1,y1,z1) along a direction (xp,...
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&medium, int &status)
Get the drift field at (x, y, z).
double WeightingPotential(const double x, const double y, const double z, const std::string &label)
Get the weighting potential at (x, y, z).
Component * GetComponent(const unsigned int i)
Retrieve the pointer to a given component.
void AddWhiteNoise(const std::string &label, const double enc, const bool poisson=true, const double q0=1.)
Add white noise to the induced signal, given a desired output ENC.
bool ComputeThresholdCrossings(const double thr, const std::string &label, int &n)
Determine the threshold crossings of the current signal curve.
void Clear()
Remove all components, electrodes and reset the sensor.
double m_yMinUser
Definition Sensor.hh:358
bool SetArea(const bool verbose=false)
Set the user area to the default.
void SetSignal(const std::string &label, const std::vector< double > &ts, const std::vector< double > &is)
Set/override the signal.
~Sensor()=default
Destructor.
bool GetThresholdCrossing(const unsigned int i, double &time, double &level, bool &rise) const
Retrieve the time and type of a given threshold crossing (after having called ComputeThresholdCrossin...
size_t GetNumberOfThresholdCrossings() const
Get the number of threshold crossings (after having called ComputeThresholdCrossings).
Definition Sensor.hh:224
void ClearElectrodes()
Remove all electrodes.
bool IsInside(const double x, const double y, const double z)
Check if a point is inside an active medium and inside the user area.
double m_fTransferSq
Definition Sensor.hh:347
Sensor()=default
Default constructor.
std::vector< std::pair< double, bool > > m_thresholdCrossings
Definition Sensor.hh:354
unsigned int m_nTimeBins
Definition Sensor.hh:334
double GetDelayedElectronSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed electron signal for a given electrode and time bin.
void EnableTransferFunctionCache(const bool on=true)
Cache integral and FFT of the transfer function instead of recomputing it at every call (default: on)...
Definition Sensor.hh:175
double CreateGPUTransferObject(SensorGPU *&sensor_gpu)
Create and initialise GPU Transfer class.
double m_xMinUser
Definition Sensor.hh:358
void SetTransferFunction(Shaper &shaper)
Set the transfer function using a Shaper object.
double m_yMaxUser
Definition Sensor.hh:359
double m_zMaxUser
Definition Sensor.hh:359
void ClearSignal()
Reset signals and induced charges of all electrodes.
double GetDelayedSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed signal for a given electrode and time bin.
void ExportSignal(const std::string &label, const std::string &filename, const bool chargeCariers=false) const
Exporting induced signal to a csv file.
std::vector< double > m_fTransferFFT
Definition Sensor.hh:349
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.
void PlotTransferFunction()
Plot the presently set transfer function.
double m_thresholdLevel
Definition Sensor.hh:355
void FFT(std::vector< double > &data, const bool inverse, const int nn)
bool InTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
Determine whether a point is in the trap radius of a wire.
double GetIonSignal(const std::string &label, const unsigned int bin)
Retrieve the ion or hole signal for a given electrode and time bin.
void EnableDebugging(const bool on=true)
Switch debugging messages on/off.
Definition Sensor.hh:82
void SetDelayedSignalAveragingOrder(const unsigned int navg)
Set the number of points to be used when averaging the delayed signal vector over a time bin (default...
Definition Sensor.hh:129
bool ConvoluteSignalFFT(const std::string &label)
double m_zMinUser
Definition Sensor.hh:358
void NewSignal()
Start a new event, when computing the average signal over multiple events.
Definition Sensor.hh:102
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
bool IsIntegrated(const std::string &label) const
Return whether the signal has been integrated/convoluted.
Class for signal processing.
Definition Shaper.hh:10
std::vector< double > delayedSignal
Definition Sensor.hh:320
std::vector< double > delayedElectronSignal
Definition Sensor.hh:323
std::vector< double > electronSignal
Definition Sensor.hh:321
std::vector< double > signal
Definition Sensor.hh:319
std::vector< double > delayedIonSignal
Definition Sensor.hh:324
std::vector< double > ionSignal
Definition Sensor.hh:322