Garfield++ 4.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>
Public Member Functions | |
AvalancheMicroscopic () | |
Constructor. | |
~AvalancheMicroscopic () | |
Destructor. | |
void | SetSensor (Sensor *sensor) |
Set the sensor. | |
void | EnablePlotting (ViewDrift *view) |
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 on calculation of induced currents (default: off). | |
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 | 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) |
Enable magnetic field in stepping algorithm (default: off). | |
void | SetCollisionSteps (const unsigned int n) |
Set number of collisions to be skipped for plotting. | |
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 |
unsigned int | GetNumberOfElectronEndpoints () const |
void | GetElectronEndpoint (const unsigned int 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 unsigned int 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 |
unsigned int | GetNumberOfElectronDriftLinePoints (const unsigned int i=0) const |
unsigned int | GetNumberOfHoleDriftLinePoints (const unsigned int i=0) const |
void | GetElectronDriftLinePoint (double &x, double &y, double &z, double &t, const int ip, const unsigned int iel=0) const |
void | GetHoleDriftLinePoint (double &x, double &y, double &z, double &t, const int ip, const unsigned int iel=0) const |
unsigned int | GetNumberOfHoleEndpoints () const |
void | GetHoleEndpoint (const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, int &status) const |
unsigned int | GetNumberOfPhotons () const |
void | GetPhoton (const unsigned int 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 x0, const double y0, const double z0, const double t0, const double e0, const double dx0=0., const double dy0=0., const double dz0=0.) |
bool | AvalancheElectron (const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0=0., const double dy0=0., const double dz0=0.) |
Calculate an avalanche initiated by a given electron. | |
bool | ResumeAvalanche () |
void | SetUserHandleStep (void(*f)(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole)) |
Set a user handling procedure. This function is 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 user handling procedure, 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.
Garfield::AvalancheMicroscopic::AvalancheMicroscopic | ( | ) |
Constructor.
Definition at line 91 of file AvalancheMicroscopic.cc.
|
inline |
bool Garfield::AvalancheMicroscopic::AvalancheElectron | ( | const double | x0, |
const double | y0, | ||
const double | z0, | ||
const double | t0, | ||
const double | e0, | ||
const double | dx0 = 0. , |
||
const double | dy0 = 0. , |
||
const double | dz0 = 0. |
||
) |
Calculate an avalanche initiated by a given electron.
Definition at line 462 of file AvalancheMicroscopic.cc.
Referenced by GarfieldPhysics::DoIt(), and main().
|
inline |
Do not apply a limit on the avalanche size.
Definition at line 113 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 235 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming | ( | ) |
Stop filling distance distribution histograms.
Definition at line 199 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming | ( | const int | type | ) |
Stop filling distance distribution histograms for a given collision type.
Definition at line 185 of file AvalancheMicroscopic.cc.
|
inline |
Stop histogramming the electron energy distribution.
Definition at line 57 of file AvalancheMicroscopic.hh.
|
inline |
Stop histogramming the hole energy distribution.
Definition at line 61 of file AvalancheMicroscopic.hh.
|
inline |
Switch off drift line plotting.
Definition at line 30 of file AvalancheMicroscopic.hh.
|
inline |
Stop histogramming the secondary electron energy distribution.
Definition at line 79 of file AvalancheMicroscopic.hh.
bool Garfield::AvalancheMicroscopic::DriftElectron | ( | const double | x0, |
const double | y0, | ||
const double | z0, | ||
const double | t0, | ||
const double | e0, | ||
const double | dx0 = 0. , |
||
const double | dy0 = 0. , |
||
const double | dz0 = 0. |
||
) |
Calculate an electron drift line.
x0,y0,z0,t0 | starting point of the electron |
e0 | initial energy of the electron |
dx0,dy0,dz0 | initial direction of the electron If the initial direction is not specified, it is sampled randomly. Secondary electrons are not transported. |
Definition at line 453 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every attachment or not.
Definition at line 36 of file AvalancheMicroscopic.hh.
|
inline |
Set a max. avalanche size (i. e. ignore ionising collisions once this size has been reached).
Definition at line 111 of file AvalancheMicroscopic.hh.
|
inline |
Switch on stepping according to band structure E(k), for semiconductors.
Definition at line 89 of file AvalancheMicroscopic.hh.
|
inline |
Switch on debugging messages.
Definition at line 234 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnableDistanceHistogramming | ( | const int | type | ) |
Fill distance distribution histograms for a given collision type.
Definition at line 163 of file AvalancheMicroscopic.cc.
|
inline |
Switch on storage of drift lines (default: off).
Definition at line 82 of file AvalancheMicroscopic.hh.
Referenced by EnablePlotting().
void Garfield::AvalancheMicroscopic::EnableElectronEnergyHistogramming | ( | TH1 * | histo | ) |
Fill a histogram with the electron energy distribution.
Definition at line 119 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every excitation or not.
Definition at line 32 of file AvalancheMicroscopic.hh.
Referenced by main().
void Garfield::AvalancheMicroscopic::EnableHoleEnergyHistogramming | ( | TH1 * | histo | ) |
Fill a histogram with the hole energy distribution.
Definition at line 129 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every ionising collision or not.
Definition at line 34 of file AvalancheMicroscopic.hh.
Referenced by main().
|
inline |
Enable magnetic field in stepping algorithm (default: off).
Definition at line 118 of file AvalancheMicroscopic.hh.
Referenced by main().
|
inline |
Switch on update of coordinates for null-collision steps (default: off).
Definition at line 94 of file AvalancheMicroscopic.hh.
|
inline |
Switch on photon transport.
Definition at line 86 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnablePlotting | ( | ViewDrift * | view | ) |
Switch on drift line plotting.
Definition at line 105 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 204 of file AvalancheMicroscopic.cc.
|
inline |
Switch on calculation of induced currents (default: off).
Definition at line 39 of file AvalancheMicroscopic.hh.
Referenced by main().
|
inline |
Integrate the weighting field over a drift line step when calculating the induced current (default: off).
Definition at line 47 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 133 of file AvalancheMicroscopic.hh.
|
inline |
Return the number of electrons and ions in the avalanche.
Definition at line 129 of file AvalancheMicroscopic.hh.
Referenced by GarfieldPhysics::DoIt().
|
inline |
Retrieve the currently set size limit.
Definition at line 115 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::GetElectronDriftLinePoint | ( | double & | x, |
double & | y, | ||
double & | z, | ||
double & | t, | ||
const int | ip, | ||
const unsigned int | iel = 0 |
||
) | const |
Definition at line 334 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::GetElectronEndpoint | ( | const unsigned int | 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 |
Definition at line 253 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::GetElectronEndpoint | ( | const unsigned int | 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 226 of file AvalancheMicroscopic.cc.
Referenced by Garfield::AvalancheGrid::GetElectronsFromAvalancheMicroscopic().
|
inline |
Retrieve the value of the energy threshold.
Definition at line 102 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::GetHoleDriftLinePoint | ( | double & | x, |
double & | y, | ||
double & | z, | ||
double & | t, | ||
const int | ip, | ||
const unsigned int | iel = 0 |
||
) | const |
Definition at line 366 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::GetHoleEndpoint | ( | const unsigned int | 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 282 of file AvalancheMicroscopic.cc.
unsigned int Garfield::AvalancheMicroscopic::GetNumberOfElectronDriftLinePoints | ( | const unsigned int | i = 0 | ) | const |
Definition at line 308 of file AvalancheMicroscopic.cc.
|
inline |
Return the number of electron trajectories in the last simulated avalanche (including captured electrons).
Definition at line 141 of file AvalancheMicroscopic.hh.
Referenced by Garfield::AvalancheGrid::GetElectronsFromAvalancheMicroscopic(), and main().
unsigned int Garfield::AvalancheMicroscopic::GetNumberOfHoleDriftLinePoints | ( | const unsigned int | i = 0 | ) | const |
Definition at line 321 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 170 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 177 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::GetPhoton | ( | const unsigned int | 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 399 of file AvalancheMicroscopic.cc.
|
inline |
Retrieve the energy threshold for transporting photons.
Definition at line 107 of file AvalancheMicroscopic.hh.
bool Garfield::AvalancheMicroscopic::ResumeAvalanche | ( | ) |
Definition at line 472 of file AvalancheMicroscopic.cc.
|
inline |
Set number of collisions to be skipped for plotting.
Definition at line 121 of file AvalancheMicroscopic.hh.
Referenced by main().
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 139 of file AvalancheMicroscopic.cc.
|
inline |
Set a (lower) energy threshold for electron transport. This can be useful for simulating delta electrons.
Definition at line 100 of file AvalancheMicroscopic.hh.
|
inline |
Set an energy threshold for photon transport.
Definition at line 105 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::SetSensor | ( | Sensor * | sensor | ) |
Set the sensor.
Definition at line 97 of file AvalancheMicroscopic.cc.
Referenced by GarfieldPhysics::DoIt(), and main().
void Garfield::AvalancheMicroscopic::SetTimeWindow | ( | const double | t0, |
const double | t1 | ||
) |
Define a time interval (only carriers inside the interval are simulated).
Definition at line 214 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleAttachment | ( | void(*)(double x, double y, double z, double t, int type, int level, Medium *m) | f | ) |
Set a user handling procedure, to be called at every attachment.
Definition at line 438 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleCollision | ( | void(*)(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) | f | ) |
Set a user handling procedure, to be called at every (real) collision.
Definition at line 431 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleInelastic | ( | void(*)(double x, double y, double z, double t, int type, int level, Medium *m) | f | ) |
Set a user handling procedure, to be called at every inelastic collision.
Definition at line 443 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleIonisation | ( | void(*)(double x, double y, double z, double t, int type, int level, Medium *m) | f | ) |
Set a user handling procedure, to be called at every ionising collision or excitation followed by Penning transfer.
Definition at line 448 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::SetUserHandleStep | ( | void(*)(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole) | f | ) |
Set a user handling procedure. This function is called at every step.
Definition at line 421 of file AvalancheMicroscopic.cc.
|
inline |
Do not restrict the time interval within which carriers are simulated.
Definition at line 126 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every attachment.
Definition at line 220 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every collision.
Definition at line 215 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every inelastic collision.
Definition at line 225 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every ionisation.
Definition at line 231 of file AvalancheMicroscopic.hh.
|
inline |
Deactivate the user handle called at every step.
Definition at line 207 of file AvalancheMicroscopic.hh.
|
inline |
Switch on calculation of the total induced charge (default: off).
Definition at line 52 of file AvalancheMicroscopic.hh.
|
inline |
Use the weighting potential (as opposed to the weighting field) for calculating the induced current.
Definition at line 42 of file AvalancheMicroscopic.hh.