Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
|
Calculate electron drift lines and avalanches using microscopic tracking. More...
#include <AvalancheMicroscopic.hh>
Classes | |
struct | Electron |
struct | Point |
Public Member Functions | |
AvalancheMicroscopic () | |
Default constructor. | |
AvalancheMicroscopic (Sensor *sensor) | |
Constructor. | |
~AvalancheMicroscopic () | |
Destructor. | |
void | SetSensor (Sensor *sensor) |
Set the sensor. | |
void | EnablePlotting (ViewDrift *view, const size_t nColl=100) |
Switch on drift line plotting. | |
void | DisablePlotting () |
Switch off drift line plotting. | |
void | EnableExcitationMarkers (const bool on=true) |
Draw a marker at every excitation or not. | |
void | EnableIonisationMarkers (const bool on=true) |
Draw a marker at every ionising collision or not. | |
void | EnableAttachmentMarkers (const bool on=true) |
Draw a marker at every attachment or not. | |
void | EnableSignalCalculation (const bool on=true) |
Switch calculation of induced currents on or off (default: enabled). | |
void | UseWeightingPotential (const bool on=true) |
void | EnableWeightingFieldIntegration (const bool on=true) |
void | UseInducedCharge (const bool on=true) |
Switch on calculation of the total induced charge (default: off). | |
void | EnablePathLengthComputation (const bool on=true) |
Compute and store the path length of each trajectory (default: off). | |
void | EnableElectronEnergyHistogramming (TH1 *histo) |
Fill a histogram with the electron energy distribution. | |
void | DisableElectronEnergyHistogramming () |
Stop histogramming the electron energy distribution. | |
void | EnableHoleEnergyHistogramming (TH1 *histo) |
Fill a histogram with the hole energy distribution. | |
void | DisableHoleEnergyHistogramming () |
Stop histogramming the hole energy distribution. | |
void | SetDistanceHistogram (TH1 *histo, const char opt='r') |
void | EnableDistanceHistogramming (const int type) |
Fill distance distribution histograms for a given collision type. | |
void | DisableDistanceHistogramming (const int type) |
Stop filling distance distribution histograms for a given collision type. | |
void | DisableDistanceHistogramming () |
Stop filling distance distribution histograms. | |
void | EnableSecondaryEnergyHistogramming (TH1 *histo) |
Fill histograms of the energy of electrons emitted in ionising collisions. | |
void | DisableSecondaryEnergyHistogramming () |
Stop histogramming the secondary electron energy distribution. | |
void | EnableDriftLines (const bool on=true) |
Switch on storage of drift lines (default: off). | |
void | EnablePhotonTransport (const bool on=true) |
void | EnableBandStructure (const bool on=true) |
Switch on stepping according to band structure E(k), for semiconductors. | |
void | EnableNullCollisionSteps (const bool on=true) |
Switch on update of coordinates for null-collision steps (default: off). | |
void | SetElectronTransportCut (const double cut) |
double | GetElectronTransportCut () const |
Retrieve the value of the energy threshold. | |
void | SetPhotonTransportCut (const double cut) |
Set an energy threshold for photon transport. | |
double | GetPhotonTransportCut () const |
Retrieve the energy threshold for transporting photons. | |
void | EnableAvalancheSizeLimit (const unsigned int size) |
void | DisableAvalancheSizeLimit () |
Do not apply a limit on the avalanche size. | |
int | GetAvalancheSizeLimit () const |
Retrieve the currently set size limit. | |
void | EnableMagneticField (const bool on=true) |
Switch on/off using the magnetic field in the stepping algorithm. | |
void | SetCollisionSteps (const unsigned int n) |
Set number of collisions to be skipped for storing drift lines. | |
void | SetTimeWindow (const double t0, const double t1) |
Define a time interval (only carriers inside the interval are simulated). | |
void | UnsetTimeWindow () |
Do not restrict the time interval within which carriers are simulated. | |
void | GetAvalancheSize (int &ne, int &ni) const |
Return the number of electrons and ions in the avalanche. | |
void | GetAvalancheSize (int &ne, int &nh, int &ni) const |
const std::vector< Electron > & | GetElectrons () const |
const std::vector< Electron > & | GetHoles () const |
size_t | GetNumberOfElectronEndpoints () const |
void | GetElectronEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, int &status) const |
void | GetElectronEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, double &dx1, double &dy1, double &dz1, int &status) const |
size_t | GetNumberOfElectronDriftLinePoints (const size_t i=0) const |
size_t | GetNumberOfHoleDriftLinePoints (const size_t i=0) const |
void | GetElectronDriftLinePoint (double &x, double &y, double &z, double &t, const size_t ip, const size_t ie=0) const |
void | GetHoleDriftLinePoint (double &x, double &y, double &z, double &t, const size_t ip, const size_t ih=0) const |
size_t | GetNumberOfHoleEndpoints () const |
void | GetHoleEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, int &status) const |
size_t | GetNumberOfPhotons () const |
void | GetPhoton (const size_t i, double &e, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const |
bool | DriftElectron (const double x, const double y, const double z, const double t, const double e, const double dx=0., const double dy=0., const double dz=0.) |
bool | AvalancheElectron (const double x, const double y, const double z, const double t, const double e, const double dx=0., const double dy=0., const double dz=0.) |
Calculate an avalanche initiated by a given electron. | |
void | AddElectron (const double x, const double y, const double z, const double t, const double e, const double dx=0., const double dy=0., const double dz=0.) |
Add an electron to the list of particles to be transported. | |
bool | ResumeAvalanche () |
Continue the avalanche simulation from the current set of electrons. | |
void | SetUserHandleStep (void(*f)(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole)) |
Set a callback function to be called at every step. | |
void | UnsetUserHandleStep () |
Deactivate the user handle called at every step. | |
void | SetUserHandleCollision (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m, double e0, double e1, double dx0, double dy0, double dz0, double dx1, double dy1, double dz1)) |
Set a callback function to be called at every (real) collision. | |
void | UnsetUserHandleCollision () |
Deactivate the user handle called at every collision. | |
void | SetUserHandleAttachment (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m)) |
Set a user handling procedure, to be called at every attachment. | |
void | UnsetUserHandleAttachment () |
Deactivate the user handle called at every attachment. | |
void | SetUserHandleInelastic (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m)) |
Set a user handling procedure, to be called at every inelastic collision. | |
void | UnsetUserHandleInelastic () |
Deactivate the user handle called at every inelastic collision. | |
void | SetUserHandleIonisation (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m)) |
void | UnsetUserHandleIonisation () |
Deactivate the user handle called at every ionisation. | |
void | EnableDebugging () |
Switch on debugging messages. | |
void | DisableDebugging () |
Calculate electron drift lines and avalanches using microscopic tracking.
Definition at line 17 of file AvalancheMicroscopic.hh.
|
inline |
Default constructor.
Definition at line 20 of file AvalancheMicroscopic.hh.
Referenced by AvalancheMicroscopic().
Garfield::AvalancheMicroscopic::AvalancheMicroscopic | ( | Sensor * | sensor | ) |
Constructor.
Definition at line 125 of file AvalancheMicroscopic.cc.
|
inline |
void Garfield::AvalancheMicroscopic::AddElectron | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t, | ||
const double | e, | ||
const double | dx = 0., | ||
const double | dy = 0., | ||
const double | dz = 0. ) |
Add an electron to the list of particles to be transported.
Definition at line 441 of file AvalancheMicroscopic.cc.
bool Garfield::AvalancheMicroscopic::AvalancheElectron | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t, | ||
const double | e, | ||
const double | dx = 0., | ||
const double | dy = 0., | ||
const double | dz = 0. ) |
Calculate an avalanche initiated by a given electron.
Definition at line 431 of file AvalancheMicroscopic.cc.
Referenced by GarfieldPhysics::DoIt(), and main().
|
inline |
Do not apply a limit on the avalanche size.
Definition at line 119 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 261 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming | ( | ) |
Stop filling distance distribution histograms.
Definition at line 230 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming | ( | const int | type | ) |
Stop filling distance distribution histograms for a given collision type.
Definition at line 216 of file AvalancheMicroscopic.cc.
|
inline |
Stop histogramming the electron energy distribution.
Definition at line 63 of file AvalancheMicroscopic.hh.
|
inline |
Stop histogramming the hole energy distribution.
Definition at line 67 of file AvalancheMicroscopic.hh.
|
inline |
Switch off drift line plotting.
Definition at line 32 of file AvalancheMicroscopic.hh.
|
inline |
Stop histogramming the secondary electron energy distribution.
Definition at line 85 of file AvalancheMicroscopic.hh.
bool Garfield::AvalancheMicroscopic::DriftElectron | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t, | ||
const double | e, | ||
const double | dx = 0., | ||
const double | dy = 0., | ||
const double | dz = 0. ) |
Calculate an electron drift line.
x,y,z,t | starting point of the electron |
e | initial energy of the electron |
dx,dy,dz | initial direction vector of the electron If the initial direction is not specified, it is sampled randomly. Secondary electrons are not transported. |
Definition at line 422 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every attachment or not.
Definition at line 38 of file AvalancheMicroscopic.hh.
|
inline |
Set a max. avalanche size (i. e. ignore ionising collisions once this size has been reached).
Definition at line 117 of file AvalancheMicroscopic.hh.
|
inline |
Switch on stepping according to band structure E(k), for semiconductors.
Definition at line 95 of file AvalancheMicroscopic.hh.
|
inline |
Switch on debugging messages.
Definition at line 260 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnableDistanceHistogramming | ( | const int | type | ) |
Fill distance distribution histograms for a given collision type.
Definition at line 194 of file AvalancheMicroscopic.cc.
|
inline |
Switch on storage of drift lines (default: off).
Definition at line 88 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnableElectronEnergyHistogramming | ( | TH1 * | histo | ) |
Fill a histogram with the electron energy distribution.
Definition at line 150 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every excitation or not.
Definition at line 34 of file AvalancheMicroscopic.hh.
Referenced by main().
void Garfield::AvalancheMicroscopic::EnableHoleEnergyHistogramming | ( | TH1 * | histo | ) |
Fill a histogram with the hole energy distribution.
Definition at line 160 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every ionising collision or not.
Definition at line 36 of file AvalancheMicroscopic.hh.
Referenced by main().
|
inline |
Switch on/off using the magnetic field in the stepping algorithm.
Definition at line 124 of file AvalancheMicroscopic.hh.
|
inline |
Switch on update of coordinates for null-collision steps (default: off).
Definition at line 100 of file AvalancheMicroscopic.hh.
|
inline |
Compute and store the path length of each trajectory (default: off).
Definition at line 57 of file AvalancheMicroscopic.hh.
|
inline |
Switch on photon transport.
Definition at line 92 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnablePlotting | ( | ViewDrift * | view, |
const size_t | nColl = 100 ) |
Switch on drift line plotting.
Definition at line 140 of file AvalancheMicroscopic.cc.
Referenced by main().
void Garfield::AvalancheMicroscopic::EnableSecondaryEnergyHistogramming | ( | TH1 * | histo | ) |
Fill histograms of the energy of electrons emitted in ionising collisions.
Definition at line 235 of file AvalancheMicroscopic.cc.
|
inline |
Switch calculation of induced currents on or off (default: enabled).
Definition at line 41 of file AvalancheMicroscopic.hh.
|
inline |
Integrate the weighting field over a drift line step when calculating the induced current (default: off).
Definition at line 49 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 142 of file AvalancheMicroscopic.hh.
|
inline |
Return the number of electrons and ions in the avalanche.
Definition at line 138 of file AvalancheMicroscopic.hh.
Referenced by GarfieldPhysics::DoIt().
|
inline |
Retrieve the currently set size limit.
Definition at line 121 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::GetElectronDriftLinePoint | ( | double & | x, |
double & | y, | ||
double & | z, | ||
double & | t, | ||
const size_t | ip, | ||
const size_t | ie = 0 ) const |
Definition at line 331 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::GetElectronEndpoint | ( | const size_t | i, |
double & | x0, | ||
double & | y0, | ||
double & | z0, | ||
double & | t0, | ||
double & | e0, | ||
double & | x1, | ||
double & | y1, | ||
double & | z1, | ||
double & | t1, | ||
double & | e1, | ||
double & | dx1, | ||
double & | dy1, | ||
double & | dz1, | ||
int & | status ) const |
void Garfield::AvalancheMicroscopic::GetElectronEndpoint | ( | const size_t | i, |
double & | x0, | ||
double & | y0, | ||
double & | z0, | ||
double & | t0, | ||
double & | e0, | ||
double & | x1, | ||
double & | y1, | ||
double & | z1, | ||
double & | t1, | ||
double & | e1, | ||
int & | status ) const |
Return the coordinates and time of start and end point of a given electron drift line.
i | index of the drift line |
x0,y0,z0,t0 | coordinates and time of the starting point |
x1,y1,z1,t1 | coordinates and time of the end point |
e0,e1 | initial and final energy |
status | status code (see GarfieldConstants.hh) |
Definition at line 257 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 162 of file AvalancheMicroscopic.hh.
Referenced by Garfield::AvalancheGrid::ImportElectronsFromAvalancheMicroscopic(), and main().
|
inline |
Retrieve the value of the energy threshold.
Definition at line 108 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::GetHoleDriftLinePoint | ( | double & | x, |
double & | y, | ||
double & | z, | ||
double & | t, | ||
const size_t | ip, | ||
const size_t | ih = 0 ) const |
Definition at line 350 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::GetHoleEndpoint | ( | const size_t | i, |
double & | x0, | ||
double & | y0, | ||
double & | z0, | ||
double & | t0, | ||
double & | e0, | ||
double & | x1, | ||
double & | y1, | ||
double & | z1, | ||
double & | t1, | ||
double & | e1, | ||
int & | status ) const |
Definition at line 284 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 163 of file AvalancheMicroscopic.hh.
size_t Garfield::AvalancheMicroscopic::GetNumberOfElectronDriftLinePoints | ( | const size_t | i = 0 | ) | const |
Definition at line 311 of file AvalancheMicroscopic.cc.
|
inline |
Return the number of electron trajectories in the last simulated avalanche (including captured electrons).
Definition at line 166 of file AvalancheMicroscopic.hh.
Referenced by main().
size_t Garfield::AvalancheMicroscopic::GetNumberOfHoleDriftLinePoints | ( | const size_t | i = 0 | ) | const |
Definition at line 321 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 193 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 198 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::GetPhoton | ( | const size_t | i, |
double & | e, | ||
double & | x0, | ||
double & | y0, | ||
double & | z0, | ||
double & | t0, | ||
double & | x1, | ||
double & | y1, | ||
double & | z1, | ||
double & | t1, | ||
int & | status ) const |
Definition at line 370 of file AvalancheMicroscopic.cc.
|
inline |
Retrieve the energy threshold for transporting photons.
Definition at line 113 of file AvalancheMicroscopic.hh.
bool Garfield::AvalancheMicroscopic::ResumeAvalanche | ( | ) |
Continue the avalanche simulation from the current set of electrons.
Definition at line 451 of file AvalancheMicroscopic.cc.
|
inline |
Set number of collisions to be skipped for storing drift lines.
Definition at line 130 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::SetDistanceHistogram | ( | TH1 * | histo, |
const char | opt = 'r' ) |
Fill histograms of the distance between successive collisions.
histo | pointer to the histogram to be filled |
opt | direction ('x', 'y', 'z', 'r') |
Definition at line 170 of file AvalancheMicroscopic.cc.
|
inline |
Set a (lower) energy threshold for electron transport. This can be useful for simulating delta electrons.
Definition at line 106 of file AvalancheMicroscopic.hh.
|
inline |
Set an energy threshold for photon transport.
Definition at line 111 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::SetSensor | ( | Sensor * | sensor | ) |
Set the sensor.
Definition at line 132 of file AvalancheMicroscopic.cc.
Referenced by GarfieldPhysics::DoIt().
void Garfield::AvalancheMicroscopic::SetTimeWindow | ( | const double | t0, |
const double | t1 ) |
Define a time interval (only carriers inside the interval are simulated).
Definition at line 245 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleAttachment | ( | void(* | f )(double x, double y, double z, double t, int type, int level, Medium *m) | ) |
Set a user handling procedure, to be called at every attachment.
Definition at line 407 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleCollision | ( | void(* | f )(double x, double y, double z, double t, int type, int level, Medium *m, double e0, double e1, double dx0, double dy0, double dz0, double dx1, double dy1, double dz1) | ) |
Set a callback function to be called at every (real) collision.
Definition at line 400 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleInelastic | ( | void(* | f )(double x, double y, double z, double t, int type, int level, Medium *m) | ) |
Set a user handling procedure, to be called at every inelastic collision.
Definition at line 412 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleIonisation | ( | void(* | f )(double x, double y, double z, double t, int type, int level, Medium *m) | ) |
Set a user handling procedure, to be called at every ionising collision or excitation followed by Penning transfer.
Definition at line 417 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleStep | ( | void(* | f )(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole) | ) |
Set a callback function to be called at every step.
Definition at line 390 of file AvalancheMicroscopic.cc.
|
inline |
Do not restrict the time interval within which carriers are simulated.
Definition at line 135 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every attachment.
Definition at line 246 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every collision.
Definition at line 241 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every inelastic collision.
Definition at line 251 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every ionisation.
Definition at line 257 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every step.
Definition at line 233 of file AvalancheMicroscopic.hh.
|
inline |
Switch on calculation of the total induced charge (default: off).
Definition at line 54 of file AvalancheMicroscopic.hh.
|
inline |
Use the weighting potential (as opposed to the weighting field) for calculating the induced current.
Definition at line 44 of file AvalancheMicroscopic.hh.