Garfield++ v2r0
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 ()
 
void EnableSignalCalculation ()
 Switch on calculation of induced currents (default: disabled).
 
void DisableSignalCalculation ()
 
void EnableInducedChargeCalculation ()
 Switch on calculation of induced charge (default: disabled).
 
void DisableInducedChargeCalculation ()
 
void EnableProjectedPathIntegration ()
 
void DisableProjectedPathIntegration ()
 
void EnableDiffusion ()
 Switch on diffusion (default: enabled)
 
void DisableDiffusion ()
 
void EnableAttachment ()
 
void DisableAttachment ()
 
void EnableTcadTraps ()
 Switch on calculating trapping with TCAD traps.
 
void DisableTcadTraps ()
 
void EnableTcadVelocity ()
 Switch on TCAD velocity maps.
 
void DisableTcadVelocity ()
 
void EnableMagneticField ()
 Enable use of magnetic field in stepping algorithm.
 
void DisableMagneticField ()
 
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 int n=100)
 
void SetTimeWindow (const double t0, const double t1)
 Define a time interval (only carriers inside the interval are drifted).
 
void UnsetTimeWindow ()
 
void SetHoles ()
 Treat positive charge carriers as holes (default: ions).
 
void SetIons ()
 
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 (int &ne, 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
 
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 with a given starting point.
 
bool DriftHole (const double x0, const double y0, const double z0, const double t0)
 
bool DriftIon (const double x0, const double y0, const double z0, const double t0)
 
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 with given starting point.
 
bool AvalancheHole (const double x0, const double y0, const double z0, const double t0, const bool electron=false)
 
bool AvalancheElectronHole (const double x0, const double y0, const double z0, const double t0)
 
void EnableDebugging ()
 Switch on debugging messages (default: off).
 
void DisableDebugging ()
 

Detailed Description

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

Definition at line 15 of file AvalancheMC.hh.

Constructor & Destructor Documentation

◆ AvalancheMC()

Garfield::AvalancheMC::AvalancheMC ( )

Constructor.

Definition at line 15 of file AvalancheMC.cc.

16 : m_sensor(NULL),
17 m_stepModel(2),
18 m_tMc(0.02),
19 m_dMc(0.001),
20 m_nMc(100),
21 m_hasTimeWindow(false),
22 m_tMin(0.),
23 m_tMax(0.),
24 m_nElectrons(0),
25 m_nHoles(0),
26 m_nIons(0),
27 m_viewer(NULL),
28 m_useSignal(false),
29 m_useInducedCharge(false),
30 m_useEquilibration(true),
31 m_useDiffusion(true),
32 m_useAttachment(true),
33 m_useBfield(false),
34 m_useIons(true),
35 m_withElectrons(true),
36 m_withHoles(true),
37 m_scaleElectronSignal(1.),
38 m_scaleHoleSignal(1.),
39 m_scaleIonSignal(1.),
40 m_useTcadTrapping(false),
41 m_useTcadVelocity(false),
42 m_debug(false) {
43
44 m_className = "AvalancheMC";
45
46 m_drift.reserve(10000);
47}

◆ ~AvalancheMC()

Garfield::AvalancheMC::~AvalancheMC ( )
inline

Destructor.

Definition at line 21 of file AvalancheMC.hh.

21{}

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 with given starting point.

Definition at line 511 of file AvalancheMC.cc.

513 {
514
515 m_withHoles = holes;
516 return Avalanche(x0, y0, z0, t0, 1, 0, 0);
517}

◆ AvalancheElectronHole()

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

Definition at line 527 of file AvalancheMC.cc.

528 {
529
530 m_withElectrons = m_withHoles = true;
531 return Avalanche(x0, y0, z0, t0, 1, 1, 0);
532}

◆ AvalancheHole()

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

Definition at line 519 of file AvalancheMC.cc.

521 {
522
523 m_withElectrons = electrons;
524 return Avalanche(x0, y0, z0, t0, 0, 1, 0);
525}

◆ DisableAttachment()

void Garfield::AvalancheMC::DisableAttachment ( )
inline

Definition at line 50 of file AvalancheMC.hh.

50{ m_useAttachment = false; }

◆ DisableDebugging()

void Garfield::AvalancheMC::DisableDebugging ( )
inline

Definition at line 151 of file AvalancheMC.hh.

151{ m_debug = false; }

◆ DisableDiffusion()

void Garfield::AvalancheMC::DisableDiffusion ( )
inline

Definition at line 45 of file AvalancheMC.hh.

45{ m_useDiffusion = false; }

◆ DisableInducedChargeCalculation()

void Garfield::AvalancheMC::DisableInducedChargeCalculation ( )
inline

Definition at line 36 of file AvalancheMC.hh.

36{ m_useInducedCharge = false; }

◆ DisableMagneticField()

void Garfield::AvalancheMC::DisableMagneticField ( )
inline

Definition at line 62 of file AvalancheMC.hh.

62{ m_useBfield = false; }

◆ DisablePlotting()

void Garfield::AvalancheMC::DisablePlotting ( )
inline

Definition at line 28 of file AvalancheMC.hh.

28{ m_viewer = NULL; }

◆ DisableProjectedPathIntegration()

void Garfield::AvalancheMC::DisableProjectedPathIntegration ( )
inline

Definition at line 41 of file AvalancheMC.hh.

41{ m_useEquilibration = false; }

◆ DisableSignalCalculation()

void Garfield::AvalancheMC::DisableSignalCalculation ( )
inline

Definition at line 32 of file AvalancheMC.hh.

32{ m_useSignal = false; }

◆ DisableTcadTraps()

void Garfield::AvalancheMC::DisableTcadTraps ( )
inline

Definition at line 54 of file AvalancheMC.hh.

54{ m_useTcadTrapping = false; }

◆ DisableTcadVelocity()

void Garfield::AvalancheMC::DisableTcadVelocity ( )
inline

Definition at line 58 of file AvalancheMC.hh.

58{ m_useTcadVelocity = false; }

◆ DriftElectron()

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

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

Definition at line 210 of file AvalancheMC.cc.

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

Referenced by GarfieldPhysics::DoIt().

◆ DriftHole()

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

Definition at line 230 of file AvalancheMC.cc.

231 {
232
233 if (!m_sensor) {
234 std::cerr << m_className << "::DriftHole:\n Sensor is not defined.\n";
235 return false;
236 }
237
238 m_endpointsElectrons.clear();
239 m_endpointsHoles.clear();
240 m_endpointsIons.clear();
241
242 m_nElectrons = 0;
243 m_nHoles = 1;
244 m_nIons = 0;
245
246 return DriftLine(x0, y0, z0, t0, 1);
247}

◆ DriftIon()

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

Definition at line 249 of file AvalancheMC.cc.

250 {
251
252 if (!m_sensor) {
253 std::cerr << m_className << "::DriftIon:\n Sensor is not defined.\n";
254 return false;
255 }
256
257 m_endpointsElectrons.clear();
258 m_endpointsHoles.clear();
259 m_endpointsIons.clear();
260
261 m_nElectrons = 0;
262 m_nHoles = 0;
263 m_nIons = 1;
264
265 return DriftLine(x0, y0, z0, t0, 2);
266}

◆ 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 49 of file AvalancheMC.hh.

49{ m_useAttachment = true; }

◆ EnableDebugging()

void Garfield::AvalancheMC::EnableDebugging ( )
inline

Switch on debugging messages (default: off).

Definition at line 150 of file AvalancheMC.hh.

150{ m_debug = true; }

◆ EnableDiffusion()

void Garfield::AvalancheMC::EnableDiffusion ( )
inline

Switch on diffusion (default: enabled)

Definition at line 44 of file AvalancheMC.hh.

44{ m_useDiffusion = true; }

◆ EnableInducedChargeCalculation()

void Garfield::AvalancheMC::EnableInducedChargeCalculation ( )
inline

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

Definition at line 35 of file AvalancheMC.hh.

35{ m_useInducedCharge = true; }

◆ EnableMagneticField()

void Garfield::AvalancheMC::EnableMagneticField ( )
inline

Enable use of magnetic field in stepping algorithm.

Definition at line 61 of file AvalancheMC.hh.

61{ m_useBfield = true; }

◆ EnablePlotting()

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

Switch on drift line plotting.

Definition at line 59 of file AvalancheMC.cc.

59 {
60
61 if (!view) {
62 std::cerr << m_className << "::EnablePlotting:\n Null pointer.\n";
63 return;
64 }
65
66 m_viewer = view;
67}

◆ EnableProjectedPathIntegration()

void Garfield::AvalancheMC::EnableProjectedPathIntegration ( )
inline

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

Definition at line 40 of file AvalancheMC.hh.

40{ m_useEquilibration = true; }

◆ EnableSignalCalculation()

void Garfield::AvalancheMC::EnableSignalCalculation ( )
inline

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

Definition at line 31 of file AvalancheMC.hh.

31{ m_useSignal = true; }

◆ EnableTcadTraps()

void Garfield::AvalancheMC::EnableTcadTraps ( )
inline

Switch on calculating trapping with TCAD traps.

Definition at line 53 of file AvalancheMC.hh.

53{ m_useTcadTrapping = true; }

◆ EnableTcadVelocity()

void Garfield::AvalancheMC::EnableTcadVelocity ( )
inline

Switch on TCAD velocity maps.

Definition at line 57 of file AvalancheMC.hh.

57{ m_useTcadVelocity = true; }

◆ GetAvalancheSize()

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

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

Definition at line 94 of file AvalancheMC.hh.

94 {
95 ne = m_nElectrons;
96 ni = m_nIons;
97 }

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

131 {
132
133 if (i >= m_drift.size()) {
134 std::cerr << m_className << "::GetDriftLinePoint:\n"
135 << " Drift line point " << i << " does not exist.\n";
136 return;
137 }
138
139 x = m_drift[i].x;
140 y = m_drift[i].y;
141 z = m_drift[i].z;
142 t = m_drift[i].t;
143}

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

191 {
192
193 if (i >= m_endpointsElectrons.size()) {
194 std::cerr << m_className << "::GetElectronEndpoint:\n"
195 << " Endpoint " << i << " does not exist.\n";
196 return;
197 }
198
199 x0 = m_endpointsElectrons[i].x0;
200 y0 = m_endpointsElectrons[i].y0;
201 z0 = m_endpointsElectrons[i].z0;
202 t0 = m_endpointsElectrons[i].t0;
203 x1 = m_endpointsElectrons[i].x1;
204 y1 = m_endpointsElectrons[i].y1;
205 z1 = m_endpointsElectrons[i].z1;
206 t1 = m_endpointsElectrons[i].t1;
207 status = m_endpointsElectrons[i].status;
208}

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

148 {
149
150 if (i >= m_endpointsHoles.size()) {
151 std::cerr << m_className << "::GetHoleEndpoint:\n"
152 << " Endpoint " << i << " does not exist.\n";
153 return;
154 }
155
156 x0 = m_endpointsHoles[i].x0;
157 y0 = m_endpointsHoles[i].y0;
158 z0 = m_endpointsHoles[i].z0;
159 t0 = m_endpointsHoles[i].t0;
160 x1 = m_endpointsHoles[i].x1;
161 y1 = m_endpointsHoles[i].y1;
162 z1 = m_endpointsHoles[i].z1;
163 t1 = m_endpointsHoles[i].t1;
164 status = m_endpointsHoles[i].status;
165}

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

169 {
170
171 if (i >= m_endpointsIons.size()) {
172 std::cerr << m_className << "::GetIonEndpoint:\n"
173 << " Endpoint " << i << " does not exist.\n";
174 return;
175 }
176
177 x0 = m_endpointsIons[i].x0;
178 y0 = m_endpointsIons[i].y0;
179 z0 = m_endpointsIons[i].z0;
180 t0 = m_endpointsIons[i].t0;
181 x1 = m_endpointsIons[i].x1;
182 y1 = m_endpointsIons[i].y1;
183 z1 = m_endpointsIons[i].z1;
184 t1 = m_endpointsIons[i].t1;
185 status = m_endpointsIons[i].status;
186}

◆ GetNumberOfDriftLinePoints()

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

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

Definition at line 100 of file AvalancheMC.hh.

100{ 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 107 of file AvalancheMC.hh.

107 {
108 return m_endpointsElectrons.size();
109 }

◆ GetNumberOfHoleEndpoints()

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

Definition at line 110 of file AvalancheMC.hh.

110 {
111 return m_endpointsHoles.size();
112 }

◆ GetNumberOfIonEndpoints()

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

Definition at line 113 of file AvalancheMC.hh.

113 {
114 return m_endpointsIons.size();
115 }

◆ SetCollisionSteps()

void Garfield::AvalancheMC::SetCollisionSteps ( const int  n = 100)

Use exponentially distributed time steps with mean equal to the specified multiple of the collision time (default model).

Definition at line 101 of file AvalancheMC.cc.

101 {
102
103 m_stepModel = 2;
104 if (n < 1) {
105 std::cerr << m_className << "::SetCollisionSteps:\n "
106 << "Number of collisions set to default value (100).\n";
107 m_nMc = 100;
108 return;
109 }
110 if (m_debug) {
111 std::cout << m_className << "::SetCollisionSteps:\n "
112 << "Number of collisions to be skipped set to " << n << ".\n";
113 }
114 m_nMc = n;
115}

◆ SetDistanceSteps()

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

Use fixed distance steps (default 10 um).

Definition at line 85 of file AvalancheMC.cc.

85 {
86
87 m_stepModel = 1;
88 if (d < Small) {
89 std::cerr << m_className << "::SetDistanceSteps:\n "
90 << "Step size is too small. Using default (10 um) instead.\n";
91 m_dMc = 0.001;
92 return;
93 }
94 if (m_debug) {
95 std::cout << m_className << "::SetDistanceSteps:\n"
96 << " Step size set to " << d << " cm.\n";
97 }
98 m_dMc = d;
99}

◆ SetElectronSignalScalingFactor()

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

Set multiplication factor for the signal induced by electrons.

Definition at line 81 of file AvalancheMC.hh.

81 {
82 m_scaleElectronSignal = scale;
83 }

◆ SetHoles()

void Garfield::AvalancheMC::SetHoles ( )
inline

Treat positive charge carriers as holes (default: ions).

Definition at line 77 of file AvalancheMC.hh.

77{ m_useIons = false; }

◆ SetHoleSignalScalingFactor()

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

Set multiplication factor for the signal induced by holes.

Definition at line 85 of file AvalancheMC.hh.

85 {
86 m_scaleHoleSignal = scale;
87 }

◆ SetIons()

void Garfield::AvalancheMC::SetIons ( )
inline

Definition at line 78 of file AvalancheMC.hh.

78{ m_useIons = true; }

◆ SetIonSignalScalingFactor()

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

Set multiplication factor for the signal induced by ions.

Definition at line 89 of file AvalancheMC.hh.

89 {
90 m_scaleIonSignal = scale;
91 }

◆ SetSensor()

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

Set the sensor.

Definition at line 49 of file AvalancheMC.cc.

49 {
50
51 if (!sensor) {
52 std::cerr << m_className << "::SetSensor:\n Null pointer.\n";
53 return;
54 }
55
56 m_sensor = sensor;
57}

Referenced by GarfieldPhysics::InitializePhysics().

◆ SetTimeSteps()

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

Use fixed-time steps (default 20 ps).

Definition at line 69 of file AvalancheMC.cc.

69 {
70
71 m_stepModel = 0;
72 if (d < Small) {
73 std::cerr << m_className << "::SetTimeSteps:\n "
74 << "Step size is too small. Using default (20 ps) instead.\n";
75 m_tMc = 0.02;
76 return;
77 }
78 if (m_debug) {
79 std::cout << m_className << "::SetTimeSteps:\n"
80 << " Step size set to " << d << " ns.\n";
81 }
82 m_tMc = d;
83}

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

117 {
118
119 if (fabs(t1 - t0) < Small) {
120 std::cerr << m_className << "::SetTimeWindow:\n"
121 << " Time interval must be greater than zero.\n";
122 return;
123 }
124
125 m_tMin = std::min(t0, t1);
126 m_tMax = std::max(t0, t1);
127 m_hasTimeWindow = true;
128}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

◆ UnsetTimeWindow()

void Garfield::AvalancheMC::UnsetTimeWindow ( )
inline

Definition at line 74 of file AvalancheMC.hh.

74{ m_hasTimeWindow = false; }

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