Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
|
#include <AvalancheMC.hh>
Classes | |
struct | EndPoint |
struct | Point |
Public Member Functions | |
AvalancheMC () | |
Default constructor. | |
AvalancheMC (Sensor *sensor) | |
Constructor. | |
~AvalancheMC () | |
Destructor. | |
void | SetSensor (Sensor *s) |
Set the sensor. | |
bool | DriftElectron (const double x, const double y, const double z, const double t) |
Simulate the drift line of an electron from a given starting point. | |
bool | DriftHole (const double x, const double y, const double z, const double t) |
Simulate the drift line of a hole from a given starting point. | |
bool | DriftIon (const double x, const double y, const double z, const double t) |
Simulate the drift line of an ion from a given starting point. | |
bool | DriftNegativeIon (const double x, const double y, const double z, const double t) |
Simulate the drift line of a negative ion from a given starting point. | |
bool | AvalancheElectron (const double x, const double y, const double z, const double t, const bool hole=false) |
bool | AvalancheHole (const double x, const double y, const double z, const double t, const bool electron=false) |
Simulate an avalanche initiated by a hole at a given starting point. | |
bool | AvalancheElectronHole (const double x, const double y, const double z, const double t) |
Simulate an avalanche initiated by an electron-hole pair. | |
void | AddElectron (const double x, const double y, const double z, const double t) |
Add an electron to the list of particles to be transported. | |
void | AddHole (const double x, const double y, const double z, const double t) |
Add a hole to the list of particles to be transported. | |
void | AddIon (const double x, const double y, const double z, const double t) |
Add an ion to the list of particles to be transported. | |
bool | ResumeAvalanche (const bool electron=true, const bool hole=true) |
Resume the simulation from the current set of charge carriers. | |
const std::vector< EndPoint > & | GetElectrons () const |
const std::vector< EndPoint > & | GetHoles () const |
const std::vector< EndPoint > & | GetIons () const |
const std::vector< EndPoint > & | GetNegativeIons () const |
size_t | GetNumberOfElectronEndpoints () const |
size_t | GetNumberOfHoleEndpoints () const |
size_t | GetNumberOfIonEndpoints () const |
Return the number of ion trajectories. | |
void | GetElectronEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const |
void | GetHoleEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const |
void | GetIonEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const |
void | EnablePlotting (ViewDrift *view) |
Switch on drift line plotting. | |
void | DisablePlotting () |
Switch off drift line plotting. | |
void | EnableDriftLines (const bool on=true) |
Switch on storage of drift lines (default: off). | |
void | EnableSignalCalculation (const bool on=true) |
Switch calculation of induced currents on or off (default: enabled). | |
void | SetSignalAveragingOrder (const unsigned int navg) |
void | UseWeightingPotential (const bool on=true) |
void | EnableInducedChargeCalculation (const bool on=true) |
Switch on calculation of induced charge (default: disabled). | |
void | EnableRKFSteps (const bool on=true) |
void | EnableProjectedPathIntegration (const bool on=true) |
void | EnableDiffusion () |
Switch on diffusion (default: enabled). | |
void | DisableDiffusion () |
Switch off diffusion. | |
void | EnableAttachment () |
void | DisableAttachment () |
Switch off attachment and multiplication. | |
void | EnableTownsendMap (const bool on=true) |
Retrieve the Townsend coefficient from the component. | |
void | EnableAttachmentMap (const bool on=true) |
Retrieve the attachment coefficient from the component. | |
void | EnableVelocityMap (const bool on=true) |
Retrieve the drift velocity from the component. | |
void | EnableAvalancheSizeLimit (const unsigned int size) |
void | DisableAvalancheSizeLimit () |
Do not limit the maximum avalanche size. | |
int | GetAvalancheSizeLimit () const |
Return the currently set avalanche size limit. | |
void | SetTimeSteps (const double d=0.02) |
Use fixed-time steps (default 20 ps). | |
void | SetDistanceSteps (const double d=0.001) |
Use fixed distance steps (default 10 um). | |
void | SetCollisionSteps (const unsigned int n=100) |
void | SetStepDistanceFunction (double(*f)(double x, double y, double z)) |
Retrieve the step distance from a user-supplied function. | |
void | SetTimeWindow (const double t0, const double t1) |
Define a time interval (only carriers inside the interval are drifted). | |
void | UnsetTimeWindow () |
Do not limit the time interval within which carriers are drifted. | |
void | SetElectronSignalScalingFactor (const double scale) |
Set multiplication factor for the signal induced by electrons. | |
void | SetHoleSignalScalingFactor (const double scale) |
Set multiplication factor for the signal induced by holes. | |
void | SetIonSignalScalingFactor (const double scale) |
Set multiplication factor for the signal induced by ions. | |
void | GetAvalancheSize (unsigned int &ne, unsigned int &ni) const |
Return the number of electrons and ions/holes in the avalanche. | |
void | EnableDebugging (const bool on=true) |
Switch debugging messages on/off (default: off). | |
Calculate drift lines and avalanches based on macroscopic transport coefficients, using Monte Carlo integration.
Definition at line 18 of file AvalancheMC.hh.
|
inline |
Garfield::AvalancheMC::AvalancheMC | ( | Sensor * | sensor | ) |
|
inline |
void Garfield::AvalancheMC::AddElectron | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Add an electron to the list of particles to be transported.
Definition at line 554 of file AvalancheMC.cc.
void Garfield::AvalancheMC::AddHole | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Add a hole to the list of particles to be transported.
Definition at line 565 of file AvalancheMC.cc.
void Garfield::AvalancheMC::AddIon | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Add an ion to the list of particles to be transported.
Definition at line 574 of file AvalancheMC.cc.
bool Garfield::AvalancheMC::AvalancheElectron | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t, | ||
const bool | hole = false ) |
Simulate an avalanche initiated by an electron at a given starting point.
x,y,z,t | coordinates and time of the initial electron |
hole | simulate the hole component of the avalanche or not |
Definition at line 526 of file AvalancheMC.cc.
bool Garfield::AvalancheMC::AvalancheElectronHole | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Simulate an avalanche initiated by an electron-hole pair.
Definition at line 544 of file AvalancheMC.cc.
bool Garfield::AvalancheMC::AvalancheHole | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t, | ||
const bool | electron = false ) |
Simulate an avalanche initiated by a hole at a given starting point.
Definition at line 535 of file AvalancheMC.cc.
|
inline |
Switch off attachment and multiplication.
Definition at line 153 of file AvalancheMC.hh.
|
inline |
Do not limit the maximum avalanche size.
Definition at line 166 of file AvalancheMC.hh.
|
inline |
|
inline |
Switch off drift line plotting.
Definition at line 111 of file AvalancheMC.hh.
bool Garfield::AvalancheMC::DriftElectron | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Simulate the drift line of an electron from a given starting point.
Definition at line 212 of file AvalancheMC.cc.
Referenced by GarfieldPhysics::DoIt(), and main().
bool Garfield::AvalancheMC::DriftHole | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Simulate the drift line of a hole from a given starting point.
Definition at line 221 of file AvalancheMC.cc.
Referenced by main().
bool Garfield::AvalancheMC::DriftIon | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Simulate the drift line of an ion from a given starting point.
Definition at line 229 of file AvalancheMC.cc.
bool Garfield::AvalancheMC::DriftNegativeIon | ( | const double | x, |
const double | y, | ||
const double | z, | ||
const double | t ) |
Simulate the drift line of a negative ion from a given starting point.
Definition at line 237 of file AvalancheMC.cc.
|
inline |
Switch on attachment (and multiplication) for drift line calculation (default: enabled). For avalanches the flag is ignored.
Definition at line 151 of file AvalancheMC.hh.
|
inline |
Retrieve the attachment coefficient from the component.
Definition at line 158 of file AvalancheMC.hh.
Referenced by main().
|
inline |
Set a maximum avalanche size (ignore further multiplication once this size has been reached).
Definition at line 164 of file AvalancheMC.hh.
|
inline |
Switch debugging messages on/off (default: off).
Definition at line 199 of file AvalancheMC.hh.
|
inline |
Switch on diffusion (default: enabled).
Definition at line 145 of file AvalancheMC.hh.
|
inline |
Switch on storage of drift lines (default: off).
Definition at line 114 of file AvalancheMC.hh.
|
inline |
Switch on calculation of induced charge (default: disabled).
Definition at line 130 of file AvalancheMC.hh.
void Garfield::AvalancheMC::EnablePlotting | ( | ViewDrift * | view | ) |
Switch on drift line plotting.
Definition at line 77 of file AvalancheMC.cc.
Referenced by main().
|
inline |
Switch on equilibration of multiplication and attachment over the drift line (default: enabled).
Definition at line 140 of file AvalancheMC.hh.
|
inline |
Switch on Runge-Kutta-Fehlberg stepping (as opposed to simple straight-line steps.
Definition at line 136 of file AvalancheMC.hh.
|
inline |
Switch calculation of induced currents on or off (default: enabled).
Definition at line 117 of file AvalancheMC.hh.
|
inline |
Retrieve the Townsend coefficient from the component.
Definition at line 156 of file AvalancheMC.hh.
|
inline |
Retrieve the drift velocity from the component.
Definition at line 160 of file AvalancheMC.hh.
|
inline |
Return the number of electrons and ions/holes in the avalanche.
Definition at line 193 of file AvalancheMC.hh.
|
inline |
Return the currently set avalanche size limit.
Definition at line 168 of file AvalancheMC.hh.
void Garfield::AvalancheMC::GetElectronEndpoint | ( | const size_t | i, |
double & | x0, | ||
double & | y0, | ||
double & | z0, | ||
double & | t0, | ||
double & | x1, | ||
double & | y1, | ||
double & | z1, | ||
double & | t1, | ||
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 |
status | status code (see GarfieldConstants.hh) |
Definition at line 192 of file AvalancheMC.cc.
Referenced by GarfieldPhysics::DoIt().
|
inline |
Definition at line 75 of file AvalancheMC.hh.
void Garfield::AvalancheMC::GetHoleEndpoint | ( | const size_t | i, |
double & | x0, | ||
double & | y0, | ||
double & | z0, | ||
double & | t0, | ||
double & | x1, | ||
double & | y1, | ||
double & | z1, | ||
double & | t1, | ||
int & | status ) const |
Definition at line 153 of file AvalancheMC.cc.
|
inline |
Definition at line 76 of file AvalancheMC.hh.
void Garfield::AvalancheMC::GetIonEndpoint | ( | const size_t | i, |
double & | x0, | ||
double & | y0, | ||
double & | z0, | ||
double & | t0, | ||
double & | x1, | ||
double & | y1, | ||
double & | z1, | ||
double & | t1, | ||
int & | status ) const |
Definition at line 173 of file AvalancheMC.cc.
|
inline |
Definition at line 77 of file AvalancheMC.hh.
|
inline |
Definition at line 78 of file AvalancheMC.hh.
|
inline |
Return the number of electron trajectories in the last simulated avalanche (including captured electrons).
Definition at line 84 of file AvalancheMC.hh.
|
inline |
Return the number of hole trajectories in the last simulated avalanche (including captured holes).
Definition at line 87 of file AvalancheMC.hh.
|
inline |
Return the number of ion trajectories.
Definition at line 89 of file AvalancheMC.hh.
bool Garfield::AvalancheMC::ResumeAvalanche | ( | const bool | electron = true, |
const bool | hole = true ) |
Resume the simulation from the current set of charge carriers.
Definition at line 584 of file AvalancheMC.cc.
void Garfield::AvalancheMC::SetCollisionSteps | ( | const unsigned int | n = 100 | ) |
Use exponentially distributed time steps with mean equal to the specified multiple of the collision time (default model).
Definition at line 116 of file AvalancheMC.cc.
void Garfield::AvalancheMC::SetDistanceSteps | ( | const double | d = 0.001 | ) |
Use fixed distance steps (default 10 um).
Definition at line 101 of file AvalancheMC.cc.
Referenced by GarfieldPhysics::DoIt(), and main().
|
inline |
Set multiplication factor for the signal induced by electrons.
Definition at line 186 of file AvalancheMC.hh.
|
inline |
Set multiplication factor for the signal induced by holes.
Definition at line 188 of file AvalancheMC.hh.
|
inline |
Set multiplication factor for the signal induced by ions.
Definition at line 190 of file AvalancheMC.hh.
void Garfield::AvalancheMC::SetSensor | ( | Sensor * | s | ) |
Set the sensor.
Definition at line 68 of file AvalancheMC.cc.
Referenced by GarfieldPhysics::DoIt(), and main().
|
inline |
Set the number of points to be used when averaging the signal vector over a time bin in the Sensor class. The averaging is done with a
Definition at line 122 of file AvalancheMC.hh.
void Garfield::AvalancheMC::SetStepDistanceFunction | ( | double(* | f )(double x, double y, double z) | ) |
Retrieve the step distance from a user-supplied function.
Definition at line 131 of file AvalancheMC.cc.
void Garfield::AvalancheMC::SetTimeSteps | ( | const double | d = 0.02 | ) |
Use fixed-time steps (default 20 ps).
Definition at line 86 of file AvalancheMC.cc.
void Garfield::AvalancheMC::SetTimeWindow | ( | const double | t0, |
const double | t1 ) |
Define a time interval (only carriers inside the interval are drifted).
Definition at line 141 of file AvalancheMC.cc.
|
inline |
Do not limit the time interval within which carriers are drifted.
Definition at line 183 of file AvalancheMC.hh.
|
inline |
Use the weighting potential (as opposed to the weighting field) for calculating the induced signal.
Definition at line 125 of file AvalancheMC.hh.