Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::Sensor Class Reference

Sensor. More...

#include <Sensor.hh>

Classes

struct  Electrode

Public Member Functions

 Sensor ()=default
 Default constructor.
 ~Sensor ()=default
 Destructor.
 Sensor (Component *comp)
 Constructor from a single component.
void AddComponent (Component *comp)
 Add a component.
size_t GetNumberOfComponents () const
 Get the number of components attached to the sensor.
ComponentGetComponent (const unsigned int i)
 Retrieve the pointer to a given component.
void EnableComponent (const unsigned int i, const bool on)
 Activate/deactivate a given component.
void EnableMagneticField (const unsigned int i, const bool on)
 Activate/deactivate use of the magnetic field of a given component.
bool HasMagneticField () const
 Does the sensor have a non-zero magnetic field?
void AddElectrode (Component *comp, const std::string &label)
 Add an electrode.
size_t GetNumberOfElectrodes () const
 Get the number of electrodes attached to the sensor.
void ClearElectrodes ()
 Remove all electrodes.
void Clear ()
 Remove all components, electrodes and reset the sensor.
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 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).
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 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 WeightingPotential (const double x, const double y, const double z, const std::string &label)
 Get the weighting potential at (x, y, z).
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).
MediumGetMedium (const double x, const double y, const double z)
 Get the medium at (x, y, z).
void EnableDebugging (const bool on=true)
 Switch debugging messages on/off.
bool SetArea (const bool verbose=false)
 Set the user area to the default.
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.
bool GetArea (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
 Return the current user area.
bool IsInArea (const double x, const double y, const double z)
 Check if a point is inside the user area.
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.
bool GetVoltageRange (double &vmin, double &vmax)
 Return the voltage range.
void NewSignal ()
 Start a new event, when computing the average signal over multiple events.
void ClearSignal ()
 Reset signals and induced charges of all electrodes.
void SetTimeWindow (const double tstart, const double tstep, const unsigned int nsteps)
 Set the time window and binning for the signal calculation.
void GetTimeWindow (double &tstart, double &tstep, unsigned int &nsteps) const
 Retrieve the time window and binning.
void EnableDelayedSignal (const bool on=true)
 Compute the component of the signal due to the delayed weighting field.
void SetDelayedSignalTimes (const std::vector< double > &ts)
 Set the points in time at which to evaluate the delayed weighting field.
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: 0).
void SetSignal (const std::string &label, const unsigned int bin, const double signal)
 Set/override the signal in a given time bin explicitly.
void SetSignal (const std::string &label, const std::vector< double > &ts, const std::vector< double > &is)
 Set/override the signal.
double GetSignal (const std::string &label, const unsigned int bin)
 Retrieve the total signal for a given electrode and time bin.
double GetSignal (const std::string &label, const unsigned int bin, const int comp)
 Retrieve the electron signal for a given electrode and time bin.
double GetPromptSignal (const std::string &label, const unsigned int bin)
 Retrieve the prompt signal for a given electrode and time bin.
double GetDelayedSignal (const std::string &label, const unsigned int bin)
 Retrieve the delayed signal for a given electrode and time bin.
double GetElectronSignal (const std::string &label, const unsigned int bin)
 Retrieve the electron signal for a given electrode and time bin.
double GetIonSignal (const std::string &label, const unsigned int bin)
 Retrieve the ion or hole signal for a given electrode and time bin.
double GetDelayedElectronSignal (const std::string &label, const unsigned int bin)
 Retrieve the delayed electron signal for a given electrode and time bin.
double GetDelayedIonSignal (const std::string &label, const unsigned int bin)
 Retrieve the delayed ion/hole signal for a given electrode and time bin.
double GetInducedCharge (const std::string &label)
 Calculated using the weighting potentials at the start and end points.
void SetTransferFunction (std::function< double(double)>)
 Set the function to be used for evaluating the transfer function.
void SetTransferFunction (const std::vector< double > &times, const std::vector< double > &values)
 Set the points to be used for interpolating the transfer function.
void SetTransferFunction (Shaper &shaper)
 Set the transfer function using a Shaper object.
double GetTransferFunction (const double t)
 Evaluate the transfer function at a given time.
void PrintTransferFunction ()
 Print some information about the presently set transfer function.
void PlotTransferFunction ()
 Plot the presently set transfer function.
void EnableTransferFunctionCache (const bool on=true)
 Cache integral and FFT of the transfer function instead of recomputing it at every call (default: on).
bool ConvoluteSignal (const std::string &label, const bool fft=false)
 Convolute the induced current on a given electrode with the transfer function.
bool ConvoluteSignals (const bool fft=false)
 Convolute all induced currents with the transfer function.
bool IntegrateSignal (const std::string &label)
 Replace the signal on a given electrode by its integral.
bool IntegrateSignals ()
 Replace the signals on all electrodes curve by their integrals.
bool IsIntegrated (const std::string &label) const
 Return whether the signal has been integrated/convoluted.
bool DelayAndSubtractFraction (const double td, const double f)
 Delay the signal and subtract an attenuated copy (modelling a constant fraction discriminator).
void SetNoiseFunction (double(*f)(double t))
 Set the function to be used for evaluating the noise component.
void AddNoise (const bool total=true, const bool electron=false, const bool ion=false)
 Add noise to the induced signal.
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.
void AddWhiteNoise (const double enc, const bool poisson=true, const double q0=1.)
 Add white noise to the induced signals on all electrodes.
bool ComputeThresholdCrossings (const double thr, const std::string &label, int &n)
 Determine the threshold crossings of the current signal curve.
size_t GetNumberOfThresholdCrossings () const
 Get the number of threshold crossings (after having called ComputeThresholdCrossings).
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 ComputeThresholdCrossings.
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 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.
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.
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 at each drift line point.
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.
void ExportSignal (const std::string &label, const std::string &filename, const bool chargeCariers=false) const
 Exporting induced signal to a csv file.
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.
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 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.
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.
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,yp,zp).
double GetTotalInducedCharge (const std::string &label)
double StepSizeHint ()
double CreateGPUTransferObject (SensorGPU *&sensor_gpu)
 Create and initialise GPU Transfer class.

Private Member Functions

bool GetBoundingBox (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
void FillSignal (Electrode &electrode, const double q, const std::vector< double > &ts, const std::vector< double > &is, const int navg, const bool delayed=false)
void FillBin (Electrode &electrode, const unsigned int bin, const double signal, const bool electron, const bool delayed)
void IntegrateSignal (Electrode &electrode)
void ConvoluteSignal (Electrode &electrode, const std::vector< double > &tab)
bool ConvoluteSignalFFT ()
bool ConvoluteSignalFFT (const std::string &label)
void ConvoluteSignalFFT (Electrode &electrode, const std::vector< double > &tab, const unsigned int nn)
double TransferFunctionSq ()
double InterpolateTransferFunctionTable (const double t) const
void MakeTransferFunctionTable (std::vector< double > &tab)
void FFT (std::vector< double > &data, const bool inverse, const int nn)

Private Attributes

std::string m_className = "Sensor"
std::mutex m_mutex
 Mutex.
std::vector< std::tuple< Component *, bool, bool > > m_components
 Components.
std::vector< Electrodem_electrodes
 Electrodes.
double m_tStart = 0.
double m_tStep = 10.
unsigned int m_nTimeBins = 200
unsigned int m_nEvents = 0
bool m_delayedSignal = false
std::vector< double > m_delayedSignalTimes
unsigned int m_nAvgDelayedSignal = 0
std::function< double(double)> m_fTransfer
Shaperm_shaper = nullptr
std::vector< std::pair< double, double > > m_fTransferTab
bool m_cacheTransferFunction = true
double m_fTransferSq = -1.
std::vector< double > m_fTransferFFT
double(* m_fNoise )(double t) = nullptr
std::vector< std::pair< double, bool > > m_thresholdCrossings
double m_thresholdLevel = 0.
double m_xMinUser = 0.
double m_yMinUser = 0.
double m_zMinUser = 0.
double m_xMaxUser = 0.
double m_yMaxUser = 0.
double m_zMaxUser = 0.
bool m_hasUserArea = false
bool m_debug = false

Detailed Description

Sensor.

Definition at line 22 of file Sensor.hh.

Constructor & Destructor Documentation

◆ Sensor() [1/2]

Garfield::Sensor::Sensor ( )
default

Default constructor.

◆ ~Sensor()

Garfield::Sensor::~Sensor ( )
default

Destructor.

◆ Sensor() [2/2]

Garfield::Sensor::Sensor ( Component * comp)

Constructor from a single component.

Member Function Documentation

◆ AddComponent()

void Garfield::Sensor::AddComponent ( Component * comp)

Add a component.

◆ AddElectrode()

void Garfield::Sensor::AddElectrode ( Component * comp,
const std::string & label )

Add an electrode.

◆ AddInducedCharge()

void Garfield::Sensor::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.

◆ AddNoise()

void Garfield::Sensor::AddNoise ( const bool total = true,
const bool electron = false,
const bool ion = false )

Add noise to the induced signal.

◆ AddSignalWeightingField() [1/2]

void Garfield::Sensor::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.

◆ AddSignalWeightingField() [2/2]

void Garfield::Sensor::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 at each drift line point.

◆ AddSignalWeightingPotential() [1/2]

void Garfield::Sensor::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.

◆ AddSignalWeightingPotential() [2/2]

void Garfield::Sensor::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.

◆ AddWhiteNoise() [1/2]

void Garfield::Sensor::AddWhiteNoise ( const double enc,
const bool poisson = true,
const double q0 = 1. )

Add white noise to the induced signals on all electrodes.

◆ AddWhiteNoise() [2/2]

void Garfield::Sensor::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.

Parameters
labelname of the electrode
encEquivalent Noise Charge, in electrons.
poissonflag to sample the number of noise pulses from a Poisson distribution. Otherwise the noise charge in each bin is sampled from a Gaussian distribution.
q0amplitude of the noise delta pulses (in electrons).

◆ Clear()

void Garfield::Sensor::Clear ( )

Remove all components, electrodes and reset the sensor.

◆ ClearElectrodes()

void Garfield::Sensor::ClearElectrodes ( )

Remove all electrodes.

◆ ClearSignal()

void Garfield::Sensor::ClearSignal ( )

Reset signals and induced charges of all electrodes.

◆ ComputeThresholdCrossings()

bool Garfield::Sensor::ComputeThresholdCrossings ( const double thr,
const std::string & label,
int & n )

Determine the threshold crossings of the current signal curve.

Parameters
thrthreshold value
labelelectrode for which to compute the threshold crossings
nnumber of threshold crossings

◆ ConvoluteSignal() [1/2]

bool Garfield::Sensor::ConvoluteSignal ( const std::string & label,
const bool fft = false )

Convolute the induced current on a given electrode with the transfer function.

◆ ConvoluteSignal() [2/2]

void Garfield::Sensor::ConvoluteSignal ( Electrode & electrode,
const std::vector< double > & tab )
private

◆ ConvoluteSignalFFT() [1/3]

bool Garfield::Sensor::ConvoluteSignalFFT ( )
private

◆ ConvoluteSignalFFT() [2/3]

bool Garfield::Sensor::ConvoluteSignalFFT ( const std::string & label)
private

◆ ConvoluteSignalFFT() [3/3]

void Garfield::Sensor::ConvoluteSignalFFT ( Electrode & electrode,
const std::vector< double > & tab,
const unsigned int nn )
private

◆ ConvoluteSignals()

bool Garfield::Sensor::ConvoluteSignals ( const bool fft = false)

Convolute all induced currents with the transfer function.

◆ CreateGPUTransferObject()

double Garfield::Sensor::CreateGPUTransferObject ( SensorGPU *& sensor_gpu)

Create and initialise GPU Transfer class.

◆ CrossedPlane()

bool Garfield::Sensor::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.

◆ CrossedWire()

bool Garfield::Sensor::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.

◆ DelayAndSubtractFraction()

bool Garfield::Sensor::DelayAndSubtractFraction ( const double td,
const double f )

Delay the signal and subtract an attenuated copy (modelling a constant fraction discriminator).

\[  v_{out} = v_{in}\left(t - t_d\right) - f v_{in}.
\]

◆ DelayedWeightingPotential()

double Garfield::Sensor::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).

◆ ElectricField() [1/2]

void Garfield::Sensor::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).

◆ ElectricField() [2/2]

void Garfield::Sensor::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).

◆ EnableComponent()

void Garfield::Sensor::EnableComponent ( const unsigned int i,
const bool on )

Activate/deactivate a given component.

◆ EnableDebugging()

void Garfield::Sensor::EnableDebugging ( const bool on = true)
inline

Switch debugging messages on/off.

Definition at line 82 of file Sensor.hh.

82{ m_debug = on; }

◆ EnableDelayedSignal()

void Garfield::Sensor::EnableDelayedSignal ( const bool on = true)
inline

Compute the component of the signal due to the delayed weighting field.

Definition at line 122 of file Sensor.hh.

122{ m_delayedSignal = on; }
bool m_delayedSignal
Definition Sensor.hh:337

◆ EnableMagneticField()

void Garfield::Sensor::EnableMagneticField ( const unsigned int i,
const bool on )

Activate/deactivate use of the magnetic field of a given component.

◆ EnableTransferFunctionCache()

void Garfield::Sensor::EnableTransferFunctionCache ( const bool on = true)
inline

Cache integral and FFT of the transfer function instead of recomputing it at every call (default: on).

Definition at line 175 of file Sensor.hh.

175 {
177 }
bool m_cacheTransferFunction
Definition Sensor.hh:345

◆ ExportSignal()

void Garfield::Sensor::ExportSignal ( const std::string & label,
const std::string & filename,
const bool chargeCariers = false ) const

Exporting induced signal to a csv file.

◆ FFT()

void Garfield::Sensor::FFT ( std::vector< double > & data,
const bool inverse,
const int nn )
private

◆ FillBin()

void Garfield::Sensor::FillBin ( Electrode & electrode,
const unsigned int bin,
const double signal,
const bool electron,
const bool delayed )
inlineprivate

Definition at line 372 of file Sensor.hh.

373 {
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 }
std::mutex m_mutex
Mutex.
Definition Sensor.hh:311

◆ FillSignal()

void Garfield::Sensor::FillSignal ( Electrode & electrode,
const double q,
const std::vector< double > & ts,
const std::vector< double > & is,
const int navg,
const bool delayed = false )
private

◆ GetArea()

bool Garfield::Sensor::GetArea ( double & xmin,
double & ymin,
double & zmin,
double & xmax,
double & ymax,
double & zmax )

Return the current user area.

◆ GetBoundingBox()

bool Garfield::Sensor::GetBoundingBox ( double & xmin,
double & ymin,
double & zmin,
double & xmax,
double & ymax,
double & zmax )
private

◆ GetComponent()

Component * Garfield::Sensor::GetComponent ( const unsigned int i)

Retrieve the pointer to a given component.

◆ GetDelayedElectronSignal()

double Garfield::Sensor::GetDelayedElectronSignal ( const std::string & label,
const unsigned int bin )

Retrieve the delayed electron signal for a given electrode and time bin.

◆ GetDelayedIonSignal()

double Garfield::Sensor::GetDelayedIonSignal ( const std::string & label,
const unsigned int bin )

Retrieve the delayed ion/hole signal for a given electrode and time bin.

◆ GetDelayedSignal()

double Garfield::Sensor::GetDelayedSignal ( const std::string & label,
const unsigned int bin )

Retrieve the delayed signal for a given electrode and time bin.

◆ GetElectronSignal()

double Garfield::Sensor::GetElectronSignal ( const std::string & label,
const unsigned int bin )

Retrieve the electron signal for a given electrode and time bin.

◆ GetInducedCharge()

double Garfield::Sensor::GetInducedCharge ( const std::string & label)

Calculated using the weighting potentials at the start and end points.

◆ GetIonSignal()

double Garfield::Sensor::GetIonSignal ( const std::string & label,
const unsigned int bin )

Retrieve the ion or hole signal for a given electrode and time bin.

◆ GetMedium()

Medium * Garfield::Sensor::GetMedium ( const double x,
const double y,
const double z )

Get the medium at (x, y, z).

◆ GetNumberOfComponents()

size_t Garfield::Sensor::GetNumberOfComponents ( ) const
inline

Get the number of components attached to the sensor.

Definition at line 34 of file Sensor.hh.

34{ return m_components.size(); }
std::vector< std::tuple< Component *, bool, bool > > m_components
Components.
Definition Sensor.hh:314

◆ GetNumberOfElectrodes()

size_t Garfield::Sensor::GetNumberOfElectrodes ( ) const
inline

Get the number of electrodes attached to the sensor.

Definition at line 47 of file Sensor.hh.

47{ return m_electrodes.size(); }
std::vector< Electrode > m_electrodes
Electrodes.
Definition Sensor.hh:329

◆ GetNumberOfThresholdCrossings()

size_t Garfield::Sensor::GetNumberOfThresholdCrossings ( ) const
inline

Get the number of threshold crossings (after having called ComputeThresholdCrossings).

Definition at line 224 of file Sensor.hh.

224 {
225 return m_thresholdCrossings.size();
226 }
std::vector< std::pair< double, bool > > m_thresholdCrossings
Definition Sensor.hh:354

◆ GetPromptSignal()

double Garfield::Sensor::GetPromptSignal ( const std::string & label,
const unsigned int bin )

Retrieve the prompt signal for a given electrode and time bin.

◆ GetSignal() [1/2]

double Garfield::Sensor::GetSignal ( const std::string & label,
const unsigned int bin )

Retrieve the total signal for a given electrode and time bin.

◆ GetSignal() [2/2]

double Garfield::Sensor::GetSignal ( const std::string & label,
const unsigned int bin,
const int comp )

Retrieve the electron signal for a given electrode and time bin.

◆ GetThresholdCrossing()

bool Garfield::Sensor::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 ComputeThresholdCrossings.

Parameters
iindex
timethreshold crossing time [ns]
levelthreshold (should correspond to the value requested).
riseflag whether the crossing is on a rising or falling slope.

◆ GetTimeWindow()

void Garfield::Sensor::GetTimeWindow ( double & tstart,
double & tstep,
unsigned int & nsteps ) const
inline

Retrieve the time window and binning.

Definition at line 114 of file Sensor.hh.

115 {
116 tstart = m_tStart;
117 tstep = m_tStep;
118 nsteps = m_nTimeBins;
119 }
double m_tStart
Definition Sensor.hh:332
unsigned int m_nTimeBins
Definition Sensor.hh:334

◆ GetTotalInducedCharge()

double Garfield::Sensor::GetTotalInducedCharge ( const std::string & label)

◆ GetTransferFunction()

double Garfield::Sensor::GetTransferFunction ( const double t)

Evaluate the transfer function at a given time.

◆ GetVoltageRange()

bool Garfield::Sensor::GetVoltageRange ( double & vmin,
double & vmax )

Return the voltage range.

◆ HasMagneticField()

bool Garfield::Sensor::HasMagneticField ( ) const

Does the sensor have a non-zero magnetic field?

◆ IntegrateFluxLine()

double Garfield::Sensor::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,yp,zp).

◆ IntegrateSignal() [1/2]

bool Garfield::Sensor::IntegrateSignal ( const std::string & label)

Replace the signal on a given electrode by its integral.

◆ IntegrateSignal() [2/2]

void Garfield::Sensor::IntegrateSignal ( Electrode & electrode)
private

◆ IntegrateSignals()

bool Garfield::Sensor::IntegrateSignals ( )

Replace the signals on all electrodes curve by their integrals.

◆ InterpolateTransferFunctionTable()

double Garfield::Sensor::InterpolateTransferFunctionTable ( const double t) const
private

◆ InTrapRadius()

bool Garfield::Sensor::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.

◆ IsInArea()

bool Garfield::Sensor::IsInArea ( const double x,
const double y,
const double z )

Check if a point is inside the user area.

◆ IsInside()

bool Garfield::Sensor::IsInside ( const double x,
const double y,
const double z )

Check if a point is inside an active medium and inside the user area.

◆ IsIntegrated()

bool Garfield::Sensor::IsIntegrated ( const std::string & label) const

Return whether the signal has been integrated/convoluted.

◆ MagneticField()

void Garfield::Sensor::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).

◆ MakeTransferFunctionTable()

void Garfield::Sensor::MakeTransferFunctionTable ( std::vector< double > & tab)
private

◆ NewSignal()

void Garfield::Sensor::NewSignal ( )
inline

Start a new event, when computing the average signal over multiple events.

Definition at line 102 of file Sensor.hh.

102{ ++m_nEvents; }
unsigned int m_nEvents
Definition Sensor.hh:335

◆ PlotSignal()

void Garfield::Sensor::PlotSignal ( const std::string & label,
TPad * pad,
const std::string optTotal = "t",
const std::string optPrompt = "",
const std::string optDelayed = "" )

Plot the induced signal.

◆ PlotTransferFunction()

void Garfield::Sensor::PlotTransferFunction ( )

Plot the presently set transfer function.

◆ PrintTransferFunction()

void Garfield::Sensor::PrintTransferFunction ( )

Print some information about the presently set transfer function.

◆ SetArea() [1/2]

bool Garfield::Sensor::SetArea ( const bool verbose = false)

Set the user area to the default.

◆ SetArea() [2/2]

bool Garfield::Sensor::SetArea ( const double xmin,
const double ymin,
const double zmin,
const double xmax,
const double ymax,
const double zmax )

Set the user area explicitly.

◆ SetDelayedSignalAveragingOrder()

void Garfield::Sensor::SetDelayedSignalAveragingOrder ( const unsigned int navg)
inline

Set the number of points to be used when averaging the delayed signal vector over a time bin (default: 0).

The averaging is done with a $2\times navg + 1$ point Newton-Raphson integration.

Definition at line 129 of file Sensor.hh.

129 {
130 m_nAvgDelayedSignal = navg;
131 }
unsigned int m_nAvgDelayedSignal
Definition Sensor.hh:339

◆ SetDelayedSignalTimes()

void Garfield::Sensor::SetDelayedSignalTimes ( const std::vector< double > & ts)

Set the points in time at which to evaluate the delayed weighting field.

◆ SetNoiseFunction()

void Garfield::Sensor::SetNoiseFunction ( double(* )(double t))

Set the function to be used for evaluating the noise component.

◆ SetSignal() [1/2]

void Garfield::Sensor::SetSignal ( const std::string & label,
const std::vector< double > & ts,
const std::vector< double > & is )

Set/override the signal.

◆ SetSignal() [2/2]

void Garfield::Sensor::SetSignal ( const std::string & label,
const unsigned int bin,
const double signal )

Set/override the signal in a given time bin explicitly.

◆ SetTimeWindow()

void Garfield::Sensor::SetTimeWindow ( const double tstart,
const double tstep,
const unsigned int nsteps )

Set the time window and binning for the signal calculation.

Parameters
tstartstart time [ns]
tstepbin width [ns]
nstepsnumber of bins

◆ SetTransferFunction() [1/3]

void Garfield::Sensor::SetTransferFunction ( const std::vector< double > & times,
const std::vector< double > & values )

Set the points to be used for interpolating the transfer function.

◆ SetTransferFunction() [2/3]

void Garfield::Sensor::SetTransferFunction ( Shaper & shaper)

Set the transfer function using a Shaper object.

◆ SetTransferFunction() [3/3]

void Garfield::Sensor::SetTransferFunction ( std::function< double(double)> )

Set the function to be used for evaluating the transfer function.

◆ StepSizeHint()

double Garfield::Sensor::StepSizeHint ( )

◆ TransferFunctionSq()

double Garfield::Sensor::TransferFunctionSq ( )
private

◆ WeightingField()

void Garfield::Sensor::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).

◆ WeightingPotential()

double Garfield::Sensor::WeightingPotential ( const double x,
const double y,
const double z,
const std::string & label )

Get the weighting potential at (x, y, z).

Member Data Documentation

◆ m_cacheTransferFunction

bool Garfield::Sensor::m_cacheTransferFunction = true
private

Definition at line 345 of file Sensor.hh.

◆ m_className

std::string Garfield::Sensor::m_className = "Sensor"
private

Definition at line 309 of file Sensor.hh.

◆ m_components

std::vector<std::tuple<Component*, bool, bool> > Garfield::Sensor::m_components
private

Components.

Definition at line 314 of file Sensor.hh.

◆ m_debug

bool Garfield::Sensor::m_debug = false
private

Definition at line 363 of file Sensor.hh.

◆ m_delayedSignal

bool Garfield::Sensor::m_delayedSignal = false
private

Definition at line 337 of file Sensor.hh.

◆ m_delayedSignalTimes

std::vector<double> Garfield::Sensor::m_delayedSignalTimes
private

Definition at line 338 of file Sensor.hh.

◆ m_electrodes

std::vector<Electrode> Garfield::Sensor::m_electrodes
private

Electrodes.

Definition at line 329 of file Sensor.hh.

◆ m_fNoise

double(* Garfield::Sensor::m_fNoise) (double t) = nullptr
private

Definition at line 352 of file Sensor.hh.

◆ m_fTransfer

std::function<double(double)> Garfield::Sensor::m_fTransfer
private

Definition at line 342 of file Sensor.hh.

◆ m_fTransferFFT

std::vector<double> Garfield::Sensor::m_fTransferFFT
private

Definition at line 349 of file Sensor.hh.

◆ m_fTransferSq

double Garfield::Sensor::m_fTransferSq = -1.
private

Definition at line 347 of file Sensor.hh.

◆ m_fTransferTab

std::vector<std::pair<double, double> > Garfield::Sensor::m_fTransferTab
private

Definition at line 344 of file Sensor.hh.

◆ m_hasUserArea

bool Garfield::Sensor::m_hasUserArea = false
private

Definition at line 360 of file Sensor.hh.

◆ m_mutex

std::mutex Garfield::Sensor::m_mutex
private

Mutex.

Definition at line 311 of file Sensor.hh.

◆ m_nAvgDelayedSignal

unsigned int Garfield::Sensor::m_nAvgDelayedSignal = 0
private

Definition at line 339 of file Sensor.hh.

◆ m_nEvents

unsigned int Garfield::Sensor::m_nEvents = 0
private

Definition at line 335 of file Sensor.hh.

◆ m_nTimeBins

unsigned int Garfield::Sensor::m_nTimeBins = 200
private

Definition at line 334 of file Sensor.hh.

◆ m_shaper

Shaper* Garfield::Sensor::m_shaper = nullptr
private

Definition at line 343 of file Sensor.hh.

◆ m_thresholdCrossings

std::vector<std::pair<double, bool> > Garfield::Sensor::m_thresholdCrossings
private

Definition at line 354 of file Sensor.hh.

◆ m_thresholdLevel

double Garfield::Sensor::m_thresholdLevel = 0.
private

Definition at line 355 of file Sensor.hh.

◆ m_tStart

double Garfield::Sensor::m_tStart = 0.
private

Definition at line 332 of file Sensor.hh.

◆ m_tStep

double Garfield::Sensor::m_tStep = 10.
private

Definition at line 333 of file Sensor.hh.

◆ m_xMaxUser

double Garfield::Sensor::m_xMaxUser = 0.
private

Definition at line 359 of file Sensor.hh.

◆ m_xMinUser

double Garfield::Sensor::m_xMinUser = 0.
private

Definition at line 358 of file Sensor.hh.

◆ m_yMaxUser

double Garfield::Sensor::m_yMaxUser = 0.
private

Definition at line 359 of file Sensor.hh.

◆ m_yMinUser

double Garfield::Sensor::m_yMinUser = 0.
private

Definition at line 358 of file Sensor.hh.

◆ m_zMaxUser

double Garfield::Sensor::m_zMaxUser = 0.
private

Definition at line 359 of file Sensor.hh.

◆ m_zMinUser

double Garfield::Sensor::m_zMinUser = 0.
private

Definition at line 358 of file Sensor.hh.


The documentation for this class was generated from the following file:
  • /builds/garfield/docs/source/Include/Garfield/Sensor.hh