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::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 ()
 Switch off drift line plotting.
 
void EnableExcitationMarkers (const bool on=true)
 Draw a marker at every excitation or not.
 
void EnableIonisationMarkers (const bool on=true)
 Draw a marker at every ionising collision or not.
 
void EnableAttachmentMarkers (const bool on=true)
 Draw a marker at every attachment or not.
 
void EnableSignalCalculation (const bool on=true)
 Switch on calculation of induced currents (default: off).
 
void UseWeightingPotential (const bool on=true)
 
void EnableWeightingFieldIntegration (const bool on=true)
 
void EnableInducedChargeCalculation (const bool on=true)
 Switch on calculation of the total induced charge (default: off).
 
void EnableElectronEnergyHistogramming (TH1 *histo)
 Fill a histogram with the electron energy distribution.
 
void DisableElectronEnergyHistogramming ()
 Stop histogramming the electron energy distribution.
 
void EnableHoleEnergyHistogramming (TH1 *histo)
 Fill a histogram with the hole energy distribution.
 
void DisableHoleEnergyHistogramming ()
 Stop histogramming the hole energy distribution.
 
void SetDistanceHistogram (TH1 *histo, const char opt='r')
 
void EnableDistanceHistogramming (const int type)
 Fill distance distribution histograms for a given collision type.
 
void DisableDistanceHistogramming (const int type)
 Stop filling distance distribution histograms for a given collision type.
 
void DisableDistanceHistogramming ()
 Stop filling distance distribution histograms.
 
void EnableSecondaryEnergyHistogramming (TH1 *histo)
 Fill histograms of the energy of electrons emitted in ionising collisions.
 
void DisableSecondaryEnergyHistogramming ()
 Stop histogramming the secondary electron energy distribution.
 
void EnableDriftLines (const bool on=true)
 Switch on storage of drift lines (default: off).
 
void EnablePhotonTransport (const bool on=true)
 
void EnableBandStructure (const bool on=true)
 Switch on stepping according to band structure E(k), for semiconductors.
 
void EnableNullCollisionSteps (const bool on=true)
 Switch on update of coordinates for null-collision steps (default: off).
 
void SetElectronTransportCut (const double cut)
 
double GetElectronTransportCut () const
 Retrieve the value of the energy threshold.
 
void SetPhotonTransportCut (const double cut)
 Set an energy threshold for photon transport.
 
double GetPhotonTransportCut () const
 Retrieve the energy threshold for transporting photons.
 
void EnableAvalancheSizeLimit (const unsigned int size)
 
void DisableAvalancheSizeLimit ()
 Do not apply a limit on the avalanche size.
 
int GetAvalancheSizeLimit () const
 Retrieve the currently set size limit.
 
void EnableMagneticField (const bool on=true)
 Enable magnetic field in stepping algorithm (default: off).
 
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 ()
 Do not restrict the time interval within which carriers are simulated.
 
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 ()
 Deactivate the user handle called at every step.
 
void SetUserHandleCollision (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m, double e0, double e1, double dx0, double dy0, double dz0, double dx1, double dy1, double dz1))
 Set a user handling procedure, to be called at every (real) collision.
 
void UnsetUserHandleCollision ()
 Deactivate the user handle called at every collision.
 
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 ()
 Deactivate the user handle called at every attachment.
 
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 ()
 Deactivate the user handle called at every inelastic collision.
 
void SetUserHandleIonisation (void(*f)(double x, double y, double z, double t, int type, int level, Medium *m))
 
void UnsetUserHandleIonisation ()
 Deactivate the user handle called at every ionisation.
 
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 91 of file AvalancheMicroscopic.cc.

91 {
92 m_endpointsElectrons.reserve(10000);
93 m_endpointsHoles.reserve(10000);
94 m_photons.reserve(1000);
95
96}

◆ ~AvalancheMicroscopic()

Garfield::AvalancheMicroscopic::~AvalancheMicroscopic ( )
inline

Destructor.

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

Calculate an avalanche initiated by a given electron.

Definition at line 469 of file AvalancheMicroscopic.cc.

473 {
474 // Clear the list of electrons, holes and photons.
475 m_endpointsElectrons.clear();
476 m_endpointsHoles.clear();
477 m_photons.clear();
478
479 // Reset the particle counters.
480 m_nElectrons = m_nHoles = m_nIons = 0;
481
482 return TransportElectron(x0, y0, z0, t0, e0, dx0, dy0, dz0, true, false);
483}

Referenced by GarfieldPhysics::DoIt(), and main().

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMicroscopic::DisableAvalancheSizeLimit ( )
inline

Do not apply a limit on the avalanche size.

Definition at line 121 of file AvalancheMicroscopic.hh.

121{ m_sizeCut = 0; }

◆ DisableDebugging()

void Garfield::AvalancheMicroscopic::DisableDebugging ( )
inline

Definition at line 241 of file AvalancheMicroscopic.hh.

241{ m_debug = false; }

◆ DisableDistanceHistogramming() [1/2]

void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming ( )

Stop filling distance distribution histograms.

Definition at line 200 of file AvalancheMicroscopic.cc.

200 {
201 m_histDistance = nullptr;
202 m_distanceHistogramType.clear();
203}

◆ DisableDistanceHistogramming() [2/2]

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

Stop filling distance distribution histograms for a given collision type.

Definition at line 186 of file AvalancheMicroscopic.cc.

186 {
187 if (std::find(m_distanceHistogramType.begin(), m_distanceHistogramType.end(),
188 type) == m_distanceHistogramType.end()) {
189 std::cerr << m_className << "::DisableDistanceHistogramming:\n"
190 << " Collision type " << type << " is not histogrammed.\n";
191 return;
192 }
193
194 m_distanceHistogramType.erase(
195 std::remove(m_distanceHistogramType.begin(),
196 m_distanceHistogramType.end(), type),
197 m_distanceHistogramType.end());
198}

◆ DisableElectronEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableElectronEnergyHistogramming ( )
inline

Stop histogramming the electron energy distribution.

Definition at line 65 of file AvalancheMicroscopic.hh.

65{ m_histElectronEnergy = nullptr; }

◆ DisableHoleEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableHoleEnergyHistogramming ( )
inline

Stop histogramming the hole energy distribution.

Definition at line 69 of file AvalancheMicroscopic.hh.

69{ m_histHoleEnergy = nullptr; }

◆ DisablePlotting()

void Garfield::AvalancheMicroscopic::DisablePlotting ( )
inline

Switch off drift line plotting.

Definition at line 30 of file AvalancheMicroscopic.hh.

30{ m_viewer = nullptr; }

◆ DisableSecondaryEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableSecondaryEnergyHistogramming ( )
inline

Stop histogramming the secondary electron energy distribution.

Definition at line 87 of file AvalancheMicroscopic.hh.

87{ m_histSecondary = nullptr; }

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

457 {
458 // Clear the list of electrons and photons.
459 m_endpointsElectrons.clear();
460 m_endpointsHoles.clear();
461 m_photons.clear();
462
463 // Reset the particle counters.
464 m_nElectrons = m_nHoles = m_nIons = 0;
465
466 return TransportElectron(x0, y0, z0, t0, e0, dx0, dy0, dz0, false, false);
467}

◆ EnableAttachmentMarkers()

void Garfield::AvalancheMicroscopic::EnableAttachmentMarkers ( const bool  on = true)
inline

Draw a marker at every attachment or not.

Definition at line 40 of file AvalancheMicroscopic.hh.

40 {
41 m_plotAttachments = on;
42 }

◆ 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 119 of file AvalancheMicroscopic.hh.

119{ m_sizeCut = size; }

◆ EnableBandStructure()

void Garfield::AvalancheMicroscopic::EnableBandStructure ( const bool  on = true)
inline

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

Definition at line 97 of file AvalancheMicroscopic.hh.

97 {
98 m_useBandStructureDefault = on;
99 }

◆ EnableDebugging()

void Garfield::AvalancheMicroscopic::EnableDebugging ( )
inline

Switch on debugging messages.

Definition at line 240 of file AvalancheMicroscopic.hh.

240{ m_debug = true; }

◆ EnableDistanceHistogramming()

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

Fill distance distribution histograms for a given collision type.

Definition at line 164 of file AvalancheMicroscopic.cc.

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 ( const bool  on = true)
inline

Switch on storage of drift lines (default: off).

Definition at line 90 of file AvalancheMicroscopic.hh.

90{ m_useDriftLines = on; }

Referenced by EnablePlotting().

◆ EnableElectronEnergyHistogramming()

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

Fill a histogram with the electron energy distribution.

Definition at line 120 of file AvalancheMicroscopic.cc.

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

◆ EnableExcitationMarkers()

void Garfield::AvalancheMicroscopic::EnableExcitationMarkers ( const bool  on = true)
inline

Draw a marker at every excitation or not.

Definition at line 32 of file AvalancheMicroscopic.hh.

32 {
33 m_plotExcitations = on;
34 }

◆ EnableHoleEnergyHistogramming()

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

Fill a histogram with the hole energy distribution.

Definition at line 130 of file AvalancheMicroscopic.cc.

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

◆ EnableInducedChargeCalculation()

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

Switch on calculation of the total induced charge (default: off).

Definition at line 58 of file AvalancheMicroscopic.hh.

58 {
59 m_doInducedCharge = on;
60 }

◆ EnableIonisationMarkers()

void Garfield::AvalancheMicroscopic::EnableIonisationMarkers ( const bool  on = true)
inline

Draw a marker at every ionising collision or not.

Definition at line 36 of file AvalancheMicroscopic.hh.

36 {
37 m_plotIonisations = on;
38 }

◆ EnableMagneticField()

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

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

Definition at line 126 of file AvalancheMicroscopic.hh.

126{ m_useBfield = on; }

◆ EnableNullCollisionSteps()

void Garfield::AvalancheMicroscopic::EnableNullCollisionSteps ( const bool  on = true)
inline

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

Definition at line 102 of file AvalancheMicroscopic.hh.

102 {
103 m_useNullCollisionSteps = on;
104 }

◆ EnablePhotonTransport()

void Garfield::AvalancheMicroscopic::EnablePhotonTransport ( const bool  on = true)
inline

Switch on photon transport.

Remarks
This feature has not been tested thoroughly.

Definition at line 94 of file AvalancheMicroscopic.hh.

94{ m_usePhotons = on; }

◆ EnablePlotting()

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

Switch on drift line plotting.

Definition at line 106 of file AvalancheMicroscopic.cc.

106 {
107 if (!view) {
108 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
109 return;
110 }
111
112 m_viewer = view;
113 if (!m_useDriftLines) {
114 std::cout << m_className << "::EnablePlotting:\n"
115 << " Enabling storage of drift line.\n";
117 }
118}
void EnableDriftLines(const bool on=true)
Switch on storage of drift lines (default: off).

Referenced by main().

◆ EnableSecondaryEnergyHistogramming()

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

Fill histograms of the energy of electrons emitted in ionising collisions.

Definition at line 205 of file AvalancheMicroscopic.cc.

205 {
206 if (!histo) {
207 std::cerr << m_className << "::EnableSecondaryEnergyHistogramming:\n"
208 << " Null pointer.\n";
209 return;
210 }
211
212 m_histSecondary = histo;
213}

◆ EnableSignalCalculation()

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

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

Definition at line 45 of file AvalancheMicroscopic.hh.

45{ m_doSignal = on; }

◆ EnableWeightingFieldIntegration()

void Garfield::AvalancheMicroscopic::EnableWeightingFieldIntegration ( const bool  on = true)
inline

Integrate the weighting field over a drift line step when calculating the induced current (default: off).

Definition at line 53 of file AvalancheMicroscopic.hh.

53 {
54 m_integrateWeightingField = on;
55 }

◆ GetAvalancheSize() [1/2]

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

Definition at line 141 of file AvalancheMicroscopic.hh.

141 {
142 ne = m_nElectrons;
143 nh = m_nHoles;
144 ni = m_nIons;
145 }

◆ 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 137 of file AvalancheMicroscopic.hh.

137 {
138 ne = m_nElectrons;
139 ni = m_nIons;
140 }

Referenced by GarfieldPhysics::DoIt().

◆ GetAvalancheSizeLimit()

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

Retrieve the currently set size limit.

Definition at line 123 of file AvalancheMicroscopic.hh.

123{ 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 335 of file AvalancheMicroscopic.cc.

337 {
338 if (iel >= m_endpointsElectrons.size()) {
339 std::cerr << m_className << "::GetElectronDriftLinePoint:\n";
340 std::cerr << " Endpoint index (" << iel << ") out of range.\n";
341 return;
342 }
343
344 if (ip <= 0) {
345 x = m_endpointsElectrons[iel].x0;
346 y = m_endpointsElectrons[iel].y0;
347 z = m_endpointsElectrons[iel].z0;
348 t = m_endpointsElectrons[iel].t0;
349 return;
350 }
351
352 const int np = m_endpointsElectrons[iel].driftLine.size();
353 if (ip > np) {
354 x = m_endpointsElectrons[iel].x;
355 y = m_endpointsElectrons[iel].y;
356 z = m_endpointsElectrons[iel].z;
357 t = m_endpointsElectrons[iel].t;
358 return;
359 }
360
361 x = m_endpointsElectrons[iel].driftLine[ip - 1].x;
362 y = m_endpointsElectrons[iel].driftLine[ip - 1].y;
363 z = m_endpointsElectrons[iel].driftLine[ip - 1].z;
364 t = m_endpointsElectrons[iel].driftLine[ip - 1].t;
365}

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

257 {
258 if (i >= m_endpointsElectrons.size()) {
259 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
260 x0 = y0 = z0 = t0 = e0 = 0.;
261 x1 = y1 = z1 = t1 = e1 = 0.;
262 dx1 = dy1 = dz1 = 0.;
263 status = 0;
264 return;
265 }
266
267 x0 = m_endpointsElectrons[i].x0;
268 y0 = m_endpointsElectrons[i].y0;
269 z0 = m_endpointsElectrons[i].z0;
270 t0 = m_endpointsElectrons[i].t0;
271 e0 = m_endpointsElectrons[i].e0;
272 x1 = m_endpointsElectrons[i].x;
273 y1 = m_endpointsElectrons[i].y;
274 z1 = m_endpointsElectrons[i].z;
275 t1 = m_endpointsElectrons[i].t;
276 e1 = m_endpointsElectrons[i].energy;
277 dx1 = m_endpointsElectrons[i].kx;
278 dy1 = m_endpointsElectrons[i].ky;
279 dz1 = m_endpointsElectrons[i].kz;
280 status = m_endpointsElectrons[i].status;
281}

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

232 {
233 if (i >= m_endpointsElectrons.size()) {
234 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
235 x0 = y0 = z0 = t0 = e0 = 0.;
236 x1 = y1 = z1 = t1 = e1 = 0.;
237 status = 0;
238 return;
239 }
240
241 x0 = m_endpointsElectrons[i].x0;
242 y0 = m_endpointsElectrons[i].y0;
243 z0 = m_endpointsElectrons[i].z0;
244 t0 = m_endpointsElectrons[i].t0;
245 e0 = m_endpointsElectrons[i].e0;
246 x1 = m_endpointsElectrons[i].x;
247 y1 = m_endpointsElectrons[i].y;
248 z1 = m_endpointsElectrons[i].z;
249 t1 = m_endpointsElectrons[i].t;
250 e1 = m_endpointsElectrons[i].energy;
251 status = m_endpointsElectrons[i].status;
252}

◆ GetElectronTransportCut()

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

Retrieve the value of the energy threshold.

Definition at line 110 of file AvalancheMicroscopic.hh.

110{ 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 367 of file AvalancheMicroscopic.cc.

370 {
371 if (ih >= m_endpointsHoles.size()) {
372 std::cerr << m_className << "::GetHoleDriftLinePoint:\n";
373 std::cerr << " Endpoint index (" << ih << ") out of range.\n";
374 return;
375 }
376
377 if (ip <= 0) {
378 x = m_endpointsHoles[ih].x0;
379 y = m_endpointsHoles[ih].y0;
380 z = m_endpointsHoles[ih].z0;
381 t = m_endpointsHoles[ih].t0;
382 return;
383 }
384
385 const int np = m_endpointsHoles[ih].driftLine.size();
386 if (ip > np) {
387 x = m_endpointsHoles[ih].x;
388 y = m_endpointsHoles[ih].y;
389 z = m_endpointsHoles[ih].z;
390 t = m_endpointsHoles[ih].t;
391 return;
392 }
393
394 x = m_endpointsHoles[ih].driftLine[ip - 1].x;
395 y = m_endpointsHoles[ih].driftLine[ip - 1].y;
396 z = m_endpointsHoles[ih].driftLine[ip - 1].z;
397 t = m_endpointsHoles[ih].driftLine[ip - 1].t;
398}

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

287 {
288 if (i >= m_endpointsHoles.size()) {
289 std::cerr << m_className << "::GetHoleEndpoint: Index out of range.\n";
290 x0 = y0 = z0 = t0 = e0 = 0.;
291 x1 = y1 = z1 = t1 = e1 = 0.;
292 status = 0;
293 return;
294 }
295
296 x0 = m_endpointsHoles[i].x0;
297 y0 = m_endpointsHoles[i].y0;
298 z0 = m_endpointsHoles[i].z0;
299 t0 = m_endpointsHoles[i].t0;
300 e0 = m_endpointsHoles[i].e0;
301 x1 = m_endpointsHoles[i].x;
302 y1 = m_endpointsHoles[i].y;
303 z1 = m_endpointsHoles[i].z;
304 t1 = m_endpointsHoles[i].t;
305 e1 = m_endpointsHoles[i].energy;
306 status = m_endpointsHoles[i].status;
307}

◆ GetNumberOfElectronDriftLinePoints()

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

Definition at line 309 of file AvalancheMicroscopic.cc.

310 {
311 if (i >= m_endpointsElectrons.size()) {
312 std::cerr << m_className << "::GetNumberOfElectronDriftLinePoints:\n";
313 std::cerr << " Endpoint index (" << i << ") out of range.\n";
314 return 0;
315 }
316
317 if (!m_useDriftLines) return 2;
318
319 return m_endpointsElectrons[i].driftLine.size() + 2;
320}

◆ 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 149 of file AvalancheMicroscopic.hh.

149 {
150 return m_endpointsElectrons.size();
151 }

Referenced by main().

◆ GetNumberOfHoleDriftLinePoints()

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

Definition at line 322 of file AvalancheMicroscopic.cc.

323 {
324 if (i >= m_endpointsHoles.size()) {
325 std::cerr << m_className << "::GetNumberOfHoleDriftLinePoints:\n";
326 std::cerr << " Endpoint index (" << i << ") out of range.\n";
327 return 0;
328 }
329
330 if (!m_useDriftLines) return 2;
331
332 return m_endpointsHoles[i].driftLine.size() + 2;
333}

◆ GetNumberOfHoleEndpoints()

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

Definition at line 178 of file AvalancheMicroscopic.hh.

178 {
179 return m_endpointsHoles.size();
180 }

◆ GetNumberOfPhotons()

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

Definition at line 185 of file AvalancheMicroscopic.hh.

185{ 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 400 of file AvalancheMicroscopic.cc.

404 {
405 if (i >= m_photons.size()) {
406 std::cerr << m_className << "::GetPhoton: Index out of range.\n";
407 return;
408 }
409
410 x0 = m_photons[i].x0;
411 x1 = m_photons[i].x1;
412 y0 = m_photons[i].y0;
413 y1 = m_photons[i].y1;
414 z0 = m_photons[i].z0;
415 z1 = m_photons[i].z1;
416 t0 = m_photons[i].t0;
417 t1 = m_photons[i].t1;
418 status = m_photons[i].status;
419 e = m_photons[i].energy;
420}

◆ GetPhotonTransportCut()

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

Retrieve the energy threshold for transporting photons.

Definition at line 115 of file AvalancheMicroscopic.hh.

115{ 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 129 of file AvalancheMicroscopic.hh.

129{ m_nCollSkip = n; }

Referenced by main().

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

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

◆ 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 108 of file AvalancheMicroscopic.hh.

108{ m_deltaCut = cut; }

◆ SetPhotonTransportCut()

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

Set an energy threshold for photon transport.

Definition at line 113 of file AvalancheMicroscopic.hh.

113{ m_gammaCut = cut; }

◆ SetSensor()

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

Set the sensor.

Definition at line 98 of file AvalancheMicroscopic.cc.

98 {
99 if (!s) {
100 std::cerr << m_className << "::SetSensor: Null pointer.\n";
101 return;
102 }
103 m_sensor = s;
104}

Referenced by GarfieldPhysics::InitializePhysics(), and main().

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

215 {
216 if (fabs(t1 - t0) < Small) {
217 std::cerr << m_className << "::SetTimeWindow:\n";
218 std::cerr << " Time interval must be greater than zero.\n";
219 return;
220 }
221
222 m_tMin = std::min(t0, t1);
223 m_tMax = std::max(t0, t1);
224 m_hasTimeWindow = true;
225}
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 439 of file AvalancheMicroscopic.cc.

440 {
441 m_userHandleAttachment = f;
442}

◆ SetUserHandleCollision()

void Garfield::AvalancheMicroscopic::SetUserHandleCollision ( void(*)(double x, double y, double z, double t, int type, int level, Medium *m, double e0, double e1, double dx0, double dy0, double dz0, double dx1, double dy1, double dz1)  f)

Set a user handling procedure, to be called at every (real) collision.

Definition at line 432 of file AvalancheMicroscopic.cc.

435 {
436 m_userHandleCollision = f;
437}

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

445 {
446 m_userHandleInelastic = f;
447}

◆ 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 or excitation followed by Penning transfer.

Definition at line 449 of file AvalancheMicroscopic.cc.

450 {
451 m_userHandleIonisation = f;
452}

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

424 {
425 if (!f) {
426 std::cerr << m_className << "::SetUserHandleStep: Null pointer.\n";
427 return;
428 }
429 m_userHandleStep = f;
430}

◆ UnsetTimeWindow()

void Garfield::AvalancheMicroscopic::UnsetTimeWindow ( )
inline

Do not restrict the time interval within which carriers are simulated.

Definition at line 134 of file AvalancheMicroscopic.hh.

134{ m_hasTimeWindow = false; }

◆ UnsetUserHandleAttachment()

void Garfield::AvalancheMicroscopic::UnsetUserHandleAttachment ( )
inline

Deactivate the user handle called at every attachment.

Definition at line 226 of file AvalancheMicroscopic.hh.

226{ m_userHandleAttachment = nullptr; }

◆ UnsetUserHandleCollision()

void Garfield::AvalancheMicroscopic::UnsetUserHandleCollision ( )
inline

Deactivate the user handle called at every collision.

Definition at line 221 of file AvalancheMicroscopic.hh.

221{ m_userHandleCollision = nullptr; }

◆ UnsetUserHandleInelastic()

void Garfield::AvalancheMicroscopic::UnsetUserHandleInelastic ( )
inline

Deactivate the user handle called at every inelastic collision.

Definition at line 231 of file AvalancheMicroscopic.hh.

231{ m_userHandleInelastic = nullptr; }

◆ UnsetUserHandleIonisation()

void Garfield::AvalancheMicroscopic::UnsetUserHandleIonisation ( )
inline

Deactivate the user handle called at every ionisation.

Definition at line 237 of file AvalancheMicroscopic.hh.

237{ m_userHandleIonisation = nullptr; }

◆ UnsetUserHandleStep()

void Garfield::AvalancheMicroscopic::UnsetUserHandleStep ( )
inline

Deactivate the user handle called at every step.

Definition at line 213 of file AvalancheMicroscopic.hh.

213{ m_userHandleStep = nullptr; }

◆ UseWeightingPotential()

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

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

Definition at line 48 of file AvalancheMicroscopic.hh.

48 {
49 m_useWeightingPotential = on;
50 }

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