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

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 ()
 

Detailed Description

Calculate electron drift lines and avalanches using microscopic tracking.

Definition at line 17 of file AvalancheMicroscopic.hh.

Constructor & Destructor Documentation

◆ AvalancheMicroscopic()

Garfield::AvalancheMicroscopic::AvalancheMicroscopic ( )

Constructor.

Definition at line 25 of file AvalancheMicroscopic.cc.

26 : m_sensor(NULL),
27 m_nElectrons(0),
28 m_nHoles(0),
29 m_nIons(0),
30 m_usePlotting(false),
31 m_viewer(NULL),
32 m_plotExcitations(true),
33 m_plotIonisations(true),
34 m_plotAttachments(true),
35 m_histElectronEnergy(NULL),
36 m_histHoleEnergy(NULL),
37 m_histDistance(NULL),
38 m_distanceOption('r'),
39 m_histSecondary(NULL),
40 m_useSignal(false),
41 m_useInducedCharge(false),
42 m_useDriftLines(false),
43 m_usePhotons(false),
44 m_useBandStructureDefault(true),
45 m_useNullCollisionSteps(false),
46 m_useBfield(false),
47 m_rb11(1.),
48 m_rb12(0.),
49 m_rb13(0.),
50 m_rb21(0.),
51 m_rb22(1.),
52 m_rb23(0.),
53 m_rb31(0.),
54 m_rb32(0.),
55 m_rb33(1.),
56 m_rx22(1.),
57 m_rx23(0.),
58 m_rx32(0.),
59 m_rx33(1.),
60 m_deltaCut(0.),
61 m_gammaCut(0.),
62 m_sizeCut(0),
63 m_nCollSkip(100),
64 m_hasTimeWindow(false),
65 m_tMin(0.),
66 m_tMax(0.),
67 m_hasUserHandleStep(false),
68 m_hasUserHandleAttachment(false),
69 m_hasUserHandleInelastic(false),
70 m_hasUserHandleIonisation(false),
71 m_userHandleStep(0),
72 m_userHandleAttachment(0),
73 m_userHandleInelastic(0),
74 m_userHandleIonisation(0),
75 m_debug(false) {
76
77 m_className = "AvalancheMicroscopic";
78
79 m_endpointsElectrons.reserve(10000);
80 m_endpointsHoles.reserve(10000);
81 m_photons.reserve(1000);
82
83}

◆ ~AvalancheMicroscopic()

Garfield::AvalancheMicroscopic::~AvalancheMicroscopic ( )
inline

Destructor.

Definition at line 23 of file AvalancheMicroscopic.hh.

23{}

Member Function Documentation

◆ AvalancheElectron()

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.

521 {
522
523 // Clear the list of electrons, holes and photons.
524 m_endpointsElectrons.clear();
525 m_endpointsHoles.clear();
526 m_photons.clear();
527
528 // Reset the particle counters.
529 m_nElectrons = m_nHoles = m_nIons = 0;
530
531 return TransportElectron(x0, y0, z0, t0, e0, dx0, dy0, dz0, true, false);
532}

Referenced by GarfieldPhysics::DoIt().

◆ DisableAttachmentMarkers()

void Garfield::AvalancheMicroscopic::DisableAttachmentMarkers ( )
inline

Definition at line 39 of file AvalancheMicroscopic.hh.

39{ m_plotAttachments = false; }

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMicroscopic::DisableAvalancheSizeLimit ( )
inline

Definition at line 99 of file AvalancheMicroscopic.hh.

99{ m_sizeCut = 0; }

◆ DisableBandStructure()

void Garfield::AvalancheMicroscopic::DisableBandStructure ( )
inline

Definition at line 81 of file AvalancheMicroscopic.hh.

81{ m_useBandStructureDefault = false; }

◆ DisableDebugging()

void Garfield::AvalancheMicroscopic::DisableDebugging ( )
inline

Definition at line 202 of file AvalancheMicroscopic.hh.

202{ m_debug = false; }

◆ DisableDistanceHistogramming() [1/2]

void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming ( )

Definition at line 206 of file AvalancheMicroscopic.cc.

206 {
207
208 m_histDistance = NULL;
209 m_distanceHistogramType.clear();
210}

◆ DisableDistanceHistogramming() [2/2]

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.

186 {
187
188 if (m_distanceHistogramType.empty()) {
189 std::cerr << m_className << "::DisableDistanceHistogramming:\n";
190 std::cerr << " Collision type " << type << " is not histogrammed.\n";
191 return;
192 }
193 const unsigned int nDistanceHistogramTypes = m_distanceHistogramType.size();
194 for (unsigned int i = 0; i < nDistanceHistogramTypes; ++i) {
195 if (m_distanceHistogramType[i] != type) continue;
196 m_distanceHistogramType.erase(m_distanceHistogramType.begin() + i);
197 std::cout << " Histogramming of collision type " << type
198 << " disabled.\n";
199 return;
200 }
201
202 std::cerr << m_className << "::DisableDistanceHistogramming:\n";
203 std::cerr << " Collision type " << type << " is not histogrammed.\n";
204}

◆ DisableDriftLines()

void Garfield::AvalancheMicroscopic::DisableDriftLines ( )
inline

Definition at line 72 of file AvalancheMicroscopic.hh.

72{ m_useDriftLines = false; }

◆ DisableElectronEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableElectronEnergyHistogramming ( )
inline

Definition at line 51 of file AvalancheMicroscopic.hh.

51{ m_histElectronEnergy = NULL; }

◆ DisableExcitationMarkers()

void Garfield::AvalancheMicroscopic::DisableExcitationMarkers ( )
inline

Definition at line 33 of file AvalancheMicroscopic.hh.

33{ m_plotExcitations = false; }

◆ DisableHoleEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableHoleEnergyHistogramming ( )
inline

Definition at line 53 of file AvalancheMicroscopic.hh.

53{ m_histHoleEnergy = NULL; }

◆ DisableInducedChargeCalculation()

void Garfield::AvalancheMicroscopic::DisableInducedChargeCalculation ( )
inline

Definition at line 47 of file AvalancheMicroscopic.hh.

47{ m_useInducedCharge = false; }

◆ DisableIonisationMarkers()

void Garfield::AvalancheMicroscopic::DisableIonisationMarkers ( )
inline

Definition at line 36 of file AvalancheMicroscopic.hh.

36{ m_plotIonisations = false; }

◆ DisableMagneticField()

void Garfield::AvalancheMicroscopic::DisableMagneticField ( )
inline

Definition at line 104 of file AvalancheMicroscopic.hh.

104{ m_useBfield = false; }

◆ DisableNullCollisionSteps()

void Garfield::AvalancheMicroscopic::DisableNullCollisionSteps ( )
inline

Definition at line 85 of file AvalancheMicroscopic.hh.

85{ m_useNullCollisionSteps = false; }

◆ DisablePhotonTransport()

void Garfield::AvalancheMicroscopic::DisablePhotonTransport ( )
inline

Definition at line 77 of file AvalancheMicroscopic.hh.

77{ m_usePhotons = false; }

◆ DisablePlotting()

void Garfield::AvalancheMicroscopic::DisablePlotting ( )

Definition at line 110 of file AvalancheMicroscopic.cc.

110 {
111
112 m_viewer = NULL;
113 m_usePlotting = false;
114}

◆ DisableSecondaryEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableSecondaryEnergyHistogramming ( )
inline

Definition at line 68 of file AvalancheMicroscopic.hh.

68{ m_histSecondary = NULL; }

◆ DisableSignalCalculation()

void Garfield::AvalancheMicroscopic::DisableSignalCalculation ( )
inline

Definition at line 43 of file AvalancheMicroscopic.hh.

43{ m_useSignal = false; }

◆ DriftElectron()

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.

Parameters
x0,y0,z0,t0starting point of the electron
e0initial energy of the electron
dx0,dy0,dz0initial 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.

504 {
505
506 // Clear the list of electrons and photons.
507 m_endpointsElectrons.clear();
508 m_endpointsHoles.clear();
509 m_photons.clear();
510
511 // Reset the particle counters.
512 m_nElectrons = m_nHoles = m_nIons = 0;
513
514 return TransportElectron(x0, y0, z0, t0, e0, dx0, dy0, dz0, false, false);
515}

◆ EnableAttachmentMarkers()

void Garfield::AvalancheMicroscopic::EnableAttachmentMarkers ( )
inline

Draw a marker at every attachment.

Definition at line 38 of file AvalancheMicroscopic.hh.

38{ m_plotAttachments = true; }

◆ EnableAvalancheSizeLimit()

void Garfield::AvalancheMicroscopic::EnableAvalancheSizeLimit ( const unsigned int  size)
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.

98{ m_sizeCut = size; }

◆ EnableBandStructure()

void Garfield::AvalancheMicroscopic::EnableBandStructure ( )
inline

Switch on stepping according to band structure E(k), for semiconductors.

Definition at line 80 of file AvalancheMicroscopic.hh.

80{ m_useBandStructureDefault = true; }

◆ EnableDebugging()

void Garfield::AvalancheMicroscopic::EnableDebugging ( )
inline

Switch on debugging messages.

Definition at line 201 of file AvalancheMicroscopic.hh.

201{ m_debug = true; }

◆ EnableDistanceHistogramming()

void Garfield::AvalancheMicroscopic::EnableDistanceHistogramming ( const int  type)

Definition at line 163 of file AvalancheMicroscopic.cc.

163 {
164
165 // Check if this type of collision is already registered
166 // for histogramming.
167 const unsigned int nDistanceHistogramTypes = m_distanceHistogramType.size();
168 if (nDistanceHistogramTypes > 0) {
169 for (unsigned int i = 0; i < nDistanceHistogramTypes; ++i) {
170 if (m_distanceHistogramType[i] != type) continue;
171 std::cout << m_className << "::EnableDistanceHistogramming:\n";
172 std::cout << " Collision type " << type
173 << " is already being histogrammed.\n";
174 return;
175 }
176 }
177
178 m_distanceHistogramType.push_back(type);
179 std::cout << m_className << "::EnableDistanceHistogramming:\n";
180 std::cout << " Histogramming of collision type " << type << " enabled.\n";
181 if (!m_histDistance) {
182 std::cout << " Don't forget to set the histogram.\n";
183 }
184}

◆ EnableDriftLines()

void Garfield::AvalancheMicroscopic::EnableDriftLines ( )
inline

Switch on storage of drift lines.

Definition at line 71 of file AvalancheMicroscopic.hh.

71{ m_useDriftLines = true; }

Referenced by EnablePlotting().

◆ EnableElectronEnergyHistogramming()

void Garfield::AvalancheMicroscopic::EnableElectronEnergyHistogramming ( TH1 *  histo)

Switch on filling histograms for electron energy distribution.

Definition at line 116 of file AvalancheMicroscopic.cc.

116 {
117
118 if (!histo) {
119 std::cerr << m_className << "::EnableElectronEnergyHistogramming:\n"
120 << " Null pointer.\n";
121 return;
122 }
123
124 m_histElectronEnergy = histo;
125}

◆ EnableExcitationMarkers()

void Garfield::AvalancheMicroscopic::EnableExcitationMarkers ( )
inline

Draw a marker at every excitation.

Definition at line 32 of file AvalancheMicroscopic.hh.

32{ m_plotExcitations = true; }

◆ EnableHoleEnergyHistogramming()

void Garfield::AvalancheMicroscopic::EnableHoleEnergyHistogramming ( TH1 *  histo)

Definition at line 127 of file AvalancheMicroscopic.cc.

127 {
128
129 if (!histo) {
130 std::cerr << m_className << "::EnableHoleEnergyHistogramming:\n"
131 << " Null pointer.\n";
132 return;
133 }
134
135 m_histHoleEnergy = histo;
136}

◆ EnableInducedChargeCalculation()

void Garfield::AvalancheMicroscopic::EnableInducedChargeCalculation ( )
inline

Switch on calculation of the total induced charge.

Definition at line 46 of file AvalancheMicroscopic.hh.

46{ m_useInducedCharge = true; }

◆ EnableIonisationMarkers()

void Garfield::AvalancheMicroscopic::EnableIonisationMarkers ( )
inline

Draw a marker at every ionising collision.

Definition at line 35 of file AvalancheMicroscopic.hh.

35{ m_plotIonisations = true; }

◆ EnableMagneticField()

void Garfield::AvalancheMicroscopic::EnableMagneticField ( )
inline

Enable magnetic field in stepping algorithm (default: off).

Definition at line 103 of file AvalancheMicroscopic.hh.

103{ m_useBfield = true; }

◆ EnableNullCollisionSteps()

void Garfield::AvalancheMicroscopic::EnableNullCollisionSteps ( )
inline

Switch on update of coordinates for null-collision steps (default: off).

Definition at line 84 of file AvalancheMicroscopic.hh.

84{ m_useNullCollisionSteps = true; }

◆ EnablePhotonTransport()

void Garfield::AvalancheMicroscopic::EnablePhotonTransport ( )
inline

Switch on photon transport.

Remarks
This feature has not been tested thoroughly.

Definition at line 76 of file AvalancheMicroscopic.hh.

76{ m_usePhotons = true; }

◆ EnablePlotting()

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

Switch on drift line plotting.

Definition at line 94 of file AvalancheMicroscopic.cc.

94 {
95
96 if (!view) {
97 std::cerr << m_className << "::EnablePlotting:\n Null pointer.\n";
98 return;
99 }
100
101 m_viewer = view;
102 m_usePlotting = true;
103 if (!m_useDriftLines) {
104 std::cout << m_className << "::EnablePlotting:\n"
105 << " Enabling storage of drift line.\n";
107 }
108}
void EnableDriftLines()
Switch on storage of drift lines.

◆ EnableSecondaryEnergyHistogramming()

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.

212 {
213
214 if (!histo) {
215 std::cerr << m_className << "::EnableSecondaryEnergyHistogramming:\n"
216 << " Null pointer.\n";
217 return;
218 }
219
220 m_histSecondary = histo;
221}

◆ EnableSignalCalculation()

void Garfield::AvalancheMicroscopic::EnableSignalCalculation ( )
inline

Switch on calculation of induced currents.

Definition at line 42 of file AvalancheMicroscopic.hh.

42{ m_useSignal = true; }

◆ GetAvalancheSize() [1/2]

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

Definition at line 118 of file AvalancheMicroscopic.hh.

118 {
119 ne = m_nElectrons;
120 nh = m_nHoles;
121 ni = m_nIons;
122 }

◆ GetAvalancheSize() [2/2]

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

Return the number of electrons and ions in the avalanche.

Definition at line 114 of file AvalancheMicroscopic.hh.

114 {
115 ne = m_nElectrons;
116 ni = m_nIons;
117 }

Referenced by GarfieldPhysics::DoIt().

◆ GetAvalancheSizeLimit()

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

Definition at line 100 of file AvalancheMicroscopic.hh.

100{ return m_sizeCut; }

◆ GetElectronDriftLinePoint()

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.

356 {
357
358 if (iel >= m_endpointsElectrons.size()) {
359 std::cerr << m_className << "::GetElectronDriftLinePoint:\n";
360 std::cerr << " Endpoint index (" << iel << ") out of range.\n";
361 return;
362 }
363
364 if (ip <= 0) {
365 x = m_endpointsElectrons[iel].x0;
366 y = m_endpointsElectrons[iel].y0;
367 z = m_endpointsElectrons[iel].z0;
368 t = m_endpointsElectrons[iel].t0;
369 return;
370 }
371
372 const int np = m_endpointsElectrons[iel].driftLine.size();
373 if (ip > np) {
374 x = m_endpointsElectrons[iel].x;
375 y = m_endpointsElectrons[iel].y;
376 z = m_endpointsElectrons[iel].z;
377 t = m_endpointsElectrons[iel].t;
378 return;
379 }
380
381 x = m_endpointsElectrons[iel].driftLine[ip - 1].x;
382 y = m_endpointsElectrons[iel].driftLine[ip - 1].y;
383 z = m_endpointsElectrons[iel].driftLine[ip - 1].z;
384 t = m_endpointsElectrons[iel].driftLine[ip - 1].t;
385}

◆ GetElectronEndpoint() [1/2]

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.

270 {
271
272 if (i >= m_endpointsElectrons.size()) {
273 std::cerr << m_className << "::GetElectronEndpoint:\n";
274 std::cerr << " Endpoint index " << i << " out of range.\n";
275 x0 = y0 = z0 = t0 = e0 = 0.;
276 x1 = y1 = z1 = t1 = e1 = 0.;
277 dx1 = dy1 = dz1 = 0.;
278 status = 0;
279 return;
280 }
281
282 x0 = m_endpointsElectrons[i].x0;
283 y0 = m_endpointsElectrons[i].y0;
284 z0 = m_endpointsElectrons[i].z0;
285 t0 = m_endpointsElectrons[i].t0;
286 e0 = m_endpointsElectrons[i].e0;
287 x1 = m_endpointsElectrons[i].x;
288 y1 = m_endpointsElectrons[i].y;
289 z1 = m_endpointsElectrons[i].z;
290 t1 = m_endpointsElectrons[i].t;
291 e1 = m_endpointsElectrons[i].energy;
292 dx1 = m_endpointsElectrons[i].kx;
293 dy1 = m_endpointsElectrons[i].ky;
294 dz1 = m_endpointsElectrons[i].kz;
295 status = m_endpointsElectrons[i].status;
296}

◆ GetElectronEndpoint() [2/2]

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.

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
e0,e1initial and final energy
statusstatus code (see GarfieldConstants.hh)

Definition at line 238 of file AvalancheMicroscopic.cc.

243 {
244
245 if (i >= m_endpointsElectrons.size()) {
246 std::cerr << m_className << "::GetElectronEndpoint:\n";
247 std::cerr << " Endpoint index " << i << " out of range.\n";
248 x0 = y0 = z0 = t0 = e0 = 0.;
249 x1 = y1 = z1 = t1 = e1 = 0.;
250 status = 0;
251 return;
252 }
253
254 x0 = m_endpointsElectrons[i].x0;
255 y0 = m_endpointsElectrons[i].y0;
256 z0 = m_endpointsElectrons[i].z0;
257 t0 = m_endpointsElectrons[i].t0;
258 e0 = m_endpointsElectrons[i].e0;
259 x1 = m_endpointsElectrons[i].x;
260 y1 = m_endpointsElectrons[i].y;
261 z1 = m_endpointsElectrons[i].z;
262 t1 = m_endpointsElectrons[i].t;
263 e1 = m_endpointsElectrons[i].energy;
264 status = m_endpointsElectrons[i].status;
265}

◆ GetElectronTransportCut()

double Garfield::AvalancheMicroscopic::GetElectronTransportCut ( ) const
inline

Definition at line 90 of file AvalancheMicroscopic.hh.

90{ return m_deltaCut; }

◆ GetHoleDriftLinePoint()

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.

390 {
391
392 if (ih >= m_endpointsHoles.size()) {
393 std::cerr << m_className << "::GetHoleDriftLinePoint:\n";
394 std::cerr << " Endpoint index (" << ih << ") out of range.\n";
395 return;
396 }
397
398 if (ip <= 0) {
399 x = m_endpointsHoles[ih].x0;
400 y = m_endpointsHoles[ih].y0;
401 z = m_endpointsHoles[ih].z0;
402 t = m_endpointsHoles[ih].t0;
403 return;
404 }
405
406 const int np = m_endpointsHoles[ih].driftLine.size();
407 if (ip > np) {
408 x = m_endpointsHoles[ih].x;
409 y = m_endpointsHoles[ih].y;
410 z = m_endpointsHoles[ih].z;
411 t = m_endpointsHoles[ih].t;
412 return;
413 }
414
415 x = m_endpointsHoles[ih].driftLine[ip - 1].x;
416 y = m_endpointsHoles[ih].driftLine[ip - 1].y;
417 z = m_endpointsHoles[ih].driftLine[ip - 1].z;
418 t = m_endpointsHoles[ih].driftLine[ip - 1].t;
419}

◆ GetHoleEndpoint()

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.

302 {
303
304 if (i >= m_endpointsHoles.size()) {
305 std::cerr << m_className << "::GetHoleEndpoint:\n";
306 std::cerr << " Endpoint index " << i << " out of range.\n";
307 x0 = y0 = z0 = t0 = e0 = 0.;
308 x1 = y1 = z1 = t1 = e1 = 0.;
309 status = 0;
310 return;
311 }
312
313 x0 = m_endpointsHoles[i].x0;
314 y0 = m_endpointsHoles[i].y0;
315 z0 = m_endpointsHoles[i].z0;
316 t0 = m_endpointsHoles[i].t0;
317 e0 = m_endpointsHoles[i].e0;
318 x1 = m_endpointsHoles[i].x;
319 y1 = m_endpointsHoles[i].y;
320 z1 = m_endpointsHoles[i].z;
321 t1 = m_endpointsHoles[i].t;
322 e1 = m_endpointsHoles[i].energy;
323 status = m_endpointsHoles[i].status;
324}

◆ GetNumberOfElectronDriftLinePoints()

unsigned int Garfield::AvalancheMicroscopic::GetNumberOfElectronDriftLinePoints ( const unsigned int  i = 0) const

Definition at line 326 of file AvalancheMicroscopic.cc.

327 {
328
329 if (i >= m_endpointsElectrons.size()) {
330 std::cerr << m_className << "::GetNumberOfElectronDriftLinePoints:\n";
331 std::cerr << " Endpoint index (" << i << ") out of range.\n";
332 return 0;
333 }
334
335 if (!m_useDriftLines) return 2;
336
337 return m_endpointsElectrons[i].driftLine.size() + 2;
338}

◆ GetNumberOfElectronEndpoints()

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

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

Definition at line 126 of file AvalancheMicroscopic.hh.

126 {
127 return m_endpointsElectrons.size();
128 }

◆ GetNumberOfHoleDriftLinePoints()

unsigned int Garfield::AvalancheMicroscopic::GetNumberOfHoleDriftLinePoints ( const unsigned int  i = 0) const

Definition at line 340 of file AvalancheMicroscopic.cc.

340 {
341
342 if (i >= m_endpointsHoles.size()) {
343 std::cerr << m_className << "::GetNumberOfHoleDriftLinePoints:\n";
344 std::cerr << " Endpoint index (" << i << ") out of range.\n";
345 return 0;
346 }
347
348 if (!m_useDriftLines) return 2;
349
350 return m_endpointsHoles[i].driftLine.size() + 2;
351}

◆ GetNumberOfHoleEndpoints()

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

Definition at line 152 of file AvalancheMicroscopic.hh.

152 {
153 return m_endpointsHoles.size();
154 }

◆ GetNumberOfPhotons()

unsigned int Garfield::AvalancheMicroscopic::GetNumberOfPhotons ( ) const
inline

Definition at line 159 of file AvalancheMicroscopic.hh.

159{ return m_photons.size(); }

◆ GetPhoton()

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.

424 {
425
426 if (i >= m_photons.size()) {
427 std::cerr << m_className << "::GetPhoton:\n Index out of range.\n";
428 return;
429 }
430
431 x0 = m_photons[i].x0;
432 x1 = m_photons[i].x1;
433 y0 = m_photons[i].y0;
434 y1 = m_photons[i].y1;
435 z0 = m_photons[i].z0;
436 z1 = m_photons[i].z1;
437 t0 = m_photons[i].t0;
438 t1 = m_photons[i].t1;
439 status = m_photons[i].status;
440 e = m_photons[i].energy;
441}

◆ GetPhotonTransportCut()

double Garfield::AvalancheMicroscopic::GetPhotonTransportCut ( ) const
inline

Definition at line 94 of file AvalancheMicroscopic.hh.

94{ return m_gammaCut; }

◆ SetCollisionSteps()

void Garfield::AvalancheMicroscopic::SetCollisionSteps ( const unsigned int  n)
inline

Set number of collisions to be skipped for plotting.

Definition at line 107 of file AvalancheMicroscopic.hh.

107{ m_nCollSkip = n; }

◆ SetDistanceHistogram()

void Garfield::AvalancheMicroscopic::SetDistanceHistogram ( TH1 *  histo,
const char  opt = 'r' 
)

Fill histograms of the distance between successive collisions.

Parameters
histopointer to the histogram to be filled
optdirection ('x', 'y', 'z', 'r')

Definition at line 138 of file AvalancheMicroscopic.cc.

138 {
139
140 if (!histo) {
141 std::cerr << m_className << "::SetDistanceHistogram:\n Null pointer.\n";
142 return;
143 }
144
145 m_histDistance = histo;
146
147 if (opt == 'x' || opt == 'y' || opt == 'z' || opt == 'r') {
148 m_distanceOption = opt;
149 } else {
150 std::cerr << m_className << "::SetDistanceHistogram:";
151 std::cerr << " Unknown option " << opt << ".\n";
152 std::cerr << " Valid options are x, y, z, r.\n";
153 std::cerr << " Using default value (r).\n";
154 m_distanceOption = 'r';
155 }
156
157 if (m_distanceHistogramType.empty()) {
158 std::cout << m_className << "::SetDistanceHistogram:\n";
159 std::cout << " Don't forget to call EnableDistanceHistogramming.\n";
160 }
161}

◆ SetElectronTransportCut()

void Garfield::AvalancheMicroscopic::SetElectronTransportCut ( const double  cut)
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.

89{ m_deltaCut = cut; }

◆ SetPhotonTransportCut()

void Garfield::AvalancheMicroscopic::SetPhotonTransportCut ( const double  cut)
inline

Set an energy threshold for photon transport.

Definition at line 93 of file AvalancheMicroscopic.hh.

93{ m_gammaCut = cut; }

◆ SetSensor()

void Garfield::AvalancheMicroscopic::SetSensor ( Sensor sensor)

Set the sensor.

Definition at line 85 of file AvalancheMicroscopic.cc.

85 {
86
87 if (!s) {
88 std::cerr << m_className << "::SetSensor:\n Null pointer.\n";
89 return;
90 }
91 m_sensor = s;
92}

Referenced by GarfieldPhysics::InitializePhysics().

◆ SetTimeWindow()

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.

223 {
224
225 if (fabs(t1 - t0) < Small) {
226 std::cerr << m_className << "::SetTimeWindow:\n";
227 std::cerr << " Time interval must be greater than zero.\n";
228 return;
229 }
230
231 m_tMin = std::min(t0, t1);
232 m_tMax = std::max(t0, t1);
233 m_hasTimeWindow = true;
234}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

◆ SetUserHandleAttachment()

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.

463 {
464
465 m_userHandleAttachment = f;
466 m_hasUserHandleAttachment = true;
467}

◆ SetUserHandleInelastic()

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.

476 {
477
478 m_userHandleInelastic = f;
479 m_hasUserHandleInelastic = true;
480}

◆ SetUserHandleIonisation()

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.

489 {
490
491 m_userHandleIonisation = f;
492 m_hasUserHandleIonisation = true;
493}

◆ SetUserHandleStep()

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.

445 {
446
447 if (!f) {
448 std::cerr << m_className << "::SetUserHandleStep:\n";
449 std::cerr << " Function pointer is a null pointer.\n";
450 return;
451 }
452 m_userHandleStep = f;
453 m_hasUserHandleStep = true;
454}

◆ UnsetTimeWindow()

void Garfield::AvalancheMicroscopic::UnsetTimeWindow ( )

Definition at line 236 of file AvalancheMicroscopic.cc.

236{ m_hasTimeWindow = false; }

◆ UnsetUserHandleAttachment()

void Garfield::AvalancheMicroscopic::UnsetUserHandleAttachment ( )

Definition at line 469 of file AvalancheMicroscopic.cc.

469 {
470
471 m_userHandleAttachment = 0;
472 m_hasUserHandleAttachment = false;
473}

◆ UnsetUserHandleInelastic()

void Garfield::AvalancheMicroscopic::UnsetUserHandleInelastic ( )

Definition at line 482 of file AvalancheMicroscopic.cc.

482 {
483
484 m_userHandleInelastic = 0;
485 m_hasUserHandleInelastic = false;
486}

◆ UnsetUserHandleIonisation()

void Garfield::AvalancheMicroscopic::UnsetUserHandleIonisation ( )

Definition at line 495 of file AvalancheMicroscopic.cc.

495 {
496
497 m_userHandleIonisation = 0;
498 m_hasUserHandleIonisation = false;
499}

◆ UnsetUserHandleStep()

void Garfield::AvalancheMicroscopic::UnsetUserHandleStep ( )

Definition at line 456 of file AvalancheMicroscopic.cc.

456 {
457
458 m_userHandleStep = 0;
459 m_hasUserHandleStep = false;
460}

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