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::AvalancheMC Class Reference

#include <AvalancheMC.hh>

Public Member Functions

 AvalancheMC ()
 Constructor.
 
 ~AvalancheMC ()
 Destructor.
 
void SetSensor (Sensor *s)
 Set the sensor.
 
void EnablePlotting (ViewDrift *view)
 Switch on drift line plotting.
 
void DisablePlotting ()
 Switch off drift line plotting.
 
void EnableSignalCalculation (const bool on=true)
 Switch on calculation of induced currents (default: disabled).
 
void SetSignalAveragingOrder (const unsigned int navg)
 
void UseWeightingPotential (const bool on=true)
 
void EnableInducedChargeCalculation (const bool on=true)
 Switch on calculation of induced charge (default: disabled).
 
void EnableRKFSteps (const bool on=true)
 
void EnableProjectedPathIntegration (const bool on=true)
 
void EnableDiffusion ()
 Switch on diffusion (default: enabled).
 
void DisableDiffusion ()
 Switch off diffusion.
 
void EnableAttachment ()
 
void DisableAttachment ()
 Switch off attachment and multiplication.
 
void EnableAttachmentMap (const bool on=true)
 Retrieve the attachment coefficient from the component.
 
void EnableVelocityMap (const bool on=true)
 Retrieve the drift velocity from the component.
 
void EnableMagneticField (const bool on=true)
 Enable use of magnetic field in stepping algorithm.
 
void EnableAvalancheSizeLimit (const unsigned int size)
 
void DisableAvalancheSizeLimit ()
 Do not limit the maximum avalanche size.
 
int GetAvalancheSizeLimit () const
 Return the currently set avalanche size limit.
 
void SetTimeSteps (const double d=0.02)
 Use fixed-time steps (default 20 ps).
 
void SetDistanceSteps (const double d=0.001)
 Use fixed distance steps (default 10 um).
 
void SetCollisionSteps (const unsigned int n=100)
 
void SetStepDistanceFunction (double(*f)(double x, double y, double z))
 Retrieve the step distance from a user-supplied function.
 
void SetTimeWindow (const double t0, const double t1)
 Define a time interval (only carriers inside the interval are drifted).
 
void UnsetTimeWindow ()
 Do not limit the time interval within which carriers are drifted.
 
void SetElectronSignalScalingFactor (const double scale)
 Set multiplication factor for the signal induced by electrons.
 
void SetHoleSignalScalingFactor (const double scale)
 Set multiplication factor for the signal induced by holes.
 
void SetIonSignalScalingFactor (const double scale)
 Set multiplication factor for the signal induced by ions.
 
void GetAvalancheSize (unsigned int &ne, unsigned int &ni) const
 Return the number of electrons and ions/holes in the avalanche.
 
unsigned int GetNumberOfDriftLinePoints () const
 Return the number of points along the last simulated drift line.
 
void GetDriftLinePoint (const unsigned int i, double &x, double &y, double &z, double &t) const
 Return the coordinates and time of a point along the last drift line.
 
unsigned int GetNumberOfElectronEndpoints () const
 
unsigned int GetNumberOfHoleEndpoints () const
 
unsigned int GetNumberOfIonEndpoints () const
 Return the number of ion trajectories.
 
void GetElectronEndpoint (const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
 
void GetHoleEndpoint (const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
 
void GetIonEndpoint (const unsigned int i, 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)
 Simulate the drift line of an electron from a given starting point.
 
bool DriftHole (const double x0, const double y0, const double z0, const double t0)
 Simulate the drift line of a hole from a given starting point.
 
bool DriftIon (const double x0, const double y0, const double z0, const double t0)
 Simulate the drift line of an ion from a given starting point.
 
bool AvalancheElectron (const double x0, const double y0, const double z0, const double t0, const bool hole=false)
 
bool AvalancheHole (const double x0, const double y0, const double z0, const double t0, const bool electron=false)
 Simulate an avalanche initiated by a hole at a given starting point.
 
bool AvalancheElectronHole (const double x0, const double y0, const double z0, const double t0)
 Simulate an avalanche initiated by an electron-hole pair.
 
bool ResumeAvalanche (const bool electron=true, const bool hole=true)
 
void EnableDebugging (const bool on=true)
 Switch debugging messages on/off (default: off).
 

Detailed Description

Calculate drift lines and avalanches based on macroscopic transport coefficients, using Monte Carlo integration.

Definition at line 17 of file AvalancheMC.hh.

Constructor & Destructor Documentation

◆ AvalancheMC()

Garfield::AvalancheMC::AvalancheMC ( )

Constructor.

Definition at line 27 of file AvalancheMC.cc.

27{ m_drift.reserve(10000); }

◆ ~AvalancheMC()

Garfield::AvalancheMC::~AvalancheMC ( )
inline

Destructor.

Definition at line 22 of file AvalancheMC.hh.

22{}

Member Function Documentation

◆ AvalancheElectron()

bool Garfield::AvalancheMC::AvalancheElectron ( const double  x0,
const double  y0,
const double  z0,
const double  t0,
const bool  hole = false 
)

Simulate an avalanche initiated by an electron at a given starting point.

Parameters
x0,y0,z0,t0coordinates and time of the initial electron
holesimulate the hole component of the avalanche or not

Definition at line 517 of file AvalancheMC.cc.

519 {
520 std::vector<DriftPoint> aval;
521 AddPoint({x0, y0, z0}, t0, Particle::Electron, 1, aval);
522 return Avalanche(aval, true, holes);
523}

◆ AvalancheElectronHole()

bool Garfield::AvalancheMC::AvalancheElectronHole ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Simulate an avalanche initiated by an electron-hole pair.

Definition at line 533 of file AvalancheMC.cc.

534 {
535 std::vector<DriftPoint> aval;
536 AddPoint({x0, y0, z0}, t0, Particle::Electron, 1, aval);
537 AddPoint({x0, y0, z0}, t0, Particle::Hole, 1, aval);
538 return Avalanche(aval, true, true);
539}

◆ AvalancheHole()

bool Garfield::AvalancheMC::AvalancheHole ( const double  x0,
const double  y0,
const double  z0,
const double  t0,
const bool  electron = false 
)

Simulate an avalanche initiated by a hole at a given starting point.

Definition at line 525 of file AvalancheMC.cc.

527 {
528 std::vector<DriftPoint> aval;
529 AddPoint({x0, y0, z0}, t0, Particle::Hole, 1, aval);
530 return Avalanche(aval, electrons, true);
531}

◆ DisableAttachment()

void Garfield::AvalancheMC::DisableAttachment ( )
inline

Switch off attachment and multiplication.

Definition at line 69 of file AvalancheMC.hh.

69{ m_useAttachment = false; }

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMC::DisableAvalancheSizeLimit ( )
inline

Do not limit the maximum avalanche size.

Definition at line 84 of file AvalancheMC.hh.

84{ m_sizeCut = 0; }

◆ DisableDiffusion()

void Garfield::AvalancheMC::DisableDiffusion ( )
inline

Switch off diffusion.

Definition at line 63 of file AvalancheMC.hh.

63{ m_useDiffusion = false; }

◆ DisablePlotting()

void Garfield::AvalancheMC::DisablePlotting ( )
inline

Switch off drift line plotting.

Definition at line 30 of file AvalancheMC.hh.

30{ m_viewer = nullptr; }

◆ DriftElectron()

bool Garfield::AvalancheMC::DriftElectron ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Simulate the drift line of an electron from a given starting point.

Definition at line 186 of file AvalancheMC.cc.

187 {
188 if (!m_sensor) {
189 std::cerr << m_className << "::DriftElectron: Sensor is not defined.\n";
190 return false;
191 }
192
193 m_endpointsElectrons.clear();
194 m_endpointsHoles.clear();
195 m_endpointsIons.clear();
196
197 m_nElectrons = 1;
198 m_nHoles = 0;
199 m_nIons = 0;
200
201 std::vector<DriftPoint> secondaries;
202 return DriftLine({x0, y0, z0}, t0, Particle::Electron, secondaries);
203}

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

◆ DriftHole()

bool Garfield::AvalancheMC::DriftHole ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Simulate the drift line of a hole from a given starting point.

Definition at line 205 of file AvalancheMC.cc.

206 {
207 if (!m_sensor) {
208 std::cerr << m_className << "::DriftHole: Sensor is not defined.\n";
209 return false;
210 }
211
212 m_endpointsElectrons.clear();
213 m_endpointsHoles.clear();
214 m_endpointsIons.clear();
215
216 m_nElectrons = 0;
217 m_nHoles = 1;
218 m_nIons = 0;
219
220 std::vector<DriftPoint> secondaries;
221 return DriftLine({x0, y0, z0}, t0, Particle::Hole, secondaries);
222}

Referenced by main().

◆ DriftIon()

bool Garfield::AvalancheMC::DriftIon ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Simulate the drift line of an ion from a given starting point.

Definition at line 224 of file AvalancheMC.cc.

225 {
226 if (!m_sensor) {
227 std::cerr << m_className << "::DriftIon: Sensor is not defined.\n";
228 return false;
229 }
230
231 m_endpointsElectrons.clear();
232 m_endpointsHoles.clear();
233 m_endpointsIons.clear();
234
235 m_nElectrons = 0;
236 m_nHoles = 0;
237 m_nIons = 1;
238
239 std::vector<DriftPoint> secondaries;
240 return DriftLine({x0, y0, z0}, t0, Particle::Ion, secondaries);
241}

◆ EnableAttachment()

void Garfield::AvalancheMC::EnableAttachment ( )
inline

Switch on attachment (and multiplication) for drift line calculation (default: enabled). For avalanches the flag is ignored.

Definition at line 67 of file AvalancheMC.hh.

67{ m_useAttachment = true; }

◆ EnableAttachmentMap()

void Garfield::AvalancheMC::EnableAttachmentMap ( const bool  on = true)
inline

Retrieve the attachment coefficient from the component.

Definition at line 72 of file AvalancheMC.hh.

72{ m_useAttachmentMap = on; }

Referenced by main().

◆ EnableAvalancheSizeLimit()

void Garfield::AvalancheMC::EnableAvalancheSizeLimit ( const unsigned int  size)
inline

Set a maximum avalanche size (ignore further multiplication once this size has been reached).

Definition at line 82 of file AvalancheMC.hh.

82{ m_sizeCut = size; }

◆ EnableDebugging()

void Garfield::AvalancheMC::EnableDebugging ( const bool  on = true)
inline

Switch debugging messages on/off (default: off).

Definition at line 178 of file AvalancheMC.hh.

178{ m_debug = on; }

◆ EnableDiffusion()

void Garfield::AvalancheMC::EnableDiffusion ( )
inline

Switch on diffusion (default: enabled).

Definition at line 61 of file AvalancheMC.hh.

61{ m_useDiffusion = true; }

◆ EnableInducedChargeCalculation()

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

Switch on calculation of induced charge (default: disabled).

Definition at line 46 of file AvalancheMC.hh.

46 {
47 m_doInducedCharge = on;
48 }

◆ EnableMagneticField()

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

Enable use of magnetic field in stepping algorithm.

Definition at line 78 of file AvalancheMC.hh.

78{ m_useBfield = on; }

◆ EnablePlotting()

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

Switch on drift line plotting.

Definition at line 38 of file AvalancheMC.cc.

38 {
39 if (!view) {
40 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
41 return;
42 }
43
44 m_viewer = view;
45}

Referenced by main().

◆ EnableProjectedPathIntegration()

void Garfield::AvalancheMC::EnableProjectedPathIntegration ( const bool  on = true)
inline

Switch on equilibration of multiplication and attachment over the drift line (default: enabled).

Definition at line 56 of file AvalancheMC.hh.

56 {
57 m_doEquilibration = on;
58 }

◆ EnableRKFSteps()

void Garfield::AvalancheMC::EnableRKFSteps ( const bool  on = true)
inline

Switch on Runge-Kutta-Fehlberg stepping (as opposed to simple straight-line steps.

Definition at line 52 of file AvalancheMC.hh.

52{ m_doRKF = on; }

◆ EnableSignalCalculation()

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

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

Definition at line 33 of file AvalancheMC.hh.

33{ m_doSignal = on; }

Referenced by main().

◆ EnableVelocityMap()

void Garfield::AvalancheMC::EnableVelocityMap ( const bool  on = true)
inline

Retrieve the drift velocity from the component.

Definition at line 75 of file AvalancheMC.hh.

75{ m_useVelocityMap = on; }

◆ GetAvalancheSize()

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

Return the number of electrons and ions/holes in the avalanche.

Definition at line 111 of file AvalancheMC.hh.

111 {
112 ne = m_nElectrons;
113 ni = std::max(m_nIons, m_nHoles);
114 }

◆ GetAvalancheSizeLimit()

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

Return the currently set avalanche size limit.

Definition at line 86 of file AvalancheMC.hh.

86{ return m_sizeCut; }

◆ GetDriftLinePoint()

void Garfield::AvalancheMC::GetDriftLinePoint ( const unsigned int  i,
double &  x,
double &  y,
double &  z,
double &  t 
) const

Return the coordinates and time of a point along the last drift line.

Definition at line 114 of file AvalancheMC.cc.

115 {
116 if (i >= m_drift.size()) {
117 std::cerr << m_className << "::GetDriftLinePoint: Index out of range.\n";
118 return;
119 }
120
121 x = m_drift[i].x[0];
122 y = m_drift[i].x[1];
123 z = m_drift[i].x[2];
124 t = m_drift[i].t;
125}

◆ GetElectronEndpoint()

void Garfield::AvalancheMC::GetElectronEndpoint ( const unsigned int  i,
double &  x0,
double &  y0,
double &  z0,
double &  t0,
double &  x1,
double &  y1,
double &  z1,
double &  t1,
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
statusstatus code (see GarfieldConstants.hh)

Definition at line 166 of file AvalancheMC.cc.

169 {
170 if (i >= m_endpointsElectrons.size()) {
171 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
172 return;
173 }
174
175 x0 = m_endpointsElectrons[i].x0[0];
176 y0 = m_endpointsElectrons[i].x0[1];
177 z0 = m_endpointsElectrons[i].x0[2];
178 t0 = m_endpointsElectrons[i].t0;
179 x1 = m_endpointsElectrons[i].x1[0];
180 y1 = m_endpointsElectrons[i].x1[1];
181 z1 = m_endpointsElectrons[i].x1[2];
182 t1 = m_endpointsElectrons[i].t1;
183 status = m_endpointsElectrons[i].status;
184}

Referenced by GarfieldPhysics::DoIt().

◆ GetHoleEndpoint()

void Garfield::AvalancheMC::GetHoleEndpoint ( const unsigned int  i,
double &  x0,
double &  y0,
double &  z0,
double &  t0,
double &  x1,
double &  y1,
double &  z1,
double &  t1,
int &  status 
) const

Definition at line 127 of file AvalancheMC.cc.

130 {
131 if (i >= m_endpointsHoles.size()) {
132 std::cerr << m_className << "::GetHoleEndpoint: Index out of range.\n";
133 return;
134 }
135
136 x0 = m_endpointsHoles[i].x0[0];
137 y0 = m_endpointsHoles[i].x0[1];
138 z0 = m_endpointsHoles[i].x0[2];
139 t0 = m_endpointsHoles[i].t0;
140 x1 = m_endpointsHoles[i].x1[0];
141 y1 = m_endpointsHoles[i].x1[1];
142 z1 = m_endpointsHoles[i].x1[2];
143 t1 = m_endpointsHoles[i].t1;
144 status = m_endpointsHoles[i].status;
145}

◆ GetIonEndpoint()

void Garfield::AvalancheMC::GetIonEndpoint ( const unsigned int  i,
double &  x0,
double &  y0,
double &  z0,
double &  t0,
double &  x1,
double &  y1,
double &  z1,
double &  t1,
int &  status 
) const

Definition at line 147 of file AvalancheMC.cc.

149 {
150 if (i >= m_endpointsIons.size()) {
151 std::cerr << m_className << "::GetIonEndpoint: Index out of range.\n";
152 return;
153 }
154
155 x0 = m_endpointsIons[i].x0[0];
156 y0 = m_endpointsIons[i].x0[1];
157 z0 = m_endpointsIons[i].x0[2];
158 t0 = m_endpointsIons[i].t0;
159 x1 = m_endpointsIons[i].x1[0];
160 y1 = m_endpointsIons[i].x1[1];
161 z1 = m_endpointsIons[i].x1[2];
162 t1 = m_endpointsIons[i].t1;
163 status = m_endpointsIons[i].status;
164}

◆ GetNumberOfDriftLinePoints()

unsigned int Garfield::AvalancheMC::GetNumberOfDriftLinePoints ( ) const
inline

Return the number of points along the last simulated drift line.

Definition at line 117 of file AvalancheMC.hh.

117{ return m_drift.size(); }

◆ GetNumberOfElectronEndpoints()

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

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

Definition at line 124 of file AvalancheMC.hh.

124 {
125 return m_endpointsElectrons.size();
126 }

◆ GetNumberOfHoleEndpoints()

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

Return the number of hole trajectories in the last simulated avalanche (including captured holes).

Definition at line 129 of file AvalancheMC.hh.

129 {
130 return m_endpointsHoles.size();
131 }

◆ GetNumberOfIonEndpoints()

unsigned int Garfield::AvalancheMC::GetNumberOfIonEndpoints ( ) const
inline

Return the number of ion trajectories.

Definition at line 133 of file AvalancheMC.hh.

133 {
134 return m_endpointsIons.size();
135 }

◆ ResumeAvalanche()

bool Garfield::AvalancheMC::ResumeAvalanche ( const bool  electron = true,
const bool  hole = true 
)

Definition at line 541 of file AvalancheMC.cc.

541 {
542
543 std::vector<DriftPoint> aval;
544 for (const auto& p: m_endpointsElectrons) {
545 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
546 AddPoint(p.x1, p.t1, Particle::Electron, 1, aval);
547 }
548 }
549 for (const auto& p: m_endpointsHoles) {
550 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
551 AddPoint(p.x1, p.t1, Particle::Hole, 1, aval);
552 }
553 }
554 for (const auto& p: m_endpointsIons) {
555 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
556 AddPoint(p.x1, p.t1, Particle::Ion, 1, aval);
557 }
558 }
559 return Avalanche(aval, electrons, holes);
560}

◆ SetCollisionSteps()

void Garfield::AvalancheMC::SetCollisionSteps ( const unsigned int  n = 100)

Use exponentially distributed time steps with mean equal to the specified multiple of the collision time (default model).

Definition at line 77 of file AvalancheMC.cc.

77 {
78 m_stepModel = StepModel::CollisionTime;
79 if (n < 1) {
80 std::cerr << m_className << "::SetCollisionSteps:\n "
81 << "Number of collisions set to default value (100).\n";
82 m_nMc = 100;
83 return;
84 }
85 if (m_debug) {
86 std::cout << m_className << "::SetCollisionSteps:\n "
87 << "Number of collisions to be skipped set to " << n << ".\n";
88 }
89 m_nMc = n;
90}

◆ SetDistanceSteps()

void Garfield::AvalancheMC::SetDistanceSteps ( const double  d = 0.001)

Use fixed distance steps (default 10 um).

Definition at line 62 of file AvalancheMC.cc.

62 {
63 m_stepModel = StepModel::FixedDistance;
64 if (d < Small) {
65 std::cerr << m_className << "::SetDistanceSteps:\n "
66 << "Step size is too small. Using default (10 um) instead.\n";
67 m_dMc = 0.001;
68 return;
69 }
70 if (m_debug) {
71 std::cout << m_className << "::SetDistanceSteps:\n"
72 << " Step size set to " << d << " cm.\n";
73 }
74 m_dMc = d;
75}

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

◆ SetElectronSignalScalingFactor()

void Garfield::AvalancheMC::SetElectronSignalScalingFactor ( const double  scale)
inline

Set multiplication factor for the signal induced by electrons.

Definition at line 104 of file AvalancheMC.hh.

104{ m_scaleE = scale; }

◆ SetHoleSignalScalingFactor()

void Garfield::AvalancheMC::SetHoleSignalScalingFactor ( const double  scale)
inline

Set multiplication factor for the signal induced by holes.

Definition at line 106 of file AvalancheMC.hh.

106{ m_scaleH = scale; }

◆ SetIonSignalScalingFactor()

void Garfield::AvalancheMC::SetIonSignalScalingFactor ( const double  scale)
inline

Set multiplication factor for the signal induced by ions.

Definition at line 108 of file AvalancheMC.hh.

108{ m_scaleI = scale; }

◆ SetSensor()

void Garfield::AvalancheMC::SetSensor ( Sensor s)

Set the sensor.

Definition at line 29 of file AvalancheMC.cc.

29 {
30 if (!sensor) {
31 std::cerr << m_className << "::SetSensor: Null pointer.\n";
32 return;
33 }
34
35 m_sensor = sensor;
36}

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

◆ SetSignalAveragingOrder()

void Garfield::AvalancheMC::SetSignalAveragingOrder ( const unsigned int  navg)
inline

Set the number of points to be used when averaging the signal vector over a time bin in the Sensor class. The averaging is done with a $2\times navg + 1$ point Newton-Raphson integration. Default: 1.

Definition at line 38 of file AvalancheMC.hh.

38{ m_navg = navg; }

◆ SetStepDistanceFunction()

void Garfield::AvalancheMC::SetStepDistanceFunction ( double(*)(double x, double y, double z)  f)

Retrieve the step distance from a user-supplied function.

Definition at line 92 of file AvalancheMC.cc.

93 {
94 if (!f) {
95 std::cerr << m_className << "::SetStepDistanceFunction: Null pointer.\n";
96 return;
97 }
98 m_fStep = f;
99 m_stepModel = StepModel::UserDistance;
100}

◆ SetTimeSteps()

void Garfield::AvalancheMC::SetTimeSteps ( const double  d = 0.02)

Use fixed-time steps (default 20 ps).

Definition at line 47 of file AvalancheMC.cc.

47 {
48 m_stepModel = StepModel::FixedTime;
49 if (d < Small) {
50 std::cerr << m_className << "::SetTimeSteps:\n "
51 << "Step size is too small. Using default (20 ps) instead.\n";
52 m_tMc = 0.02;
53 return;
54 }
55 if (m_debug) {
56 std::cout << m_className << "::SetTimeSteps:\n"
57 << " Step size set to " << d << " ns.\n";
58 }
59 m_tMc = d;
60}

◆ SetTimeWindow()

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

Define a time interval (only carriers inside the interval are drifted).

Definition at line 102 of file AvalancheMC.cc.

102 {
103 if (fabs(t1 - t0) < Small) {
104 std::cerr << m_className << "::SetTimeWindow:\n"
105 << " Time interval must be greater than zero.\n";
106 return;
107 }
108
109 m_tMin = std::min(t0, t1);
110 m_tMax = std::max(t0, t1);
111 m_hasTimeWindow = true;
112}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

◆ UnsetTimeWindow()

void Garfield::AvalancheMC::UnsetTimeWindow ( )
inline

Do not limit the time interval within which carriers are drifted.

Definition at line 101 of file AvalancheMC.hh.

101{ m_hasTimeWindow = false; }

◆ UseWeightingPotential()

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

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

Definition at line 41 of file AvalancheMC.hh.

41 {
42 m_useWeightingPotential = on;
43 }

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