Garfield++ v2r0
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 <vector>
5
6#include "ComponentBase.hh"
7
8namespace Garfield {
9
10/// %Sensor
11
12class Sensor {
13
14 public:
15 /// Constructor
16 Sensor();
17 /// Destructor
19
20 /// Add a component.
21 void AddComponent(ComponentBase* comp);
22 unsigned int GetNumberOfComponents() const { return m_components.size(); }
23 ComponentBase* GetComponent(const unsigned int componentNumber);
24
25 /// Add an electrode.
26 void AddElectrode(ComponentBase* comp, const std::string& label);
27 unsigned int GetNumberOfElectrodes() const { return m_electrodes.size(); }
28 /// Remove all components, electrodes and reset the sensor.
29 void Clear();
30
31 /// Get the drift field and potential at (x, y, z).
32 void ElectricField(const double x, const double y, const double z, double& ex,
33 double& ey, double& ez, double& v, Medium*& medium,
34 int& status);
35 /// Get the drift field at (x, y, z).
36 void ElectricField(const double x, const double y, const double z, double& ex,
37 double& ey, double& ez, Medium*& medium, int& status);
38
39 /// Get the magnetic field at (x, y, z).
40 void MagneticField(const double x, const double y, const double z, double& bx,
41 double& by, double& bz, int& status);
42
43 /// Get the weighting field at (x, y, z).
44 void WeightingField(const double x, const double y, const double z,
45 double& wx, double& wy, double& wz,
46 const std::string& label);
47 /// Get the weighting potential at (x, y, z).
48 double WeightingPotential(const double x, const double y, const double z,
49 const std::string& label);
50
51 /// Get the medium at (x, y, z).
52 bool GetMedium(const double x, const double y, const double z,
53 Medium*& medium);
54
55 /// Set the user area to the default.
56 bool SetArea();
57 /// Set the user area explicitly.
58 bool SetArea(const double xmin, const double ymin, const double zmin,
59 const double xmax, const double ymax, const double zmax);
60 /// Return the current user area.
61 bool GetArea(double& xmin, double& ymin, double& zmin, double& xmax,
62 double& ymax, double& zmax);
63 /// Check if a point is inside the user area.
64 bool IsInArea(const double x, const double y, const double z);
65
66 bool IsWireCrossed(const double x0, const double y0, const double z0,
67 const double x1, const double y1, const double z1,
68 double& xc, double& yc, double& zc);
69
70 bool IsInTrapRadius(const double q0, const double x0, const double y0,
71 const double z0, double& xw, double& yw, double& rw);
72
73 /// Return the voltage range.
74 bool GetVoltageRange(double& vmin, double& vmax);
75
76 /// Start a new event, when computing the average signal over multiple events.
77 void NewSignal() { ++m_nEvents; }
78 /// Reset signals and induced charges of all electrodes.
79 void ClearSignal();
80 void AddSignal(const double q, const double t, const double dt,
81 const double x, const double y, const double z,
82 const double vx, const double vy, const double vz);
83 void AddInducedCharge(const double q, const double x0, const double y0,
84 const double z0, const double x1, const double y1,
85 const double z1);
86 // Set/get the time window and binning for the signal calculation
87 void SetTimeWindow(const double tstart, const double tstep,
88 const unsigned int nsteps);
89 void GetTimeWindow(double& tstart, double& tstep, unsigned int& nsteps) {
90 tstart = m_tStart;
91 tstep = m_tStep;
92 nsteps = m_nTimeBins;
93 }
94 double GetSignal(const std::string& label, const unsigned int bin);
95 double GetElectronSignal(const std::string& label, const unsigned int bin);
96 double GetIonSignal(const std::string& label, const unsigned int bin);
97 double GetInducedCharge(const std::string& label);
98 void SetTransferFunction(double (*f)(double t));
99 void SetTransferFunction(const std::vector<double>& times,
100 const std::vector<double>& values);
101 double GetTransferFunction(const double t);
102 bool ConvoluteSignal();
103 bool IntegrateSignal();
104 void SetNoiseFunction(double (*f)(double t));
105 void AddNoise(const bool total = true, const bool electron = false,
106 const bool ion = false);
107 bool ComputeThresholdCrossings(const double thr, const std::string& label,
108 int& n);
109 unsigned int GetNumberOfThresholdCrossings() const {
110 return m_thresholdCrossings.size();
111 }
112 bool GetThresholdCrossing(const unsigned int i, double& time, double& level,
113 bool& rise) const;
114
115 /// Switch on debugging messages
116 void EnableDebugging() { m_debug = true; }
117 void DisableDebugging() { m_debug = false; }
118
119 private:
120 std::string m_className;
121
122 // Components
123 struct component {
124 ComponentBase* comp;
125 };
126 std::vector<component> m_components;
127 int m_lastComponent;
128
129 // Electrodes
130 struct electrode {
131 ComponentBase* comp;
132 std::string label;
133 std::vector<double> signal;
134 std::vector<double> electronsignal;
135 std::vector<double> ionsignal;
136 double charge;
137 };
138 std::vector<electrode> m_electrodes;
139
140 // Time window for signals
141 unsigned int m_nTimeBins;
142 double m_tStart, m_tStep;
143 unsigned int m_nEvents;
144 static double m_signalConversion;
145
146 // Transfer function
147 bool m_hasTransferFunction;
148 double (*m_fTransfer)(double t);
149 std::vector<double> m_transferFunctionTimes;
150 std::vector<double> m_transferFunctionValues;
151
152 // Noise
153 bool m_hasNoiseFunction;
154 double (*m_fNoise)(double t);
155
156 struct thresholdCrossing {
157 double time;
158 bool rise;
159 };
160 std::vector<thresholdCrossing> m_thresholdCrossings;
161 double m_thresholdLevel;
162
163 // User bounding box
164 bool m_hasUserArea;
165 double m_xMinUser, m_yMinUser, m_zMinUser;
166 double m_xMaxUser, m_yMaxUser, m_zMaxUser;
167
168 // Switch on/off debugging messages
169 bool m_debug;
170
171 // Return the current sensor size
172 bool GetBoundingBox(double& xmin, double& ymin, double& zmin, double& xmax,
173 double& ymax, double& zmax);
174
175 double InterpolateTransferFunctionTable(double t);
176};
177}
178
179#endif
Abstract base class for components.
Abstract base class for media.
Definition: Medium.hh:11
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:101
double GetElectronSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:563
bool GetVoltageRange(double &vmin, double &vmax)
Return the voltage range.
Definition: Sensor.cc:369
void SetTimeWindow(const double tstart, const double tstep, const unsigned int nsteps)
Definition: Sensor.cc:527
bool IntegrateSignal()
Definition: Sensor.cc:760
bool IsInTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
Definition: Sensor.cc:301
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:117
double GetTransferFunction(const double t)
Definition: Sensor.cc:695
~Sensor()
Destructor.
Definition: Sensor.hh:18
Sensor()
Constructor.
Definition: Sensor.cc:15
unsigned int GetNumberOfElectrodes() const
Definition: Sensor.hh:27
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:48
bool SetArea()
Set the user area to the default.
Definition: Sensor.cc:175
double GetSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:601
void SetNoiseFunction(double(*f)(double t))
Definition: Sensor.cc:785
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
Definition: Sensor.cc:264
void AddElectrode(ComponentBase *comp, const std::string &label)
Add an electrode.
Definition: Sensor.cc:327
double GetInducedCharge(const std::string &label)
Definition: Sensor.cc:619
void AddNoise(const bool total=true, const bool electron=false, const bool ion=false)
Definition: Sensor.cc:796
void AddSignal(const double q, const double t, const double dt, const double x, const double y, const double z, const double vx, const double vy, const double vz)
Definition: Sensor.cc:417
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:136
bool ComputeThresholdCrossings(const double thr, const std::string &label, int &n)
Definition: Sensor.cc:817
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:150
ComponentBase * GetComponent(const unsigned int componentNumber)
Definition: Sensor.cc:38
void Clear()
Remove all components, electrodes and reset the sensor.
Definition: Sensor.cc:357
void SetTransferFunction(double(*f)(double t))
Definition: Sensor.cc:636
bool GetThresholdCrossing(const unsigned int i, double &time, double &level, bool &rise) const
Definition: Sensor.cc:960
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)
Definition: Sensor.cc:288
bool ConvoluteSignal()
Definition: Sensor.cc:702
void AddComponent(ComponentBase *comp)
Add a component.
Definition: Sensor.cc:314
unsigned int GetNumberOfThresholdCrossings() const
Definition: Sensor.hh:109
void DisableDebugging()
Definition: Sensor.hh:117
void EnableDebugging()
Switch on debugging messages.
Definition: Sensor.hh:116
void ClearSignal()
Reset signals and induced charges of all electrodes.
Definition: Sensor.cc:403
void GetTimeWindow(double &tstart, double &tstep, unsigned int &nsteps)
Definition: Sensor.hh:89
void AddInducedCharge(const double q, const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Definition: Sensor.cc:501
unsigned int GetNumberOfComponents() const
Definition: Sensor.hh:22
double GetIonSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:583
void NewSignal()
Start a new event, when computing the average signal over multiple events.
Definition: Sensor.hh:77
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
Definition: Sensor.cc:237