Garfield++ v1r0
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

#include <AvalancheMicroscopic.hh>

Public Member Functions

 AvalancheMicroscopic ()
 
 ~AvalancheMicroscopic ()
 
void SetSensor (Sensor *sensor)
 
void EnablePlotting (ViewDrift *view)
 
void DisablePlotting ()
 
void EnableExcitationMarkers ()
 
void DisableExcitationMarkers ()
 
void EnableIonisationMarkers ()
 
void DisableIonisationMarkers ()
 
void EnableAttachmentMarkers ()
 
void DisableAttachmentMarkers ()
 
void EnableSignalCalculation ()
 
void DisableSignalCalculation ()
 
void EnableInducedChargeCalculation ()
 
void DisableInducedChargeCalculation ()
 
void EnableElectronEnergyHistogramming (TH1 *histo)
 
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)
 
void DisableDistanceHistogramming ()
 
void EnableSecondaryEnergyHistogramming (TH1 *histo)
 
void DisableSecondaryEnergyHistogramming ()
 
void EnableDriftLines ()
 
void DisableDriftLines ()
 
void EnablePhotonTransport ()
 
void DisablePhotonTransport ()
 
void EnableBandStructure ()
 
void DisableBandStructure ()
 
void EnableNullCollisionSteps ()
 
void DisableNullCollisionSteps ()
 
void SetElectronTransportCut (const double cut)
 
double GetElectronTransportCut () const
 
void SetPhotonTransportCut (const double cut)
 
double GetPhotonTransportCut () const
 
void EnableAvalancheSizeLimit (const int size)
 
void DisableAvalancheSizeLimit ()
 
int GetAvalancheSizeLimit () const
 
void EnableMagneticField ()
 
void DisableMagneticField ()
 
void SetCollisionSteps (const int n)
 
void SetTimeWindow (const double t0, const double t1)
 
void UnsetTimeWindow ()
 
void GetAvalancheSize (int &ne, int &ni) const
 
void GetAvalancheSize (int &ne, int &nh, int &ni) const
 
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
 
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
 
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.)
 
void SetUserHandleStep (void(*f)(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole))
 
void UnsetUserHandleStep ()
 
void SetUserHandleAttachment (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m))
 
void UnsetUserHandleAttachment ()
 
void SetUserHandleInelastic (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m))
 
void UnsetUserHandleInelastic ()
 
void SetUserHandleIonisation (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m))
 
void UnsetUserHandleIonisation ()
 
void EnableDebugging ()
 
void DisableDebugging ()
 

Detailed Description

Definition at line 16 of file AvalancheMicroscopic.hh.

Constructor & Destructor Documentation

◆ AvalancheMicroscopic()

Garfield::AvalancheMicroscopic::AvalancheMicroscopic ( )

Definition at line 12 of file AvalancheMicroscopic.cc.

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

◆ ~AvalancheMicroscopic()

Garfield::AvalancheMicroscopic::~AvalancheMicroscopic ( )
inline

Definition at line 22 of file AvalancheMicroscopic.hh.

22{}

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

Definition at line 551 of file AvalancheMicroscopic.cc.

555 {
556
557 // Clear the list of electrons, holes and photons.
558 m_endpointsElectrons.clear();
559 m_endpointsHoles.clear();
560 m_photons.clear();
561
562 // Reset the particle counters.
563 m_nPhotons = m_nElectrons = m_nHoles = m_nIons = 0;
564 m_nElectronEndpoints = m_nHoleEndpoints = 0;
565
566 return TransportElectron(x0, y0, z0, t0, e0, dx0, dy0, dz0, true, false);
567}

Referenced by GarfieldPhysics::DoIt().

◆ DisableAttachmentMarkers()

void Garfield::AvalancheMicroscopic::DisableAttachmentMarkers ( )
inline

Definition at line 34 of file AvalancheMicroscopic.hh.

34{ m_plotAttachments = false; }

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMicroscopic::DisableAvalancheSizeLimit ( )
inline

Definition at line 86 of file AvalancheMicroscopic.hh.

86{ m_sizeCut = -1; }

◆ DisableBandStructure()

void Garfield::AvalancheMicroscopic::DisableBandStructure ( )
inline

Definition at line 69 of file AvalancheMicroscopic.hh.

69{ m_useBandStructureDefault = false; }

◆ DisableDebugging()

void Garfield::AvalancheMicroscopic::DisableDebugging ( )
inline

Definition at line 168 of file AvalancheMicroscopic.hh.

168{ m_debug = false; }

◆ DisableDistanceHistogramming() [1/2]

void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming ( )

Definition at line 219 of file AvalancheMicroscopic.cc.

219 {
220
221 m_hasDistanceHistogram = false;
222 m_distanceHistogramType.clear();
223}

◆ DisableDistanceHistogramming() [2/2]

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

Definition at line 198 of file AvalancheMicroscopic.cc.

198 {
199
200 if (m_distanceHistogramType.empty()) {
201 std::cerr << m_className << "::DisableDistanceHistogramming:\n";
202 std::cerr << " Collision type " << type << " is not histogrammed.\n";
203 return;
204 }
205 const unsigned int nDistanceHistogramTypes = m_distanceHistogramType.size();
206 for (int i = nDistanceHistogramTypes; i--;) {
207 if (m_distanceHistogramType[i] == type) {
208 m_distanceHistogramType.erase(m_distanceHistogramType.begin() + i);
209 std::cout << " Histogramming of collision type " << type
210 << " disabled.\n";
211 return;
212 }
213 }
214
215 std::cerr << m_className << "::DisableDistanceHistogramming:\n";
216 std::cerr << " Collision type " << type << " is not histogrammed.\n";
217}

◆ DisableDriftLines()

void Garfield::AvalancheMicroscopic::DisableDriftLines ( )
inline

Definition at line 61 of file AvalancheMicroscopic.hh.

61{ m_useDriftLines = false; }

◆ DisableElectronEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableElectronEnergyHistogramming ( )

Definition at line 125 of file AvalancheMicroscopic.cc.

125 {
126
127 m_hasElectronEnergyHistogram = false;
128}

◆ DisableExcitationMarkers()

void Garfield::AvalancheMicroscopic::DisableExcitationMarkers ( )
inline

Definition at line 30 of file AvalancheMicroscopic.hh.

30{ m_plotExcitations = false; }

◆ DisableHoleEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableHoleEnergyHistogramming ( )

Definition at line 142 of file AvalancheMicroscopic.cc.

142 {
143
144 m_hasHoleEnergyHistogram = false;
145}

◆ DisableInducedChargeCalculation()

void Garfield::AvalancheMicroscopic::DisableInducedChargeCalculation ( )
inline

Definition at line 42 of file AvalancheMicroscopic.hh.

42{ m_useInducedCharge = false; }

◆ DisableIonisationMarkers()

void Garfield::AvalancheMicroscopic::DisableIonisationMarkers ( )
inline

Definition at line 32 of file AvalancheMicroscopic.hh.

32{ m_plotIonisations = false; }

◆ DisableMagneticField()

void Garfield::AvalancheMicroscopic::DisableMagneticField ( )
inline

Definition at line 91 of file AvalancheMicroscopic.hh.

91{ m_useBfield = false; }

◆ DisableNullCollisionSteps()

void Garfield::AvalancheMicroscopic::DisableNullCollisionSteps ( )
inline

Definition at line 73 of file AvalancheMicroscopic.hh.

73{ m_useNullCollisionSteps = false; }

◆ DisablePhotonTransport()

void Garfield::AvalancheMicroscopic::DisablePhotonTransport ( )
inline

Definition at line 65 of file AvalancheMicroscopic.hh.

65{ m_usePhotons = false; }

◆ DisablePlotting()

void Garfield::AvalancheMicroscopic::DisablePlotting ( )

Definition at line 107 of file AvalancheMicroscopic.cc.

107 {
108
109 m_viewer = NULL;
110 m_usePlotting = false;
111}

◆ DisableSecondaryEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableSecondaryEnergyHistogramming ( )

Definition at line 237 of file AvalancheMicroscopic.cc.

237 {
238
239 m_hasSecondaryHistogram = false;
240}

◆ DisableSignalCalculation()

void Garfield::AvalancheMicroscopic::DisableSignalCalculation ( )
inline

Definition at line 38 of file AvalancheMicroscopic.hh.

38{ 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. 
)

Definition at line 534 of file AvalancheMicroscopic.cc.

537 {
538
539 // Clear the list of electrons and photons.
540 m_endpointsElectrons.clear();
541 m_endpointsHoles.clear();
542 m_photons.clear();
543
544 // Reset the particle counters.
545 m_nPhotons = m_nElectrons = m_nHoles = m_nIons = 0;
546 m_nElectronEndpoints = m_nHoleEndpoints = 0;
547
548 return TransportElectron(x0, y0, z0, t0, e0, dx0, dy0, dz0, false, false);
549}

◆ EnableAttachmentMarkers()

void Garfield::AvalancheMicroscopic::EnableAttachmentMarkers ( )
inline

Definition at line 33 of file AvalancheMicroscopic.hh.

33{ m_plotAttachments = true; }

◆ EnableAvalancheSizeLimit()

void Garfield::AvalancheMicroscopic::EnableAvalancheSizeLimit ( const int  size)
inline

Definition at line 85 of file AvalancheMicroscopic.hh.

85{ m_sizeCut = size; }

◆ EnableBandStructure()

void Garfield::AvalancheMicroscopic::EnableBandStructure ( )
inline

Definition at line 68 of file AvalancheMicroscopic.hh.

68{ m_useBandStructureDefault = true; }

◆ EnableDebugging()

void Garfield::AvalancheMicroscopic::EnableDebugging ( )
inline

Definition at line 167 of file AvalancheMicroscopic.hh.

167{ m_debug = true; }

◆ EnableDistanceHistogramming()

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

Definition at line 174 of file AvalancheMicroscopic.cc.

174 {
175
176 // Check if this type of collision is already registered
177 // for histogramming.
178 const unsigned int nDistanceHistogramTypes = m_distanceHistogramType.size();
179 if (nDistanceHistogramTypes > 0) {
180 for (int i = nDistanceHistogramTypes; i--;) {
181 if (m_distanceHistogramType[i] == type) {
182 std::cout << m_className << "::EnableDistanceHistogramming:\n";
183 std::cout << " Collision type " << type
184 << " is already histogrammed.\n";
185 return;
186 }
187 }
188 }
189
190 m_distanceHistogramType.push_back(type);
191 std::cout << m_className << "::EnableDistanceHistogramming:\n";
192 std::cout << " Histogramming of collision type " << type << " enabled.\n";
193 if (!m_hasDistanceHistogram) {
194 std::cout << " Don't forget to set the histogram.\n";
195 }
196}

◆ EnableDriftLines()

void Garfield::AvalancheMicroscopic::EnableDriftLines ( )
inline

Definition at line 60 of file AvalancheMicroscopic.hh.

60{ m_useDriftLines = true; }

Referenced by EnablePlotting().

◆ EnableElectronEnergyHistogramming()

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

Definition at line 113 of file AvalancheMicroscopic.cc.

113 {
114
115 if (!histo) {
116 std::cerr << m_className << "::EnableElectronEnergyHistogramming:\n";
117 std::cerr << " Histogram pointer is a null pointer.\n";
118 return;
119 }
120
121 m_histElectronEnergy = histo;
122 m_hasElectronEnergyHistogram = true;
123}

◆ EnableExcitationMarkers()

void Garfield::AvalancheMicroscopic::EnableExcitationMarkers ( )
inline

Definition at line 29 of file AvalancheMicroscopic.hh.

29{ m_plotExcitations = true; }

◆ EnableHoleEnergyHistogramming()

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

Definition at line 130 of file AvalancheMicroscopic.cc.

130 {
131
132 if (!histo) {
133 std::cerr << m_className << "::EnableHoleEnergyHistogramming:\n";
134 std::cerr << " Histogram pointer is a null pointer.\n";
135 return;
136 }
137
138 m_histHoleEnergy = histo;
139 m_hasHoleEnergyHistogram = true;
140}

◆ EnableInducedChargeCalculation()

void Garfield::AvalancheMicroscopic::EnableInducedChargeCalculation ( )
inline

Definition at line 41 of file AvalancheMicroscopic.hh.

41{ m_useInducedCharge = true; }

◆ EnableIonisationMarkers()

void Garfield::AvalancheMicroscopic::EnableIonisationMarkers ( )
inline

Definition at line 31 of file AvalancheMicroscopic.hh.

31{ m_plotIonisations = true; }

◆ EnableMagneticField()

void Garfield::AvalancheMicroscopic::EnableMagneticField ( )
inline

Definition at line 90 of file AvalancheMicroscopic.hh.

90{ m_useBfield = true; }

◆ EnableNullCollisionSteps()

void Garfield::AvalancheMicroscopic::EnableNullCollisionSteps ( )
inline

Definition at line 72 of file AvalancheMicroscopic.hh.

72{ m_useNullCollisionSteps = true; }

◆ EnablePhotonTransport()

void Garfield::AvalancheMicroscopic::EnablePhotonTransport ( )
inline

Definition at line 64 of file AvalancheMicroscopic.hh.

64{ m_usePhotons = true; }

◆ EnablePlotting()

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

Definition at line 90 of file AvalancheMicroscopic.cc.

90 {
91
92 if (!view) {
93 std::cerr << m_className << "::EnablePlotting:\n";
94 std::cerr << " Viewer pointer is a null pointer.\n";
95 return;
96 }
97
98 m_viewer = view;
99 m_usePlotting = true;
100 if (!m_useDriftLines) {
101 std::cout << m_className << "::EnablePlotting:\n";
102 std::cout << " Enabling storage of drift line.\n";
104 }
105}

◆ EnableSecondaryEnergyHistogramming()

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

Definition at line 225 of file AvalancheMicroscopic.cc.

225 {
226
227 if (!histo) {
228 std::cerr << m_className << "::EnableSecondaryEnergyHistogramming:\n";
229 std::cerr << " Histogram pointer is a null pointer.\n";
230 return;
231 }
232
233 m_histSecondary = histo;
234 m_hasSecondaryHistogram = true;
235}

◆ EnableSignalCalculation()

void Garfield::AvalancheMicroscopic::EnableSignalCalculation ( )
inline

Definition at line 37 of file AvalancheMicroscopic.hh.

37{ m_useSignal = true; }

◆ GetAvalancheSize() [1/2]

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

Definition at line 103 of file AvalancheMicroscopic.hh.

103 {
104 ne = m_nElectrons;
105 nh = m_nHoles;
106 ni = m_nIons;
107 }

◆ GetAvalancheSize() [2/2]

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

Definition at line 99 of file AvalancheMicroscopic.hh.

99 {
100 ne = m_nElectrons;
101 ni = m_nIons;
102 }

Referenced by GarfieldPhysics::DoIt().

◆ GetAvalancheSizeLimit()

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

Definition at line 87 of file AvalancheMicroscopic.hh.

87{ 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 385 of file AvalancheMicroscopic.cc.

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

◆ 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 299 of file AvalancheMicroscopic.cc.

302 {
303
304 if (i >= m_nElectronEndpoints) {
305 std::cerr << m_className << "::GetElectronEndpoint:\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 dx1 = dy1 = dz1 = 0.;
310 status = 0;
311 return;
312 }
313
314 x0 = m_endpointsElectrons[i].x0;
315 y0 = m_endpointsElectrons[i].y0;
316 z0 = m_endpointsElectrons[i].z0;
317 t0 = m_endpointsElectrons[i].t0;
318 e0 = m_endpointsElectrons[i].e0;
319 x1 = m_endpointsElectrons[i].x;
320 y1 = m_endpointsElectrons[i].y;
321 z1 = m_endpointsElectrons[i].z;
322 t1 = m_endpointsElectrons[i].t;
323 e1 = m_endpointsElectrons[i].energy;
324 dx1 = m_endpointsElectrons[i].kx;
325 dy1 = m_endpointsElectrons[i].ky;
326 dz1 = m_endpointsElectrons[i].kz;
327 status = m_endpointsElectrons[i].status;
328}

◆ 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

Definition at line 270 of file AvalancheMicroscopic.cc.

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

◆ GetElectronTransportCut()

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

Definition at line 78 of file AvalancheMicroscopic.hh.

78{ 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 419 of file AvalancheMicroscopic.cc.

422 {
423
424 if (ih >= m_nHoleEndpoints) {
425 std::cerr << m_className << "::GetHoleDriftLinePoint:\n";
426 std::cerr << " Endpoint index (" << ih << ") out of range.\n";
427 return;
428 }
429
430 if (ip <= 0) {
431 x = m_endpointsHoles[ih].x0;
432 y = m_endpointsHoles[ih].y0;
433 z = m_endpointsHoles[ih].z0;
434 t = m_endpointsHoles[ih].t0;
435 return;
436 }
437
438 const int np = m_endpointsHoles[ih].driftLine.size();
439 if (ip > np) {
440 x = m_endpointsHoles[ih].x;
441 y = m_endpointsHoles[ih].y;
442 z = m_endpointsHoles[ih].z;
443 t = m_endpointsHoles[ih].t;
444 return;
445 }
446
447 x = m_endpointsHoles[ih].driftLine[ip - 1].x;
448 y = m_endpointsHoles[ih].driftLine[ip - 1].y;
449 z = m_endpointsHoles[ih].driftLine[ip - 1].z;
450 t = m_endpointsHoles[ih].driftLine[ip - 1].t;
451}

◆ 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 330 of file AvalancheMicroscopic.cc.

334 {
335
336 if (i >= m_nHoleEndpoints) {
337 std::cerr << m_className << "::GetHoleEndpoint:\n";
338 std::cerr << " Endpoint index " << i << " out of range.\n";
339 x0 = y0 = z0 = t0 = e0 = 0.;
340 x1 = y1 = z1 = t1 = e1 = 0.;
341 status = 0;
342 return;
343 }
344
345 x0 = m_endpointsHoles[i].x0;
346 y0 = m_endpointsHoles[i].y0;
347 z0 = m_endpointsHoles[i].z0;
348 t0 = m_endpointsHoles[i].t0;
349 e0 = m_endpointsHoles[i].e0;
350 x1 = m_endpointsHoles[i].x;
351 y1 = m_endpointsHoles[i].y;
352 z1 = m_endpointsHoles[i].z;
353 t1 = m_endpointsHoles[i].t;
354 e1 = m_endpointsHoles[i].energy;
355 status = m_endpointsHoles[i].status;
356}

◆ GetNumberOfElectronDriftLinePoints()

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

Definition at line 358 of file AvalancheMicroscopic.cc.

359 {
360
361 if (i >= m_nElectronEndpoints) {
362 std::cerr << m_className << "::GetNumberOfElectronDriftLinePoints:\n";
363 std::cerr << " Endpoint index (" << i << ") out of range.\n";
364 return 0;
365 }
366
367 if (!m_useDriftLines) return 2;
368
369 return m_endpointsElectrons[i].driftLine.size() + 2;
370}

◆ GetNumberOfElectronEndpoints()

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

Definition at line 109 of file AvalancheMicroscopic.hh.

109{ return m_nElectronEndpoints; }

◆ GetNumberOfHoleDriftLinePoints()

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

Definition at line 372 of file AvalancheMicroscopic.cc.

372 {
373
374 if (i >= m_nHoleEndpoints) {
375 std::cerr << m_className << "::GetNumberOfHoleDriftLinePoints:\n";
376 std::cerr << " Endpoint index (" << i << ") out of range.\n";
377 return 0;
378 }
379
380 if (!m_useDriftLines) return 2;
381
382 return m_endpointsHoles[i].driftLine.size() + 2;
383}

◆ GetNumberOfHoleEndpoints()

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

Definition at line 125 of file AvalancheMicroscopic.hh.

125{ return m_nHoleEndpoints; }

◆ GetNumberOfPhotons()

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

Definition at line 130 of file AvalancheMicroscopic.hh.

130{ return m_usePhotons ? m_nPhotons : 0; }

◆ 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 453 of file AvalancheMicroscopic.cc.

456 {
457
458 if (i >= m_nPhotons) {
459 std::cerr << m_className << "::GetPhoton:\n";
460 std::cerr << " Photon " << i << " does not exist.\n";
461 return;
462 }
463
464 x0 = m_photons[i].x0;
465 x1 = m_photons[i].x1;
466 y0 = m_photons[i].y0;
467 y1 = m_photons[i].y1;
468 z0 = m_photons[i].z0;
469 z1 = m_photons[i].z1;
470 t0 = m_photons[i].t0;
471 t1 = m_photons[i].t1;
472 status = m_photons[i].status;
473 e = m_photons[i].energy;
474}

◆ GetPhotonTransportCut()

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

Definition at line 82 of file AvalancheMicroscopic.hh.

82{ return m_gammaCut; }

◆ SetCollisionSteps()

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

Definition at line 242 of file AvalancheMicroscopic.cc.

242 {
243
244 if (n <= 0) {
245 std::cerr << m_className << "::SetCollisionSteps:\n";
246 std::cerr << " Number of collisions to be skipped set to"
247 << " default value (100).\n";
248 m_nCollSkip = 100;
249 return;
250 }
251
252 m_nCollSkip = n;
253}

◆ SetDistanceHistogram()

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

Definition at line 147 of file AvalancheMicroscopic.cc.

147 {
148
149 if (!histo) {
150 std::cerr << m_className << "::SetDistanceHistogram:\n";
151 std::cerr << " Histogram pointer is a null pointer.\n";
152 return;
153 }
154
155 m_histDistance = histo;
156 m_hasDistanceHistogram = true;
157
158 if (opt == 'x' || opt == 'y' || opt == 'z' || opt == 'r') {
159 m_distanceOption = opt;
160 } else {
161 std::cerr << m_className << "::SetDistanceHistogram:";
162 std::cerr << " Unknown option " << opt << ".\n";
163 std::cerr << " Valid options are x, y, z, r.\n";
164 std::cerr << " Using default value (r).\n";
165 m_distanceOption = 'r';
166 }
167
168 if (m_distanceHistogramType.empty()) {
169 std::cout << m_className << "::SetDistanceHistogram:\n";
170 std::cout << " Don't forget to call EnableDistanceHistogramming.\n";
171 }
172}

◆ SetElectronTransportCut()

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

Definition at line 77 of file AvalancheMicroscopic.hh.

77{ m_deltaCut = cut; }

◆ SetPhotonTransportCut()

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

Definition at line 81 of file AvalancheMicroscopic.hh.

81{ m_gammaCut = cut; }

◆ SetSensor()

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

Definition at line 80 of file AvalancheMicroscopic.cc.

80 {
81
82 if (!s) {
83 std::cerr << m_className << "::SetSensor:\n";
84 std::cerr << " Sensor pointer is a null pointer.\n";
85 return;
86 }
87 m_sensor = s;
88}

Referenced by GarfieldPhysics::InitializePhysics().

◆ SetTimeWindow()

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

Definition at line 255 of file AvalancheMicroscopic.cc.

255 {
256
257 if (fabs(t1 - t0) < Small) {
258 std::cerr << m_className << "::SetTimeWindow:\n";
259 std::cerr << " Time interval must be greater than zero.\n";
260 return;
261 }
262
263 m_tMin = std::min(t0, t1);
264 m_tMax = std::max(t0, t1);
265 m_hasTimeWindow = true;
266}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:616

◆ SetUserHandleAttachment()

void Garfield::AvalancheMicroscopic::SetUserHandleAttachment ( void(*)(double x, double y, double z, double t, int type, int level, Medium *m)  f)

Definition at line 495 of file AvalancheMicroscopic.cc.

496 {
497
498 m_userHandleAttachment = f;
499 m_hasUserHandleAttachment = true;
500}

◆ SetUserHandleInelastic()

void Garfield::AvalancheMicroscopic::SetUserHandleInelastic ( void(*)(double x, double y, double z, double t, int type, int level, Medium *m)  f)

Definition at line 508 of file AvalancheMicroscopic.cc.

509 {
510
511 m_userHandleInelastic = f;
512 m_hasUserHandleInelastic = true;
513}

◆ SetUserHandleIonisation()

void Garfield::AvalancheMicroscopic::SetUserHandleIonisation ( void(*)(double x, double y, double z, double t, int type, int level, Medium *m)  f)

Definition at line 521 of file AvalancheMicroscopic.cc.

522 {
523
524 m_userHandleIonisation = f;
525 m_hasUserHandleIonisation = true;
526}

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

Definition at line 476 of file AvalancheMicroscopic.cc.

478 {
479
480 if (!f) {
481 std::cerr << m_className << "::SetUserHandleStep:\n";
482 std::cerr << " Function pointer is a null pointer.\n";
483 return;
484 }
485 m_userHandleStep = f;
486 m_hasUserHandleStep = true;
487}

◆ UnsetTimeWindow()

void Garfield::AvalancheMicroscopic::UnsetTimeWindow ( )

Definition at line 268 of file AvalancheMicroscopic.cc.

268{ m_hasTimeWindow = false; }

◆ UnsetUserHandleAttachment()

void Garfield::AvalancheMicroscopic::UnsetUserHandleAttachment ( )

Definition at line 502 of file AvalancheMicroscopic.cc.

502 {
503
504 m_userHandleAttachment = 0;
505 m_hasUserHandleAttachment = false;
506}

◆ UnsetUserHandleInelastic()

void Garfield::AvalancheMicroscopic::UnsetUserHandleInelastic ( )

Definition at line 515 of file AvalancheMicroscopic.cc.

515 {
516
517 m_userHandleInelastic = 0;
518 m_hasUserHandleInelastic = false;
519}

◆ UnsetUserHandleIonisation()

void Garfield::AvalancheMicroscopic::UnsetUserHandleIonisation ( )

Definition at line 528 of file AvalancheMicroscopic.cc.

528 {
529
530 m_userHandleIonisation = 0;
531 m_hasUserHandleIonisation = false;
532}

◆ UnsetUserHandleStep()

void Garfield::AvalancheMicroscopic::UnsetUserHandleStep ( )

Definition at line 489 of file AvalancheMicroscopic.cc.

489 {
490
491 m_userHandleStep = 0;
492 m_hasUserHandleStep = false;
493}

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