Garfield++ v2r0
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 () |
void | EnableExcitationMarkers () |
Draw a marker at every excitation. | |
void | DisableExcitationMarkers () |
void | EnableIonisationMarkers () |
Draw a marker at every ionising collision. | |
void | DisableIonisationMarkers () |
void | EnableAttachmentMarkers () |
Draw a marker at every attachment. | |
void | DisableAttachmentMarkers () |
void | EnableSignalCalculation () |
Switch on calculation of induced currents. | |
void | DisableSignalCalculation () |
void | EnableInducedChargeCalculation () |
Switch on calculation of the total induced charge. | |
void | DisableInducedChargeCalculation () |
void | EnableElectronEnergyHistogramming (TH1 *histo) |
Switch on filling histograms for electron energy distribution. | |
void | DisableElectronEnergyHistogramming () |
void | EnableHoleEnergyHistogramming (TH1 *histo) |
void | DisableHoleEnergyHistogramming () |
void | SetDistanceHistogram (TH1 *histo, const char opt='r') |
void | EnableDistanceHistogramming (const int type) |
void | DisableDistanceHistogramming (const int type) |
Switch on distance distribution histograms for a given collision type. | |
void | DisableDistanceHistogramming () |
void | EnableSecondaryEnergyHistogramming (TH1 *histo) |
Fill histograms of the energy of electrons emitted in ionising collisions. | |
void | DisableSecondaryEnergyHistogramming () |
void | EnableDriftLines () |
Switch on storage of drift lines. | |
void | DisableDriftLines () |
void | EnablePhotonTransport () |
void | DisablePhotonTransport () |
void | EnableBandStructure () |
Switch on stepping according to band structure E(k), for semiconductors. | |
void | DisableBandStructure () |
void | EnableNullCollisionSteps () |
Switch on update of coordinates for null-collision steps (default: off). | |
void | DisableNullCollisionSteps () |
void | SetElectronTransportCut (const double cut) |
double | GetElectronTransportCut () const |
void | SetPhotonTransportCut (const double cut) |
Set an energy threshold for photon transport. | |
double | GetPhotonTransportCut () const |
void | EnableAvalancheSizeLimit (const unsigned int size) |
void | DisableAvalancheSizeLimit () |
int | GetAvalancheSizeLimit () const |
void | EnableMagneticField () |
Enable magnetic field in stepping algorithm (default: off). | |
void | DisableMagneticField () |
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 () |
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. | |
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 () |
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 () |
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 () |
void | 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. | |
void | UnsetUserHandleIonisation () |
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 25 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 517 of file AvalancheMicroscopic.cc.
Referenced by GarfieldPhysics::DoIt().
|
inline |
Definition at line 39 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 99 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 81 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 202 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming | ( | ) |
Definition at line 206 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming | ( | const int | type | ) |
Switch on distance distribution histograms for a given collision type.
Definition at line 186 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 72 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 51 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 33 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 53 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 47 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 36 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 104 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 85 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 77 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::DisablePlotting | ( | ) |
Definition at line 110 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 68 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 43 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 501 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every attachment.
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 98 of file AvalancheMicroscopic.hh.
|
inline |
Switch on stepping according to band structure E(k), for semiconductors.
Definition at line 80 of file AvalancheMicroscopic.hh.
|
inline |
Switch on debugging messages.
Definition at line 201 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnableDistanceHistogramming | ( | const int | type | ) |
Definition at line 163 of file AvalancheMicroscopic.cc.
|
inline |
Switch on storage of drift lines.
Definition at line 71 of file AvalancheMicroscopic.hh.
Referenced by EnablePlotting().
void Garfield::AvalancheMicroscopic::EnableElectronEnergyHistogramming | ( | TH1 * | histo | ) |
Switch on filling histograms for electron energy distribution.
Definition at line 116 of file AvalancheMicroscopic.cc.
|
inline |
Draw a marker at every excitation.
Definition at line 32 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnableHoleEnergyHistogramming | ( | TH1 * | histo | ) |
Definition at line 127 of file AvalancheMicroscopic.cc.
|
inline |
Switch on calculation of the total induced charge.
Definition at line 46 of file AvalancheMicroscopic.hh.
|
inline |
Draw a marker at every ionising collision.
Definition at line 35 of file AvalancheMicroscopic.hh.
|
inline |
Enable magnetic field in stepping algorithm (default: off).
Definition at line 103 of file AvalancheMicroscopic.hh.
|
inline |
Switch on update of coordinates for null-collision steps (default: off).
Definition at line 84 of file AvalancheMicroscopic.hh.
|
inline |
Switch on photon transport.
Definition at line 76 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::EnablePlotting | ( | ViewDrift * | view | ) |
Switch on drift line plotting.
Definition at line 94 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::EnableSecondaryEnergyHistogramming | ( | TH1 * | histo | ) |
Fill histograms of the energy of electrons emitted in ionising collisions.
Definition at line 212 of file AvalancheMicroscopic.cc.
|
inline |
Switch on calculation of induced currents.
Definition at line 42 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 118 of file AvalancheMicroscopic.hh.
|
inline |
Return the number of electrons and ions in the avalanche.
Definition at line 114 of file AvalancheMicroscopic.hh.
Referenced by GarfieldPhysics::DoIt().
|
inline |
Definition at line 100 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 353 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 267 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 238 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 90 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 387 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 298 of file AvalancheMicroscopic.cc.
unsigned int Garfield::AvalancheMicroscopic::GetNumberOfElectronDriftLinePoints | ( | const unsigned int | i = 0 | ) | const |
Definition at line 326 of file AvalancheMicroscopic.cc.
|
inline |
Return the number of electron trajectories in the last simulated avalanche (including captured electrons).
Definition at line 126 of file AvalancheMicroscopic.hh.
unsigned int Garfield::AvalancheMicroscopic::GetNumberOfHoleDriftLinePoints | ( | const unsigned int | i = 0 | ) | const |
Definition at line 340 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 152 of file AvalancheMicroscopic.hh.
|
inline |
Definition at line 159 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 421 of file AvalancheMicroscopic.cc.
|
inline |
Definition at line 94 of file AvalancheMicroscopic.hh.
|
inline |
Set number of collisions to be skipped for plotting.
Definition at line 107 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 138 of file AvalancheMicroscopic.cc.
|
inline |
Set a (lower) energy threshold for electron transport. This can be useful for simulating delta electrons.
Definition at line 89 of file AvalancheMicroscopic.hh.
|
inline |
Set an energy threshold for photon transport.
Definition at line 93 of file AvalancheMicroscopic.hh.
void Garfield::AvalancheMicroscopic::SetSensor | ( | Sensor * | sensor | ) |
Set the sensor.
Definition at line 85 of file AvalancheMicroscopic.cc.
Referenced by GarfieldPhysics::InitializePhysics().
void Garfield::AvalancheMicroscopic::SetTimeWindow | ( | const double | t0, |
const double | t1 | ||
) |
Define a time interval (only carriers inside the interval are simulated).
Definition at line 223 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 462 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 475 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.
Definition at line 488 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 443 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::UnsetTimeWindow | ( | ) |
Definition at line 236 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::UnsetUserHandleAttachment | ( | ) |
Definition at line 469 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::UnsetUserHandleInelastic | ( | ) |
Definition at line 482 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::UnsetUserHandleIonisation | ( | ) |
Definition at line 495 of file AvalancheMicroscopic.cc.
void Garfield::AvalancheMicroscopic::UnsetUserHandleStep | ( | ) |
Definition at line 456 of file AvalancheMicroscopic.cc.