Garfield++ 4.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 UseInducedCharge (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.
 
bool ResumeAvalanche ()
 
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}

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

466 {
467 std::vector<Electron> stack;
468 AddToStack(x0, y0, z0, t0, e0, dx0, dy0, dz0, 0, false, stack);
469 return TransportElectrons(stack, true);
470}

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

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMicroscopic::DisableAvalancheSizeLimit ( )
inline

Do not apply a limit on the avalanche size.

Definition at line 113 of file AvalancheMicroscopic.hh.

113{ m_sizeCut = 0; }

◆ DisableDebugging()

void Garfield::AvalancheMicroscopic::DisableDebugging ( )
inline

Definition at line 235 of file AvalancheMicroscopic.hh.

235{ m_debug = false; }

◆ DisableDistanceHistogramming() [1/2]

void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming ( )

Stop filling distance distribution histograms.

Definition at line 199 of file AvalancheMicroscopic.cc.

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

◆ DisableDistanceHistogramming() [2/2]

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

Stop filling distance distribution histograms for a given collision type.

Definition at line 185 of file AvalancheMicroscopic.cc.

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

◆ DisableElectronEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableElectronEnergyHistogramming ( )
inline

Stop histogramming the electron energy distribution.

Definition at line 57 of file AvalancheMicroscopic.hh.

57{ m_histElectronEnergy = nullptr; }

◆ DisableHoleEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableHoleEnergyHistogramming ( )
inline

Stop histogramming the hole energy distribution.

Definition at line 61 of file AvalancheMicroscopic.hh.

61{ 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 79 of file AvalancheMicroscopic.hh.

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

456 {
457 std::vector<Electron> stack;
458 AddToStack(x0, y0, z0, t0, e0, dx0, dy0, dz0, 0, false, stack);
459 return TransportElectrons(stack, false);
460}

◆ EnableAttachmentMarkers()

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

Draw a marker at every attachment or not.

Definition at line 36 of file AvalancheMicroscopic.hh.

36{ m_plotAttachments = on; }

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

111{ 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 89 of file AvalancheMicroscopic.hh.

89 {
90 m_useBandStructure = on;
91 }

◆ EnableDebugging()

void Garfield::AvalancheMicroscopic::EnableDebugging ( )
inline

Switch on debugging messages.

Definition at line 234 of file AvalancheMicroscopic.hh.

234{ m_debug = true; }

◆ EnableDistanceHistogramming()

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

Fill distance distribution histograms for a given collision type.

Definition at line 163 of file AvalancheMicroscopic.cc.

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

◆ EnableDriftLines()

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

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

Definition at line 82 of file AvalancheMicroscopic.hh.

82{ m_storeDriftLines = on; }

Referenced by EnablePlotting().

◆ EnableElectronEnergyHistogramming()

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

Fill a histogram with the electron energy distribution.

Definition at line 119 of file AvalancheMicroscopic.cc.

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

◆ 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{ m_plotExcitations = on; }

Referenced by main().

◆ EnableHoleEnergyHistogramming()

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

Fill a histogram with the hole energy distribution.

Definition at line 129 of file AvalancheMicroscopic.cc.

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

◆ EnableIonisationMarkers()

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

Draw a marker at every ionising collision or not.

Definition at line 34 of file AvalancheMicroscopic.hh.

34{ m_plotIonisations = on; }

Referenced by main().

◆ EnableMagneticField()

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

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

Definition at line 118 of file AvalancheMicroscopic.hh.

118{ m_useBfield = on; }

Referenced by main().

◆ EnableNullCollisionSteps()

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

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

Definition at line 94 of file AvalancheMicroscopic.hh.

94 {
95 m_useNullCollisionSteps = on;
96 }

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

86{ m_usePhotons = on; }

◆ EnablePlotting()

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

Switch on drift line plotting.

Definition at line 105 of file AvalancheMicroscopic.cc.

105 {
106 if (!view) {
107 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
108 return;
109 }
110
111 m_viewer = view;
112 if (!m_storeDriftLines) {
113 std::cout << m_className << "::EnablePlotting:\n"
114 << " Enabling storage of drift line.\n";
116 }
117}
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 204 of file AvalancheMicroscopic.cc.

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

◆ EnableSignalCalculation()

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

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

Definition at line 39 of file AvalancheMicroscopic.hh.

39{ m_doSignal = on; }

Referenced by main().

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

47 {
48 m_integrateWeightingField = on;
49 }

◆ GetAvalancheSize() [1/2]

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

Definition at line 133 of file AvalancheMicroscopic.hh.

133 {
134 ne = m_nElectrons;
135 nh = m_nHoles;
136 ni = m_nIons;
137 }

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

129 {
130 ne = m_nElectrons;
131 ni = m_nIons;
132 }

Referenced by GarfieldPhysics::DoIt().

◆ GetAvalancheSizeLimit()

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

Retrieve the currently set size limit.

Definition at line 115 of file AvalancheMicroscopic.hh.

115{ 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 334 of file AvalancheMicroscopic.cc.

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

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

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

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

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

Referenced by Garfield::AvalancheGrid::GetElectronsFromAvalancheMicroscopic().

◆ GetElectronTransportCut()

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

Retrieve the value of the energy threshold.

Definition at line 102 of file AvalancheMicroscopic.hh.

102{ 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 366 of file AvalancheMicroscopic.cc.

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

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

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

◆ GetNumberOfElectronDriftLinePoints()

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

Definition at line 308 of file AvalancheMicroscopic.cc.

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

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

141 {
142 return m_endpointsElectrons.size();
143 }

Referenced by Garfield::AvalancheGrid::GetElectronsFromAvalancheMicroscopic(), and main().

◆ GetNumberOfHoleDriftLinePoints()

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

Definition at line 321 of file AvalancheMicroscopic.cc.

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

◆ GetNumberOfHoleEndpoints()

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

Definition at line 170 of file AvalancheMicroscopic.hh.

170 {
171 return m_endpointsHoles.size();
172 }

◆ GetNumberOfPhotons()

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

Definition at line 177 of file AvalancheMicroscopic.hh.

177{ 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 399 of file AvalancheMicroscopic.cc.

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

◆ GetPhotonTransportCut()

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

Retrieve the energy threshold for transporting photons.

Definition at line 107 of file AvalancheMicroscopic.hh.

107{ return m_gammaCut; }

◆ ResumeAvalanche()

bool Garfield::AvalancheMicroscopic::ResumeAvalanche ( )

Definition at line 472 of file AvalancheMicroscopic.cc.

472 {
473 std::vector<Electron> stack;
474 for (const auto& p : m_endpointsElectrons) {
475 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
476 AddToStack(p.x, p.y, p.z, p.t, p.energy,
477 p.kx, p.ky, p.kz, p.band, false, stack);
478 }
479 }
480 for (const auto& p : m_endpointsHoles) {
481 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
482 AddToStack(p.x, p.y, p.z, p.t, p.energy,
483 p.kx, p.ky, p.kz, p.band, true, stack);
484 }
485 }
486 return TransportElectrons(stack, true);
487}

◆ SetCollisionSteps()

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

Set number of collisions to be skipped for plotting.

Definition at line 121 of file AvalancheMicroscopic.hh.

121{ 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 139 of file AvalancheMicroscopic.cc.

139 {
140 if (!histo) {
141 std::cerr << m_className << "::SetDistanceHistogram: 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 100 of file AvalancheMicroscopic.hh.

100{ m_deltaCut = cut; }

◆ SetPhotonTransportCut()

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

Set an energy threshold for photon transport.

Definition at line 105 of file AvalancheMicroscopic.hh.

105{ m_gammaCut = cut; }

◆ SetSensor()

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

Set the sensor.

Definition at line 97 of file AvalancheMicroscopic.cc.

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

Referenced by GarfieldPhysics::DoIt(), 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 214 of file AvalancheMicroscopic.cc.

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

439 {
440 m_userHandleAttachment = f;
441}

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

434 {
435 m_userHandleCollision = f;
436}

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

444 {
445 m_userHandleInelastic = f;
446}

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

449 {
450 m_userHandleIonisation = f;
451}

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

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

◆ UnsetTimeWindow()

void Garfield::AvalancheMicroscopic::UnsetTimeWindow ( )
inline

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

Definition at line 126 of file AvalancheMicroscopic.hh.

126{ m_hasTimeWindow = false; }

◆ UnsetUserHandleAttachment()

void Garfield::AvalancheMicroscopic::UnsetUserHandleAttachment ( )
inline

Deactivate the user handle called at every attachment.

Definition at line 220 of file AvalancheMicroscopic.hh.

220{ m_userHandleAttachment = nullptr; }

◆ UnsetUserHandleCollision()

void Garfield::AvalancheMicroscopic::UnsetUserHandleCollision ( )
inline

Deactivate the user handle called at every collision.

Definition at line 215 of file AvalancheMicroscopic.hh.

215{ m_userHandleCollision = nullptr; }

◆ UnsetUserHandleInelastic()

void Garfield::AvalancheMicroscopic::UnsetUserHandleInelastic ( )
inline

Deactivate the user handle called at every inelastic collision.

Definition at line 225 of file AvalancheMicroscopic.hh.

225{ m_userHandleInelastic = nullptr; }

◆ UnsetUserHandleIonisation()

void Garfield::AvalancheMicroscopic::UnsetUserHandleIonisation ( )
inline

Deactivate the user handle called at every ionisation.

Definition at line 231 of file AvalancheMicroscopic.hh.

231{ m_userHandleIonisation = nullptr; }

◆ UnsetUserHandleStep()

void Garfield::AvalancheMicroscopic::UnsetUserHandleStep ( )
inline

Deactivate the user handle called at every step.

Definition at line 207 of file AvalancheMicroscopic.hh.

207{ m_userHandleStep = nullptr; }

◆ UseInducedCharge()

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

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

Definition at line 52 of file AvalancheMicroscopic.hh.

52{ m_doInducedCharge = on; }

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

42 {
43 m_useWeightingPotential = on;
44 }

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