Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::AvalancheMC Class Reference

#include <AvalancheMC.hh>

Public Member Functions

 AvalancheMC ()
 Constructor.
 
 ~AvalancheMC ()
 Destructor.
 
void SetSensor (Sensor *s)
 Set the sensor.
 
void EnablePlotting (ViewDrift *view)
 Switch on drift line plotting.
 
void DisablePlotting ()
 Switch off drift line plotting.
 
void EnableSignalCalculation (const bool on=true)
 Switch on calculation of induced currents (default: disabled).
 
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 EnableTcadTraps (const bool on=true)
 Switch on calculating trapping with TCAD traps.
 
void EnableTcadVelocity (const bool on=true)
 Switch on TCAD velocity maps.
 
void EnableMagneticField (const bool on=true)
 Enable use of magnetic field in stepping algorithm.
 
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.
 
unsigned int GetNumberOfDriftLinePoints () const
 Return the number of points along the last simulated drift line.
 
void GetDriftLinePoint (const unsigned int i, double &x, double &y, double &z, double &t) const
 Return the coordinates and time of a point along the last drift line.
 
unsigned int GetNumberOfElectronEndpoints () const
 
unsigned int GetNumberOfHoleEndpoints () const
 
unsigned int GetNumberOfIonEndpoints () const
 Return the number of ion trajectories.
 
void GetElectronEndpoint (const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
 
void GetHoleEndpoint (const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
 
void GetIonEndpoint (const unsigned int i, 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)
 Simulate the drift line of an electron from a given starting point.
 
bool DriftHole (const double x0, const double y0, const double z0, const double t0)
 Simulate the drift line of a hole from a given starting point.
 
bool DriftIon (const double x0, const double y0, const double z0, const double t0)
 Simulate the drift line of an ion from a given starting point.
 
bool AvalancheElectron (const double x0, const double y0, const double z0, const double t0, const bool hole=false)
 Simulate an avalanche initiated by an electron at a given starting point.
 
bool AvalancheHole (const double x0, const double y0, const double z0, const double t0, const bool electron=false)
 Simulate an avalanche initiated by a hole at a given starting point.
 
bool AvalancheElectronHole (const double x0, const double y0, const double z0, const double t0)
 
void EnableDebugging ()
 Switch on debugging messages (default: off).
 
void DisableDebugging ()
 Switch off debugging messages.
 

Detailed Description

Calculate drift lines and avalanches based on macroscopic transport coefficients, using Monte Carlo integration.

Definition at line 17 of file AvalancheMC.hh.

Constructor & Destructor Documentation

◆ AvalancheMC()

Garfield::AvalancheMC::AvalancheMC ( )

Constructor.

Definition at line 29 of file AvalancheMC.cc.

29{ m_drift.reserve(10000); }

◆ ~AvalancheMC()

Garfield::AvalancheMC::~AvalancheMC ( )
inline

Destructor.

Definition at line 22 of file AvalancheMC.hh.

22{}

Member Function Documentation

◆ AvalancheElectron()

bool Garfield::AvalancheMC::AvalancheElectron ( const double  x0,
const double  y0,
const double  z0,
const double  t0,
const bool  hole = false 
)

Simulate an avalanche initiated by an electron at a given starting point.

Definition at line 523 of file AvalancheMC.cc.

525 {
526 return Avalanche(x0, y0, z0, t0, 1, 0, 0, true, holes);
527}

◆ AvalancheElectronHole()

bool Garfield::AvalancheMC::AvalancheElectronHole ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Definition at line 535 of file AvalancheMC.cc.

536 {
537 return Avalanche(x0, y0, z0, t0, 1, 1, 0, true, true);
538}

◆ AvalancheHole()

bool Garfield::AvalancheMC::AvalancheHole ( const double  x0,
const double  y0,
const double  z0,
const double  t0,
const bool  electron = false 
)

Simulate an avalanche initiated by a hole at a given starting point.

Definition at line 529 of file AvalancheMC.cc.

531 {
532 return Avalanche(x0, y0, z0, t0, 0, 1, 0, electrons, true);
533}

◆ DisableAttachment()

void Garfield::AvalancheMC::DisableAttachment ( )
inline

Switch off attachment and multiplication.

Definition at line 69 of file AvalancheMC.hh.

69{ m_useAttachment = false; }

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMC::DisableAvalancheSizeLimit ( )
inline

Do not limit the maximum avalanche size.

Definition at line 83 of file AvalancheMC.hh.

83{ m_sizeCut = 0; }

◆ DisableDebugging()

void Garfield::AvalancheMC::DisableDebugging ( )
inline

Switch off debugging messages.

Definition at line 174 of file AvalancheMC.hh.

174{ m_debug = false; }

◆ DisableDiffusion()

void Garfield::AvalancheMC::DisableDiffusion ( )
inline

Switch off diffusion.

Definition at line 63 of file AvalancheMC.hh.

63{ m_useDiffusion = false; }

◆ DisablePlotting()

void Garfield::AvalancheMC::DisablePlotting ( )
inline

Switch off drift line plotting.

Definition at line 30 of file AvalancheMC.hh.

30{ m_viewer = nullptr; }

◆ DriftElectron()

bool Garfield::AvalancheMC::DriftElectron ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Simulate the drift line of an electron from a given starting point.

Definition at line 189 of file AvalancheMC.cc.

190 {
191 if (!m_sensor) {
192 std::cerr << m_className << "::DriftElectron: Sensor is not defined.\n";
193 return false;
194 }
195
196 m_endpointsElectrons.clear();
197 m_endpointsHoles.clear();
198 m_endpointsIons.clear();
199
200 m_nElectrons = 1;
201 m_nHoles = 0;
202 m_nIons = 0;
203
204 return DriftLine(x0, y0, z0, t0, Particle::Electron);
205}

Referenced by GarfieldPhysics::DoIt().

◆ DriftHole()

bool Garfield::AvalancheMC::DriftHole ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Simulate the drift line of a hole from a given starting point.

Definition at line 207 of file AvalancheMC.cc.

208 {
209 if (!m_sensor) {
210 std::cerr << m_className << "::DriftHole: Sensor is not defined.\n";
211 return false;
212 }
213
214 m_endpointsElectrons.clear();
215 m_endpointsHoles.clear();
216 m_endpointsIons.clear();
217
218 m_nElectrons = 0;
219 m_nHoles = 1;
220 m_nIons = 0;
221
222 return DriftLine(x0, y0, z0, t0, Particle::Hole);
223}

◆ DriftIon()

bool Garfield::AvalancheMC::DriftIon ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Simulate the drift line of an ion from a given starting point.

Definition at line 225 of file AvalancheMC.cc.

226 {
227 if (!m_sensor) {
228 std::cerr << m_className << "::DriftIon: Sensor is not defined.\n";
229 return false;
230 }
231
232 m_endpointsElectrons.clear();
233 m_endpointsHoles.clear();
234 m_endpointsIons.clear();
235
236 m_nElectrons = 0;
237 m_nHoles = 0;
238 m_nIons = 1;
239
240 return DriftLine(x0, y0, z0, t0, Particle::Ion);
241}

◆ EnableAttachment()

void Garfield::AvalancheMC::EnableAttachment ( )
inline

Switch on attachment (and multiplication) for drift line calculation (default: enabled). For avalanches the flag is ignored.

Definition at line 67 of file AvalancheMC.hh.

67{ m_useAttachment = true; }

◆ EnableAvalancheSizeLimit()

void Garfield::AvalancheMC::EnableAvalancheSizeLimit ( const unsigned int  size)
inline

Set a maximum avalanche size (ignore further multiplication once this size has been reached).

Definition at line 81 of file AvalancheMC.hh.

81{ m_sizeCut = size; }

◆ EnableDebugging()

void Garfield::AvalancheMC::EnableDebugging ( )
inline

Switch on debugging messages (default: off).

Definition at line 172 of file AvalancheMC.hh.

172{ m_debug = true; }

◆ EnableDiffusion()

void Garfield::AvalancheMC::EnableDiffusion ( )
inline

Switch on diffusion (default: enabled).

Definition at line 61 of file AvalancheMC.hh.

61{ m_useDiffusion = true; }

◆ EnableInducedChargeCalculation()

void Garfield::AvalancheMC::EnableInducedChargeCalculation ( const bool  on = true)
inline

Switch on calculation of induced charge (default: disabled).

Definition at line 46 of file AvalancheMC.hh.

46 {
47 m_doInducedCharge = on;
48 }

◆ EnableMagneticField()

void Garfield::AvalancheMC::EnableMagneticField ( const bool  on = true)
inline

Enable use of magnetic field in stepping algorithm.

Definition at line 77 of file AvalancheMC.hh.

77{ m_useBfield = on; }

◆ EnablePlotting()

void Garfield::AvalancheMC::EnablePlotting ( ViewDrift view)

Switch on drift line plotting.

Definition at line 40 of file AvalancheMC.cc.

40 {
41 if (!view) {
42 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
43 return;
44 }
45
46 m_viewer = view;
47}

◆ EnableProjectedPathIntegration()

void Garfield::AvalancheMC::EnableProjectedPathIntegration ( const bool  on = true)
inline

Switch on equilibration of multiplication and attachment over the drift line (default: enabled).

Definition at line 56 of file AvalancheMC.hh.

56 {
57 m_doEquilibration = on;
58 }

◆ EnableRKFSteps()

void Garfield::AvalancheMC::EnableRKFSteps ( const bool  on = true)
inline

Switch on Runge-Kutta-Fehlberg stepping (as opposed to simple straight-line steps.

Definition at line 52 of file AvalancheMC.hh.

52{ m_doRKF = on; }

◆ EnableSignalCalculation()

void Garfield::AvalancheMC::EnableSignalCalculation ( const bool  on = true)
inline

Switch on calculation of induced currents (default: disabled).

Definition at line 33 of file AvalancheMC.hh.

33{ m_doSignal = on; }

◆ EnableTcadTraps()

void Garfield::AvalancheMC::EnableTcadTraps ( const bool  on = true)
inline

Switch on calculating trapping with TCAD traps.

Definition at line 72 of file AvalancheMC.hh.

72{ m_useTcadTrapping = on; }

◆ EnableTcadVelocity()

void Garfield::AvalancheMC::EnableTcadVelocity ( const bool  on = true)
inline

Switch on TCAD velocity maps.

Definition at line 74 of file AvalancheMC.hh.

74{ m_useTcadVelocity = on; }

◆ GetAvalancheSize()

void Garfield::AvalancheMC::GetAvalancheSize ( unsigned int &  ne,
unsigned int &  ni 
) const
inline

Return the number of electrons and ions/holes in the avalanche.

Definition at line 110 of file AvalancheMC.hh.

110 {
111 ne = m_nElectrons;
112 ni = m_nIons;
113 }

◆ GetAvalancheSizeLimit()

int Garfield::AvalancheMC::GetAvalancheSizeLimit ( ) const
inline

Return the currently set avalanche size limit.

Definition at line 85 of file AvalancheMC.hh.

85{ return m_sizeCut; }

◆ GetDriftLinePoint()

void Garfield::AvalancheMC::GetDriftLinePoint ( const unsigned int  i,
double &  x,
double &  y,
double &  z,
double &  t 
) const

Return the coordinates and time of a point along the last drift line.

Definition at line 117 of file AvalancheMC.cc.

118 {
119 if (i >= m_drift.size()) {
120 std::cerr << m_className << "::GetDriftLinePoint: Index out of range.\n";
121 return;
122 }
123
124 x = m_drift[i].x[0];
125 y = m_drift[i].x[1];
126 z = m_drift[i].x[2];
127 t = m_drift[i].t;
128}

◆ GetElectronEndpoint()

void Garfield::AvalancheMC::GetElectronEndpoint ( const unsigned int  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.

Parameters
iindex of the drift line
x0,y0,z0,t0coordinates and time of the starting point
x1,y1,z1,t1coordinates and time of the end point
statusstatus code (see GarfieldConstants.hh)

Definition at line 169 of file AvalancheMC.cc.

172 {
173 if (i >= m_endpointsElectrons.size()) {
174 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
175 return;
176 }
177
178 x0 = m_endpointsElectrons[i].x0[0];
179 y0 = m_endpointsElectrons[i].x0[1];
180 z0 = m_endpointsElectrons[i].x0[2];
181 t0 = m_endpointsElectrons[i].t0;
182 x1 = m_endpointsElectrons[i].x1[0];
183 y1 = m_endpointsElectrons[i].x1[1];
184 z1 = m_endpointsElectrons[i].x1[2];
185 t1 = m_endpointsElectrons[i].t1;
186 status = m_endpointsElectrons[i].status;
187}

Referenced by GarfieldPhysics::DoIt().

◆ GetHoleEndpoint()

void Garfield::AvalancheMC::GetHoleEndpoint ( const unsigned int  i,
double &  x0,
double &  y0,
double &  z0,
double &  t0,
double &  x1,
double &  y1,
double &  z1,
double &  t1,
int &  status 
) const

Definition at line 130 of file AvalancheMC.cc.

133 {
134 if (i >= m_endpointsHoles.size()) {
135 std::cerr << m_className << "::GetHoleEndpoint: Index out of range.\n";
136 return;
137 }
138
139 x0 = m_endpointsHoles[i].x0[0];
140 y0 = m_endpointsHoles[i].x0[1];
141 z0 = m_endpointsHoles[i].x0[2];
142 t0 = m_endpointsHoles[i].t0;
143 x1 = m_endpointsHoles[i].x1[0];
144 y1 = m_endpointsHoles[i].x1[1];
145 z1 = m_endpointsHoles[i].x1[2];
146 t1 = m_endpointsHoles[i].t1;
147 status = m_endpointsHoles[i].status;
148}

◆ GetIonEndpoint()

void Garfield::AvalancheMC::GetIonEndpoint ( const unsigned int  i,
double &  x0,
double &  y0,
double &  z0,
double &  t0,
double &  x1,
double &  y1,
double &  z1,
double &  t1,
int &  status 
) const

Definition at line 150 of file AvalancheMC.cc.

152 {
153 if (i >= m_endpointsIons.size()) {
154 std::cerr << m_className << "::GetIonEndpoint: Index out of range.\n";
155 return;
156 }
157
158 x0 = m_endpointsIons[i].x0[0];
159 y0 = m_endpointsIons[i].x0[1];
160 z0 = m_endpointsIons[i].x0[2];
161 t0 = m_endpointsIons[i].t0;
162 x1 = m_endpointsIons[i].x1[0];
163 y1 = m_endpointsIons[i].x1[1];
164 z1 = m_endpointsIons[i].x1[2];
165 t1 = m_endpointsIons[i].t1;
166 status = m_endpointsIons[i].status;
167}

◆ GetNumberOfDriftLinePoints()

unsigned int Garfield::AvalancheMC::GetNumberOfDriftLinePoints ( ) const
inline

Return the number of points along the last simulated drift line.

Definition at line 116 of file AvalancheMC.hh.

116{ return m_drift.size(); }

◆ GetNumberOfElectronEndpoints()

unsigned int Garfield::AvalancheMC::GetNumberOfElectronEndpoints ( ) const
inline

Return the number of electron trajectories in the last simulated avalanche (including captured electrons).

Definition at line 123 of file AvalancheMC.hh.

123 {
124 return m_endpointsElectrons.size();
125 }

◆ GetNumberOfHoleEndpoints()

unsigned int Garfield::AvalancheMC::GetNumberOfHoleEndpoints ( ) const
inline

Return the number of hole trajectories in the last simulated avalanche (including captured holes).

Definition at line 128 of file AvalancheMC.hh.

128 {
129 return m_endpointsHoles.size();
130 }

◆ GetNumberOfIonEndpoints()

unsigned int Garfield::AvalancheMC::GetNumberOfIonEndpoints ( ) const
inline

Return the number of ion trajectories.

Definition at line 132 of file AvalancheMC.hh.

132 {
133 return m_endpointsIons.size();
134 }

◆ SetCollisionSteps()

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 79 of file AvalancheMC.cc.

79 {
80 m_stepModel = StepModel::CollisionTime;
81 if (n < 1) {
82 std::cerr << m_className << "::SetCollisionSteps:\n "
83 << "Number of collisions set to default value (100).\n";
84 m_nMc = 100;
85 return;
86 }
87 if (m_debug) {
88 std::cout << m_className << "::SetCollisionSteps:\n "
89 << "Number of collisions to be skipped set to " << n << ".\n";
90 }
91 m_nMc = n;
92}

◆ SetDistanceSteps()

void Garfield::AvalancheMC::SetDistanceSteps ( const double  d = 0.001)

Use fixed distance steps (default 10 um).

Definition at line 64 of file AvalancheMC.cc.

64 {
65 m_stepModel = StepModel::FixedDistance;
66 if (d < Small) {
67 std::cerr << m_className << "::SetDistanceSteps:\n "
68 << "Step size is too small. Using default (10 um) instead.\n";
69 m_dMc = 0.001;
70 return;
71 }
72 if (m_debug) {
73 std::cout << m_className << "::SetDistanceSteps:\n"
74 << " Step size set to " << d << " cm.\n";
75 }
76 m_dMc = d;
77}

◆ SetElectronSignalScalingFactor()

void Garfield::AvalancheMC::SetElectronSignalScalingFactor ( const double  scale)
inline

Set multiplication factor for the signal induced by electrons.

Definition at line 103 of file AvalancheMC.hh.

103{ m_scaleE = scale; }

◆ SetHoleSignalScalingFactor()

void Garfield::AvalancheMC::SetHoleSignalScalingFactor ( const double  scale)
inline

Set multiplication factor for the signal induced by holes.

Definition at line 105 of file AvalancheMC.hh.

105{ m_scaleH = scale; }

◆ SetIonSignalScalingFactor()

void Garfield::AvalancheMC::SetIonSignalScalingFactor ( const double  scale)
inline

Set multiplication factor for the signal induced by ions.

Definition at line 107 of file AvalancheMC.hh.

107{ m_scaleI = scale; }

◆ SetSensor()

void Garfield::AvalancheMC::SetSensor ( Sensor s)

Set the sensor.

Definition at line 31 of file AvalancheMC.cc.

31 {
32 if (!sensor) {
33 std::cerr << m_className << "::SetSensor: Null pointer.\n";
34 return;
35 }
36
37 m_sensor = sensor;
38}

Referenced by GarfieldPhysics::InitializePhysics().

◆ SetSignalAveragingOrder()

void Garfield::AvalancheMC::SetSignalAveragingOrder ( const unsigned int  navg)
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 $2\times navg + 1$ point Newton-Raphson integration. Default: 1.

Definition at line 38 of file AvalancheMC.hh.

38{ m_navg = navg; }

◆ SetStepDistanceFunction()

void Garfield::AvalancheMC::SetStepDistanceFunction ( double(*)(double x, double y, double z)  f)

Retrieve the step distance from a user-supplied function.

Definition at line 94 of file AvalancheMC.cc.

95 {
96
97 if (!f) {
98 std::cerr << m_className << "::SetStepDistanceFunction: Null pointer.\n";
99 return;
100 }
101 m_fStep = f;
102 m_stepModel = StepModel::UserDistance;
103}

◆ SetTimeSteps()

void Garfield::AvalancheMC::SetTimeSteps ( const double  d = 0.02)

Use fixed-time steps (default 20 ps).

Definition at line 49 of file AvalancheMC.cc.

49 {
50 m_stepModel = StepModel::FixedTime;
51 if (d < Small) {
52 std::cerr << m_className << "::SetTimeSteps:\n "
53 << "Step size is too small. Using default (20 ps) instead.\n";
54 m_tMc = 0.02;
55 return;
56 }
57 if (m_debug) {
58 std::cout << m_className << "::SetTimeSteps:\n"
59 << " Step size set to " << d << " ns.\n";
60 }
61 m_tMc = d;
62}

◆ SetTimeWindow()

void Garfield::AvalancheMC::SetTimeWindow ( const double  t0,
const double  t1 
)

Define a time interval (only carriers inside the interval are drifted).

Definition at line 105 of file AvalancheMC.cc.

105 {
106 if (fabs(t1 - t0) < Small) {
107 std::cerr << m_className << "::SetTimeWindow:\n"
108 << " Time interval must be greater than zero.\n";
109 return;
110 }
111
112 m_tMin = std::min(t0, t1);
113 m_tMax = std::max(t0, t1);
114 m_hasTimeWindow = true;
115}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

◆ UnsetTimeWindow()

void Garfield::AvalancheMC::UnsetTimeWindow ( )
inline

Do not limit the time interval within which carriers are drifted.

Definition at line 100 of file AvalancheMC.hh.

100{ m_hasTimeWindow = false; }

◆ UseWeightingPotential()

void Garfield::AvalancheMC::UseWeightingPotential ( const bool  on = true)
inline

Use the weighting potential (as opposed to the weighting field) for calculating the induced signal.

Definition at line 41 of file AvalancheMC.hh.

41 {
42 m_useWeightingPotential = on;
43 }

The documentation for this class was generated from the following files: