43 void ElectricField(
const double x,
const double y,
const double z,
double& ex,
44 double& ey,
double& ez,
double& v,
Medium*& medium,
47 void ElectricField(
const double x,
const double y,
const double z,
double& ex,
48 double& ey,
double& ez,
Medium*& medium,
int& status);
51 void MagneticField(
const double x,
const double y,
const double z,
double& bx,
52 double& by,
double& bz,
int& status);
55 void WeightingField(
const double x,
const double y,
const double z,
56 double& wx,
double& wy,
double& wz,
57 const std::string& label);
60 const std::string& label);
63 bool GetMedium(
const double x,
const double y,
const double z,
72 bool SetArea(
const double xmin,
const double ymin,
const double zmin,
73 const double xmax,
const double ymax,
const double zmax);
75 bool GetArea(
double& xmin,
double& ymin,
double& zmin,
double& xmax,
76 double& ymax,
double& zmax);
78 bool IsInArea(
const double x,
const double y,
const double z);
94 const unsigned int nsteps);
97 unsigned int& nsteps)
const {
100 nsteps = m_nTimeBins;
112 m_nAvgDelayedSignal = navg;
116 void SetSignal(
const std::string& label,
const unsigned int bin,
117 const double signal);
119 double GetSignal(
const std::string& label,
const unsigned int bin);
121 double GetSignal(
const std::string& label,
const unsigned int bin,
124 double GetPromptSignal(
const std::string& label,
const unsigned int bin);
130 double GetIonSignal(
const std::string& label,
const unsigned int bin);
133 const unsigned int bin);
143 const std::vector<double>& values);
153 m_cacheTransferFunction = on;
157 bool ConvoluteSignal(
const std::string& label,
const bool fft =
false);
176 void AddNoise(
const bool total =
true,
const bool electron =
false,
177 const bool ion =
false);
186 void AddWhiteNoise(
const std::string& label,
const double enc,
187 const bool poisson =
true,
const double q0 = 1.);
189 void AddWhiteNoise(
const double enc,
const bool poisson =
true,
190 const double q0 = 1.);
202 return m_thresholdCrossings.size();
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);
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);
229 void ExportSignal(
const std::string& label,
const std::string& filename);
233 const double z0,
const double x1,
const double y1,
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,
242 bool IsInTrapRadius(
const double q0,
const double x0,
const double y0,
243 const double z0,
double& xw,
double& yw,
double& rw);
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);
254 std::string m_className =
"Sensor";
259 std::vector<std::tuple<Component*, bool, bool> > m_components;
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;
275 std::vector<Electrode> m_electrodes;
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;
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;
292 double m_fTransferSq = -1.;
294 std::vector<double> m_fTransferFFT;
297 double (*m_fNoise)(
double t) =
nullptr;
299 std::vector<std::pair<double, bool> > m_thresholdCrossings;
300 double m_thresholdLevel = 0.;
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;
308 bool m_debug =
false;
311 bool GetBoundingBox(
double& xmin,
double& ymin,
double& zmin,
double& xmax,
312 double& ymax,
double& zmax);
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;
323 electrode.electronSignal[bin] += signal;
324 if (delayed) electrode.delayedElectronSignal[bin] += signal;
326 electrode.ionSignal[bin] += signal;
327 if (delayed) electrode.delayedIonSignal[bin] += signal;
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);
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);
Abstract base class for components.
Abstract base class for media.
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).
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.
bool ConvoluteSignal(const std::string &label, const bool fft=false)
double GetTotalInducedCharge(const std::string label)
bool GetVoltageRange(double &vmin, double &vmax)
Return the voltage range.
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)
void ExportSignal(const std::string &label, const std::string &filename)
Exporting induced signal to a csv file.
void SetSignal(const std::string &label, const unsigned int bin, const double signal)
Set/override the signal in a given time bin explicitly.
void EnableDelayedSignal(const bool on=true)
Compute the component of the signal due to the delayed weighting field.
double GetDelayedIonSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed ion/hole signal for a given electrode and time bin.
bool IsInTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
bool DelayAndSubtractFraction(const double td, const double f)
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).
double GetTransferFunction(const double t)
Evaluate the transfer function at a given time.
void EnableComponent(const unsigned int i, const bool on)
Activate/deactivate a given component.
void GetTimeWindow(double &tstart, double &tstep, unsigned int &nsteps) const
Retrieve the time window and binning.
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 PrintTransferFunction()
Print some information about the presently set transfer function.
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).
bool SetArea()
Set the user area to the default.
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.
void SetNoiseFunction(double(*f)(double t))
Set the function to be used for evaluating the noise component.
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 ConvoluteSignals(const bool fft=false)
Convolute all induced currents with the transfer function.
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
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.
size_t GetNumberOfElectrodes() const
Get the number of electrodes attached to the sensor.
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)
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.)
bool ComputeThresholdCrossings(const double thr, const std::string &label, int &n)
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
void Clear()
Remove all components, electrodes and reset the sensor.
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.
void SetTransferFunction(double(*f)(double t))
Set the function to be used for evaluating the transfer function.
bool GetThresholdCrossing(const unsigned int i, double &time, double &level, bool &rise) const
size_t GetNumberOfThresholdCrossings() const
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.
void EnableTransferFunctionCache(const bool on=true)
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)
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 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.
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.
void SetDelayedSignalAveragingOrder(const unsigned int navg)
void NewSignal()
Start a new event, when computing the average signal over multiple events.
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.