![]() |
Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
|
#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. | |
Component * | GetComponent (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). | |
Medium * | GetMedium (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 > ×, 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< Electrode > | m_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 |
Shaper * | m_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 |
|
default |
Default constructor.
|
default |
Destructor.
Garfield::Sensor::Sensor | ( | Component * | comp | ) |
Constructor from a single component.
void Garfield::Sensor::AddComponent | ( | Component * | comp | ) |
Add a component.
void Garfield::Sensor::AddElectrode | ( | Component * | comp, |
const std::string & | label ) |
Add an electrode.
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.
void Garfield::Sensor::AddNoise | ( | const bool | total = true, |
const bool | electron = false, | ||
const bool | ion = false ) |
Add noise to the induced signal.
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.
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.
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.
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.
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.
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.
label | name of the electrode |
enc | Equivalent Noise Charge, in electrons. |
poisson | flag to sample the number of noise pulses from a Poisson distribution. Otherwise the noise charge in each bin is sampled from a Gaussian distribution. |
q0 | amplitude of the noise delta pulses (in electrons). |
void Garfield::Sensor::Clear | ( | ) |
Remove all components, electrodes and reset the sensor.
void Garfield::Sensor::ClearElectrodes | ( | ) |
Remove all electrodes.
void Garfield::Sensor::ClearSignal | ( | ) |
Reset signals and induced charges of all electrodes.
bool Garfield::Sensor::ComputeThresholdCrossings | ( | const double | thr, |
const std::string & | label, | ||
int & | n ) |
Determine the threshold crossings of the current signal curve.
thr | threshold value |
label | electrode for which to compute the threshold crossings |
n | number of threshold crossings |
bool Garfield::Sensor::ConvoluteSignal | ( | const std::string & | label, |
const bool | fft = false ) |
Convolute the induced current on a given electrode with the transfer function.
|
private |
|
private |
|
private |
|
private |
bool Garfield::Sensor::ConvoluteSignals | ( | const bool | fft = false | ) |
Convolute all induced currents with the transfer function.
double Garfield::Sensor::CreateGPUTransferObject | ( | SensorGPU *& | sensor_gpu | ) |
Create and initialise GPU Transfer class.
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.
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.
bool Garfield::Sensor::DelayAndSubtractFraction | ( | const double | td, |
const double | f ) |
Delay the signal and subtract an attenuated copy (modelling a constant fraction discriminator).
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).
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).
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).
void Garfield::Sensor::EnableComponent | ( | const unsigned int | i, |
const bool | on ) |
Activate/deactivate a given component.
|
inline |
|
inline |
Compute the component of the signal due to the delayed weighting field.
Definition at line 122 of file Sensor.hh.
void Garfield::Sensor::EnableMagneticField | ( | const unsigned int | i, |
const bool | on ) |
Activate/deactivate use of the magnetic field of a given component.
|
inline |
Cache integral and FFT of the transfer function instead of recomputing it at every call (default: on).
void Garfield::Sensor::ExportSignal | ( | const std::string & | label, |
const std::string & | filename, | ||
const bool | chargeCariers = false ) const |
Exporting induced signal to a csv file.
|
private |
|
inlineprivate |
Definition at line 372 of file Sensor.hh.
|
private |
bool Garfield::Sensor::GetArea | ( | double & | xmin, |
double & | ymin, | ||
double & | zmin, | ||
double & | xmax, | ||
double & | ymax, | ||
double & | zmax ) |
Return the current user area.
|
private |
Component * Garfield::Sensor::GetComponent | ( | const unsigned int | i | ) |
Retrieve the pointer to a given component.
double Garfield::Sensor::GetDelayedElectronSignal | ( | const std::string & | label, |
const unsigned int | bin ) |
Retrieve the delayed electron signal for a given electrode and time bin.
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.
double Garfield::Sensor::GetDelayedSignal | ( | const std::string & | label, |
const unsigned int | bin ) |
Retrieve the delayed signal for a given electrode and time bin.
double Garfield::Sensor::GetElectronSignal | ( | const std::string & | label, |
const unsigned int | bin ) |
Retrieve the electron signal for a given electrode and time bin.
double Garfield::Sensor::GetInducedCharge | ( | const std::string & | label | ) |
Calculated using the weighting potentials at the start and end points.
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.
Medium * Garfield::Sensor::GetMedium | ( | const double | x, |
const double | y, | ||
const double | z ) |
Get the medium at (x, y, z).
|
inline |
Get the number of components attached to the sensor.
Definition at line 34 of file Sensor.hh.
|
inline |
Get the number of electrodes attached to the sensor.
Definition at line 47 of file Sensor.hh.
|
inline |
Get the number of threshold crossings (after having called ComputeThresholdCrossings).
double Garfield::Sensor::GetPromptSignal | ( | const std::string & | label, |
const unsigned int | bin ) |
Retrieve the prompt signal for a given electrode and time bin.
double Garfield::Sensor::GetSignal | ( | const std::string & | label, |
const unsigned int | bin ) |
Retrieve the total signal for a given electrode and time bin.
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.
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.
i | index |
time | threshold crossing time [ns] |
level | threshold (should correspond to the value requested). |
rise | flag whether the crossing is on a rising or falling slope. |
|
inline |
Retrieve the time window and binning.
double Garfield::Sensor::GetTotalInducedCharge | ( | const std::string & | label | ) |
double Garfield::Sensor::GetTransferFunction | ( | const double | t | ) |
Evaluate the transfer function at a given time.
bool Garfield::Sensor::GetVoltageRange | ( | double & | vmin, |
double & | vmax ) |
Return the voltage range.
bool Garfield::Sensor::HasMagneticField | ( | ) | const |
Does the sensor have a non-zero magnetic field?
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).
bool Garfield::Sensor::IntegrateSignal | ( | const std::string & | label | ) |
Replace the signal on a given electrode by its integral.
|
private |
bool Garfield::Sensor::IntegrateSignals | ( | ) |
Replace the signals on all electrodes curve by their integrals.
|
private |
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.
bool Garfield::Sensor::IsInArea | ( | const double | x, |
const double | y, | ||
const double | z ) |
Check if a point is inside the user area.
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.
bool Garfield::Sensor::IsIntegrated | ( | const std::string & | label | ) | const |
Return whether the signal has been integrated/convoluted.
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).
|
private |
|
inline |
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.
void Garfield::Sensor::PlotTransferFunction | ( | ) |
Plot the presently set transfer function.
void Garfield::Sensor::PrintTransferFunction | ( | ) |
Print some information about the presently set transfer function.
bool Garfield::Sensor::SetArea | ( | const bool | verbose = false | ) |
Set the user area to the default.
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.
|
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
void Garfield::Sensor::SetDelayedSignalTimes | ( | const std::vector< double > & | ts | ) |
Set the points in time at which to evaluate the delayed weighting field.
void Garfield::Sensor::SetNoiseFunction | ( | double(* | f )(double t) | ) |
Set the function to be used for evaluating the noise component.
void Garfield::Sensor::SetSignal | ( | const std::string & | label, |
const std::vector< double > & | ts, | ||
const std::vector< double > & | is ) |
Set/override the signal.
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.
void Garfield::Sensor::SetTimeWindow | ( | const double | tstart, |
const double | tstep, | ||
const unsigned int | nsteps ) |
Set the time window and binning for the signal calculation.
tstart | start time [ns] |
tstep | bin width [ns] |
nsteps | number of bins |
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.
void Garfield::Sensor::SetTransferFunction | ( | Shaper & | shaper | ) |
Set the transfer function using a Shaper object.
void Garfield::Sensor::SetTransferFunction | ( | std::function< double(double)> | ) |
Set the function to be used for evaluating the transfer function.
double Garfield::Sensor::StepSizeHint | ( | ) |
|
private |
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).
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).
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |