Garfield++ 5.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>

Classes

struct  Electron
 
struct  Point
 

Public Member Functions

 AvalancheMicroscopic ()
 Default constructor.
 
 AvalancheMicroscopic (Sensor *sensor)
 Constructor.
 
 ~AvalancheMicroscopic ()
 Destructor.
 
void SetSensor (Sensor *sensor)
 Set the sensor.
 
void EnablePlotting (ViewDrift *view, const size_t nColl=100)
 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 calculation of induced currents on or off (default: enabled).
 
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 EnablePathLengthComputation (const bool on=true)
 Compute and store the path length of each trajectory (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)
 Switch on/off using the magnetic field in the stepping algorithm.
 
void SetCollisionSteps (const unsigned int n)
 Set number of collisions to be skipped for storing drift lines.
 
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
 
const std::vector< Electron > & GetElectrons () const
 
const std::vector< Electron > & GetHoles () const
 
size_t GetNumberOfElectronEndpoints () const
 
void GetElectronEndpoint (const size_t 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 size_t 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
 
size_t GetNumberOfElectronDriftLinePoints (const size_t i=0) const
 
size_t GetNumberOfHoleDriftLinePoints (const size_t i=0) const
 
void GetElectronDriftLinePoint (double &x, double &y, double &z, double &t, const size_t ip, const size_t ie=0) const
 
void GetHoleDriftLinePoint (double &x, double &y, double &z, double &t, const size_t ip, const size_t ih=0) const
 
size_t GetNumberOfHoleEndpoints () const
 
void GetHoleEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, int &status) const
 
size_t GetNumberOfPhotons () const
 
void GetPhoton (const size_t 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 x, const double y, const double z, const double t, const double e, const double dx=0., const double dy=0., const double dz=0.)
 
bool AvalancheElectron (const double x, const double y, const double z, const double t, const double e, const double dx=0., const double dy=0., const double dz=0.)
 Calculate an avalanche initiated by a given electron.
 
void AddElectron (const double x, const double y, const double z, const double t, const double e, const double dx=0., const double dy=0., const double dz=0.)
 Add an electron to the list of particles to be transported.
 
bool ResumeAvalanche ()
 Continue the avalanche simulation from the current set of electrons.
 
void SetUserHandleStep (void(*f)(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole))
 Set a callback function to be 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 callback function 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() [1/2]

Garfield::AvalancheMicroscopic::AvalancheMicroscopic ( )
inline

Default constructor.

Definition at line 20 of file AvalancheMicroscopic.hh.

20: AvalancheMicroscopic(nullptr) {}
AvalancheMicroscopic()
Default constructor.

Referenced by AvalancheMicroscopic().

◆ AvalancheMicroscopic() [2/2]

Garfield::AvalancheMicroscopic::AvalancheMicroscopic ( Sensor * sensor)

Constructor.

Definition at line 125 of file AvalancheMicroscopic.cc.

125 :
126 m_sensor(sensor) {
127 m_electrons.reserve(10000);
128 m_holes.reserve(10000);
129 m_photons.reserve(1000);
130}

◆ ~AvalancheMicroscopic()

Garfield::AvalancheMicroscopic::~AvalancheMicroscopic ( )
inline

Destructor.

Definition at line 24 of file AvalancheMicroscopic.hh.

24{}

Member Function Documentation

◆ AddElectron()

void Garfield::AvalancheMicroscopic::AddElectron ( const double x,
const double y,
const double z,
const double t,
const double e,
const double dx = 0.,
const double dy = 0.,
const double dz = 0. )

Add an electron to the list of particles to be transported.

Definition at line 441 of file AvalancheMicroscopic.cc.

443 {
444
445 Electron electron;
446 electron.status = StatusAlive;
447 electron.path.emplace_back(MakePoint(x, y, z, t, e, dx, dy, dz, 0));
448 m_electrons.push_back(std::move(electron));
449}

◆ AvalancheElectron()

bool Garfield::AvalancheMicroscopic::AvalancheElectron ( const double x,
const double y,
const double z,
const double t,
const double e,
const double dx = 0.,
const double dy = 0.,
const double dz = 0. )

Calculate an avalanche initiated by a given electron.

Definition at line 431 of file AvalancheMicroscopic.cc.

433 {
434
435 std::vector<std::pair<Point, bool> > particles;
436 Point p = MakePoint(x, y, z, t, e, dx, dy, dz, 0);
437 particles.emplace_back(std::make_pair(std::move(p), false));
438 return TransportElectrons(particles, true);
439}

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

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMicroscopic::DisableAvalancheSizeLimit ( )
inline

Do not apply a limit on the avalanche size.

Definition at line 119 of file AvalancheMicroscopic.hh.

119{ m_sizeCut = 0; }

◆ DisableDebugging()

void Garfield::AvalancheMicroscopic::DisableDebugging ( )
inline

Definition at line 261 of file AvalancheMicroscopic.hh.

261{ m_debug = false; }

◆ DisableDistanceHistogramming() [1/2]

void Garfield::AvalancheMicroscopic::DisableDistanceHistogramming ( )

Stop filling distance distribution histograms.

Definition at line 230 of file AvalancheMicroscopic.cc.

230 {
231 m_histDistance = nullptr;
232 m_distanceHistogramType.clear();
233}

◆ DisableDistanceHistogramming() [2/2]

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

Stop filling distance distribution histograms for a given collision type.

Definition at line 216 of file AvalancheMicroscopic.cc.

216 {
217 if (std::find(m_distanceHistogramType.begin(), m_distanceHistogramType.end(),
218 type) == m_distanceHistogramType.end()) {
219 std::cerr << m_className << "::DisableDistanceHistogramming:\n"
220 << " Collision type " << type << " is not histogrammed.\n";
221 return;
222 }
223
224 m_distanceHistogramType.erase(
225 std::remove(m_distanceHistogramType.begin(),
226 m_distanceHistogramType.end(), type),
227 m_distanceHistogramType.end());
228}

◆ DisableElectronEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableElectronEnergyHistogramming ( )
inline

Stop histogramming the electron energy distribution.

Definition at line 63 of file AvalancheMicroscopic.hh.

63{ m_histElectronEnergy = nullptr; }

◆ DisableHoleEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableHoleEnergyHistogramming ( )
inline

Stop histogramming the hole energy distribution.

Definition at line 67 of file AvalancheMicroscopic.hh.

67{ m_histHoleEnergy = nullptr; }

◆ DisablePlotting()

void Garfield::AvalancheMicroscopic::DisablePlotting ( )
inline

Switch off drift line plotting.

Definition at line 32 of file AvalancheMicroscopic.hh.

32{ m_viewer = nullptr; }

◆ DisableSecondaryEnergyHistogramming()

void Garfield::AvalancheMicroscopic::DisableSecondaryEnergyHistogramming ( )
inline

Stop histogramming the secondary electron energy distribution.

Definition at line 85 of file AvalancheMicroscopic.hh.

85{ m_histSecondary = nullptr; }

◆ DriftElectron()

bool Garfield::AvalancheMicroscopic::DriftElectron ( const double x,
const double y,
const double z,
const double t,
const double e,
const double dx = 0.,
const double dy = 0.,
const double dz = 0. )

Calculate an electron drift line.

Parameters
x,y,z,tstarting point of the electron
einitial energy of the electron
dx,dy,dzinitial direction vector of the electron If the initial direction is not specified, it is sampled randomly. Secondary electrons are not transported.

Definition at line 422 of file AvalancheMicroscopic.cc.

424 {
425 std::vector<std::pair<Point, bool> > particles;
426 Point p = MakePoint(x, y, z, t, e, dx, dy, dz, 0);
427 particles.emplace_back(std::make_pair(std::move(p), false));
428 return TransportElectrons(particles, false);
429}

◆ EnableAttachmentMarkers()

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

Draw a marker at every attachment or not.

Definition at line 38 of file AvalancheMicroscopic.hh.

38{ 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 117 of file AvalancheMicroscopic.hh.

117{ 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 95 of file AvalancheMicroscopic.hh.

95 {
96 m_useBandStructure = on;
97 }

◆ EnableDebugging()

void Garfield::AvalancheMicroscopic::EnableDebugging ( )
inline

Switch on debugging messages.

Definition at line 260 of file AvalancheMicroscopic.hh.

260{ m_debug = true; }

◆ EnableDistanceHistogramming()

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

Fill distance distribution histograms for a given collision type.

Definition at line 194 of file AvalancheMicroscopic.cc.

194 {
195 // Check if this type of collision is already registered
196 // for histogramming.
197 const unsigned int nDistanceHistogramTypes = m_distanceHistogramType.size();
198 if (nDistanceHistogramTypes > 0) {
199 for (unsigned int i = 0; i < nDistanceHistogramTypes; ++i) {
200 if (m_distanceHistogramType[i] != type) continue;
201 std::cout << m_className << "::EnableDistanceHistogramming:\n";
202 std::cout << " Collision type " << type
203 << " is already being histogrammed.\n";
204 return;
205 }
206 }
207
208 m_distanceHistogramType.push_back(type);
209 std::cout << m_className << "::EnableDistanceHistogramming:\n";
210 std::cout << " Histogramming of collision type " << type << " enabled.\n";
211 if (!m_histDistance) {
212 std::cout << " Don't forget to set the histogram.\n";
213 }
214}

◆ EnableDriftLines()

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

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

Definition at line 88 of file AvalancheMicroscopic.hh.

88{ m_storeDriftLines = on; }

◆ EnableElectronEnergyHistogramming()

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

Fill a histogram with the electron energy distribution.

Definition at line 150 of file AvalancheMicroscopic.cc.

150 {
151 if (!histo) {
152 std::cerr << m_className << "::EnableElectronEnergyHistogramming:\n"
153 << " Null pointer.\n";
154 return;
155 }
156
157 m_histElectronEnergy = histo;
158}

◆ EnableExcitationMarkers()

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

Draw a marker at every excitation or not.

Definition at line 34 of file AvalancheMicroscopic.hh.

34{ m_plotExcitations = on; }

Referenced by main().

◆ EnableHoleEnergyHistogramming()

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

Fill a histogram with the hole energy distribution.

Definition at line 160 of file AvalancheMicroscopic.cc.

160 {
161 if (!histo) {
162 std::cerr << m_className << "::EnableHoleEnergyHistogramming:\n"
163 << " Null pointer.\n";
164 return;
165 }
166
167 m_histHoleEnergy = histo;
168}

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

Referenced by main().

◆ EnableMagneticField()

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

Switch on/off using the magnetic field in the stepping algorithm.

Definition at line 124 of file AvalancheMicroscopic.hh.

124 {
125 m_useBfieldAuto = false;
126 m_useBfield = on;
127 }

◆ EnableNullCollisionSteps()

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

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

Definition at line 100 of file AvalancheMicroscopic.hh.

100 {
101 m_useNullCollisionSteps = on;
102 }

◆ EnablePathLengthComputation()

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

Compute and store the path length of each trajectory (default: off).

Definition at line 57 of file AvalancheMicroscopic.hh.

57 {
58 m_computePathLength = on;
59 }

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

92{ m_usePhotons = on; }

◆ EnablePlotting()

void Garfield::AvalancheMicroscopic::EnablePlotting ( ViewDrift * view,
const size_t nColl = 100 )

Switch on drift line plotting.

Definition at line 140 of file AvalancheMicroscopic.cc.

141 {
142 if (!view) {
143 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
144 return;
145 }
146 m_viewer = view;
147 m_nCollPlot = std::max(nColl, 1ul);
148}

Referenced by main().

◆ EnableSecondaryEnergyHistogramming()

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

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

Definition at line 235 of file AvalancheMicroscopic.cc.

235 {
236 if (!histo) {
237 std::cerr << m_className << "::EnableSecondaryEnergyHistogramming:\n"
238 << " Null pointer.\n";
239 return;
240 }
241
242 m_histSecondary = histo;
243}

◆ EnableSignalCalculation()

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

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

Definition at line 41 of file AvalancheMicroscopic.hh.

41{ 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 49 of file AvalancheMicroscopic.hh.

49 {
50 m_integrateWeightingField = on;
51 }

◆ GetAvalancheSize() [1/2]

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

Definition at line 142 of file AvalancheMicroscopic.hh.

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

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

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

Referenced by GarfieldPhysics::DoIt().

◆ GetAvalancheSizeLimit()

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

Retrieve the currently set size limit.

Definition at line 121 of file AvalancheMicroscopic.hh.

121{ return m_sizeCut; }

◆ GetElectronDriftLinePoint()

void Garfield::AvalancheMicroscopic::GetElectronDriftLinePoint ( double & x,
double & y,
double & z,
double & t,
const size_t ip,
const size_t ie = 0 ) const

Definition at line 331 of file AvalancheMicroscopic.cc.

333 {
334 if (ie >= m_electrons.size()) {
335 std::cerr << m_className << "::GetElectronDriftLinePoint:\n"
336 << " Endpoint index (" << ie << ") out of range.\n";
337 return;
338 }
339 if (ip >= m_electrons[ie].path.size()) {
340 std::cerr << m_className << "::GetElectronDriftLinePoint:\n"
341 << " Drift line point index (" << ip << ") out of range.\n";
342 return;
343 }
344 x = m_electrons[ie].path[ip].x;
345 y = m_electrons[ie].path[ip].y;
346 z = m_electrons[ie].path[ip].z;
347 t = m_electrons[ie].path[ip].t;
348}

◆ GetElectronEndpoint() [1/2]

void Garfield::AvalancheMicroscopic::GetElectronEndpoint ( const size_t 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

◆ GetElectronEndpoint() [2/2]

void Garfield::AvalancheMicroscopic::GetElectronEndpoint ( const size_t 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 257 of file AvalancheMicroscopic.cc.

260 {
261 if (i >= m_electrons.size()) {
262 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
263 status = -3;
264 return;
265 }
266 if (m_electrons[i].path.empty()) {
267 std::cerr << m_className << "::GetElectronEndpoint: Empty drift line.\n";
268 status = -3;
269 return;
270 }
271 x0 = m_electrons[i].path[0].x;
272 y0 = m_electrons[i].path[0].y;
273 z0 = m_electrons[i].path[0].z;
274 t0 = m_electrons[i].path[0].t;
275 e0 = m_electrons[i].path[0].energy;
276 x1 = m_electrons[i].path.back().x;
277 y1 = m_electrons[i].path.back().y;
278 z1 = m_electrons[i].path.back().z;
279 t1 = m_electrons[i].path.back().t;
280 e1 = m_electrons[i].path.back().energy;
281 status = m_electrons[i].status;
282}

◆ GetElectrons()

const std::vector< Electron > & Garfield::AvalancheMicroscopic::GetElectrons ( ) const
inline

Definition at line 162 of file AvalancheMicroscopic.hh.

162{ return m_electrons; }

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

◆ GetElectronTransportCut()

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

Retrieve the value of the energy threshold.

Definition at line 108 of file AvalancheMicroscopic.hh.

108{ return m_deltaCut; }

◆ GetHoleDriftLinePoint()

void Garfield::AvalancheMicroscopic::GetHoleDriftLinePoint ( double & x,
double & y,
double & z,
double & t,
const size_t ip,
const size_t ih = 0 ) const

Definition at line 350 of file AvalancheMicroscopic.cc.

353 {
354 if (ih >= m_holes.size()) {
355 std::cerr << m_className << "::GetHoleDriftLinePoint:\n"
356 << " Endpoint index (" << ih << ") out of range.\n";
357 return;
358 }
359 if (ip >= m_holes[ih].path.size()) {
360 std::cerr << m_className << "::GetHoleDriftLinePoint:\n"
361 << " Drift line point index (" << ip << ") out of range.\n";
362 return;
363 }
364 x = m_holes[ih].path[ip].x;
365 y = m_holes[ih].path[ip].y;
366 z = m_holes[ih].path[ip].z;
367 t = m_holes[ih].path[ip].t;
368}

◆ GetHoleEndpoint()

void Garfield::AvalancheMicroscopic::GetHoleEndpoint ( const size_t 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 284 of file AvalancheMicroscopic.cc.

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

◆ GetHoles()

const std::vector< Electron > & Garfield::AvalancheMicroscopic::GetHoles ( ) const
inline

Definition at line 163 of file AvalancheMicroscopic.hh.

163{ return m_holes; }

◆ GetNumberOfElectronDriftLinePoints()

size_t Garfield::AvalancheMicroscopic::GetNumberOfElectronDriftLinePoints ( const size_t i = 0) const

Definition at line 311 of file AvalancheMicroscopic.cc.

312 {
313 if (i >= m_electrons.size()) {
314 std::cerr << m_className << "::GetNumberOfElectronDriftLinePoints: "
315 << "Index out of range.\n";
316 return 0;
317 }
318 return m_electrons[i].path.size();
319}

◆ GetNumberOfElectronEndpoints()

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

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

Definition at line 166 of file AvalancheMicroscopic.hh.

166 {
167 return m_electrons.size();
168 }

Referenced by main().

◆ GetNumberOfHoleDriftLinePoints()

size_t Garfield::AvalancheMicroscopic::GetNumberOfHoleDriftLinePoints ( const size_t i = 0) const

Definition at line 321 of file AvalancheMicroscopic.cc.

322 {
323 if (i >= m_holes.size()) {
324 std::cerr << m_className << "::GetNumberOfHoleDriftLinePoints: "
325 << "Index out of range.\n";
326 return 0;
327 }
328 return m_holes[i].path.size();
329}

◆ GetNumberOfHoleEndpoints()

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

Definition at line 193 of file AvalancheMicroscopic.hh.

193{ return m_holes.size(); }

◆ GetNumberOfPhotons()

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

Definition at line 198 of file AvalancheMicroscopic.hh.

198{ return m_photons.size(); }

◆ GetPhoton()

void Garfield::AvalancheMicroscopic::GetPhoton ( const size_t 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 370 of file AvalancheMicroscopic.cc.

372 {
373 if (i >= m_photons.size()) {
374 std::cerr << m_className << "::GetPhoton: Index out of range.\n";
375 return;
376 }
377
378 x0 = m_photons[i].x0;
379 x1 = m_photons[i].x1;
380 y0 = m_photons[i].y0;
381 y1 = m_photons[i].y1;
382 z0 = m_photons[i].z0;
383 z1 = m_photons[i].z1;
384 t0 = m_photons[i].t0;
385 t1 = m_photons[i].t1;
386 status = m_photons[i].status;
387 e = m_photons[i].energy;
388}

◆ GetPhotonTransportCut()

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

Retrieve the energy threshold for transporting photons.

Definition at line 113 of file AvalancheMicroscopic.hh.

113{ return m_gammaCut; }

◆ ResumeAvalanche()

bool Garfield::AvalancheMicroscopic::ResumeAvalanche ( )

Continue the avalanche simulation from the current set of electrons.

Definition at line 451 of file AvalancheMicroscopic.cc.

451 {
452 std::vector<std::pair<Point, bool> > particles;
453 for (const auto& p : m_electrons) {
454 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
455 particles.emplace_back(std::make_pair(p.path.back(), false));
456 }
457 }
458 for (const auto& p : m_holes) {
459 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
460 particles.emplace_back(std::make_pair(p.path.back(), true));
461 }
462 }
463 return TransportElectrons(particles, true);
464}

◆ SetCollisionSteps()

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

Set number of collisions to be skipped for storing drift lines.

Definition at line 130 of file AvalancheMicroscopic.hh.

130{ m_nCollSkip = n; }

◆ SetDistanceHistogram()

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

Fill histograms of the distance between successive collisions.

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

Definition at line 170 of file AvalancheMicroscopic.cc.

170 {
171 if (!histo) {
172 std::cerr << m_className << "::SetDistanceHistogram: Null pointer.\n";
173 return;
174 }
175
176 m_histDistance = histo;
177
178 if (opt == 'x' || opt == 'y' || opt == 'z' || opt == 'r') {
179 m_distanceOption = opt;
180 } else {
181 std::cerr << m_className << "::SetDistanceHistogram:";
182 std::cerr << " Unknown option " << opt << ".\n";
183 std::cerr << " Valid options are x, y, z, r.\n";
184 std::cerr << " Using default value (r).\n";
185 m_distanceOption = 'r';
186 }
187
188 if (m_distanceHistogramType.empty()) {
189 std::cout << m_className << "::SetDistanceHistogram:\n";
190 std::cout << " Don't forget to call EnableDistanceHistogramming.\n";
191 }
192}

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

106{ m_deltaCut = cut; }

◆ SetPhotonTransportCut()

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

Set an energy threshold for photon transport.

Definition at line 111 of file AvalancheMicroscopic.hh.

111{ m_gammaCut = cut; }

◆ SetSensor()

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

Set the sensor.

Definition at line 132 of file AvalancheMicroscopic.cc.

132 {
133 if (!s) {
134 std::cerr << m_className << "::SetSensor: Null pointer.\n";
135 return;
136 }
137 m_sensor = s;
138}

Referenced by GarfieldPhysics::DoIt().

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

245 {
246 if (fabs(t1 - t0) < Small) {
247 std::cerr << m_className << "::SetTimeWindow:\n";
248 std::cerr << " Time interval must be greater than zero.\n";
249 return;
250 }
251
252 m_tMin = std::min(t0, t1);
253 m_tMax = std::max(t0, t1);
254 m_hasTimeWindow = true;
255}
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))

Set a user handling procedure, to be called at every attachment.

Definition at line 407 of file AvalancheMicroscopic.cc.

408 {
409 m_userHandleAttachment = f;
410}

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

Set a callback function to be called at every (real) collision.

Definition at line 400 of file AvalancheMicroscopic.cc.

403 {
404 m_userHandleCollision = f;
405}

◆ SetUserHandleInelastic()

void Garfield::AvalancheMicroscopic::SetUserHandleInelastic ( void(* )(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.

Definition at line 412 of file AvalancheMicroscopic.cc.

413 {
414 m_userHandleInelastic = f;
415}

◆ SetUserHandleIonisation()

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

Set a user handling procedure, to be called at every ionising collision or excitation followed by Penning transfer.

Definition at line 417 of file AvalancheMicroscopic.cc.

418 {
419 m_userHandleIonisation = f;
420}

◆ SetUserHandleStep()

void Garfield::AvalancheMicroscopic::SetUserHandleStep ( void(* )(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole))

Set a callback function to be called at every step.

Definition at line 390 of file AvalancheMicroscopic.cc.

392 {
393 if (!f) {
394 std::cerr << m_className << "::SetUserHandleStep: Null pointer.\n";
395 return;
396 }
397 m_userHandleStep = f;
398}

◆ UnsetTimeWindow()

void Garfield::AvalancheMicroscopic::UnsetTimeWindow ( )
inline

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

Definition at line 135 of file AvalancheMicroscopic.hh.

135{ m_hasTimeWindow = false; }

◆ UnsetUserHandleAttachment()

void Garfield::AvalancheMicroscopic::UnsetUserHandleAttachment ( )
inline

Deactivate the user handle called at every attachment.

Definition at line 246 of file AvalancheMicroscopic.hh.

246{ m_userHandleAttachment = nullptr; }

◆ UnsetUserHandleCollision()

void Garfield::AvalancheMicroscopic::UnsetUserHandleCollision ( )
inline

Deactivate the user handle called at every collision.

Definition at line 241 of file AvalancheMicroscopic.hh.

241{ m_userHandleCollision = nullptr; }

◆ UnsetUserHandleInelastic()

void Garfield::AvalancheMicroscopic::UnsetUserHandleInelastic ( )
inline

Deactivate the user handle called at every inelastic collision.

Definition at line 251 of file AvalancheMicroscopic.hh.

251{ m_userHandleInelastic = nullptr; }

◆ UnsetUserHandleIonisation()

void Garfield::AvalancheMicroscopic::UnsetUserHandleIonisation ( )
inline

Deactivate the user handle called at every ionisation.

Definition at line 257 of file AvalancheMicroscopic.hh.

257{ m_userHandleIonisation = nullptr; }

◆ UnsetUserHandleStep()

void Garfield::AvalancheMicroscopic::UnsetUserHandleStep ( )
inline

Deactivate the user handle called at every step.

Definition at line 233 of file AvalancheMicroscopic.hh.

233{ 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 54 of file AvalancheMicroscopic.hh.

54{ 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 44 of file AvalancheMicroscopic.hh.

44 {
45 m_useWeightingPotential = on;
46 }

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