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

#include <AvalancheMC.hh>

Classes

struct  EndPoint
 
struct  Point
 

Public Member Functions

 AvalancheMC ()
 Default constructor.
 
 AvalancheMC (Sensor *sensor)
 Constructor.
 
 ~AvalancheMC ()
 Destructor.
 
void SetSensor (Sensor *s)
 Set the sensor.
 
bool DriftElectron (const double x, const double y, const double z, const double t)
 Simulate the drift line of an electron from a given starting point.
 
bool DriftHole (const double x, const double y, const double z, const double t)
 Simulate the drift line of a hole from a given starting point.
 
bool DriftIon (const double x, const double y, const double z, const double t)
 Simulate the drift line of an ion from a given starting point.
 
bool DriftNegativeIon (const double x, const double y, const double z, const double t)
 Simulate the drift line of a negative ion from a given starting point.
 
bool AvalancheElectron (const double x, const double y, const double z, const double t, const bool hole=false)
 
bool AvalancheHole (const double x, const double y, const double z, const double t, const bool electron=false)
 Simulate an avalanche initiated by a hole at a given starting point.
 
bool AvalancheElectronHole (const double x, const double y, const double z, const double t)
 Simulate an avalanche initiated by an electron-hole pair.
 
void AddElectron (const double x, const double y, const double z, const double t)
 Add an electron to the list of particles to be transported.
 
void AddHole (const double x, const double y, const double z, const double t)
 Add a hole to the list of particles to be transported.
 
void AddIon (const double x, const double y, const double z, const double t)
 Add an ion to the list of particles to be transported.
 
bool ResumeAvalanche (const bool electron=true, const bool hole=true)
 Resume the simulation from the current set of charge carriers.
 
const std::vector< EndPoint > & GetElectrons () const
 
const std::vector< EndPoint > & GetHoles () const
 
const std::vector< EndPoint > & GetIons () const
 
const std::vector< EndPoint > & GetNegativeIons () const
 
size_t GetNumberOfElectronEndpoints () const
 
size_t GetNumberOfHoleEndpoints () const
 
size_t GetNumberOfIonEndpoints () const
 Return the number of ion trajectories.
 
void GetElectronEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
 
void GetHoleEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
 
void GetIonEndpoint (const size_t i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
 
void EnablePlotting (ViewDrift *view)
 Switch on drift line plotting.
 
void DisablePlotting ()
 Switch off drift line plotting.
 
void EnableDriftLines (const bool on=true)
 Switch on storage of drift lines (default: off).
 
void EnableSignalCalculation (const bool on=true)
 Switch calculation of induced currents on or off (default: enabled).
 
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 EnableTownsendMap (const bool on=true)
 Retrieve the Townsend coefficient from the component.
 
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 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.
 
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 18 of file AvalancheMC.hh.

Constructor & Destructor Documentation

◆ AvalancheMC() [1/2]

Garfield::AvalancheMC::AvalancheMC ( )
inline

Default constructor.

Definition at line 21 of file AvalancheMC.hh.

21: AvalancheMC(nullptr) {}
AvalancheMC()
Default constructor.

Referenced by AvalancheMC().

◆ AvalancheMC() [2/2]

Garfield::AvalancheMC::AvalancheMC ( Sensor * sensor)

Constructor.

Definition at line 66 of file AvalancheMC.cc.

66: m_sensor(sensor) { }

◆ ~AvalancheMC()

Garfield::AvalancheMC::~AvalancheMC ( )
inline

Destructor.

Definition at line 25 of file AvalancheMC.hh.

25{}

Member Function Documentation

◆ AddElectron()

void Garfield::AvalancheMC::AddElectron ( const double x,
const double y,
const double z,
const double t )

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

Definition at line 554 of file AvalancheMC.cc.

555 {
556
557 EndPoint p;
558 p.status = StatusAlive;
559 p.path = {MakePoint(x, y, z, t)};
560 m_electrons.push_back(std::move(p));
561 // TODO
562 ++m_nElectrons;
563}

◆ AddHole()

void Garfield::AvalancheMC::AddHole ( const double x,
const double y,
const double z,
const double t )

Add a hole to the list of particles to be transported.

Definition at line 565 of file AvalancheMC.cc.

566 {
567 EndPoint p;
568 p.status = StatusAlive;
569 p.path = {MakePoint(x, y, z, t)};
570 m_holes.push_back(std::move(p));
571 // TODO
572 ++m_nHoles;
573}

◆ AddIon()

void Garfield::AvalancheMC::AddIon ( const double x,
const double y,
const double z,
const double t )

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

Definition at line 574 of file AvalancheMC.cc.

575 {
576 EndPoint p;
577 p.status = StatusAlive;
578 p.path = {MakePoint(x, y, z, t)};
579 m_ions.push_back(std::move(p));
580 // TODO
581 ++m_nIons;
582}

◆ AvalancheElectron()

bool Garfield::AvalancheMC::AvalancheElectron ( const double x,
const double y,
const double z,
const double t,
const bool hole = false )

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

Parameters
x,y,z,tcoordinates and time of the initial electron
holesimulate the hole component of the avalanche or not

Definition at line 526 of file AvalancheMC.cc.

528 {
529 std::vector<std::pair<Point, Particle> > particles;
530 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
532 return TransportParticles(particles, true, holes, true);
533}

◆ AvalancheElectronHole()

bool Garfield::AvalancheMC::AvalancheElectronHole ( const double x,
const double y,
const double z,
const double t )

Simulate an avalanche initiated by an electron-hole pair.

Definition at line 544 of file AvalancheMC.cc.

545 {
546 std::vector<std::pair<Point, Particle> > particles;
547 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
549 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
551 return TransportParticles(particles, true, true, true);
552}

◆ AvalancheHole()

bool Garfield::AvalancheMC::AvalancheHole ( const double x,
const double y,
const double z,
const double t,
const bool electron = false )

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

Definition at line 535 of file AvalancheMC.cc.

537 {
538 std::vector<std::pair<Point, Particle> > particles;
539 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
541 return TransportParticles(particles, electrons, true, true);
542}

◆ DisableAttachment()

void Garfield::AvalancheMC::DisableAttachment ( )
inline

Switch off attachment and multiplication.

Definition at line 153 of file AvalancheMC.hh.

153{ m_useAttachment = false; }

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMC::DisableAvalancheSizeLimit ( )
inline

Do not limit the maximum avalanche size.

Definition at line 166 of file AvalancheMC.hh.

166{ m_sizeCut = 0; }

◆ DisableDiffusion()

void Garfield::AvalancheMC::DisableDiffusion ( )
inline

Switch off diffusion.

Definition at line 147 of file AvalancheMC.hh.

147{ m_useDiffusion = false; }

◆ DisablePlotting()

void Garfield::AvalancheMC::DisablePlotting ( )
inline

Switch off drift line plotting.

Definition at line 111 of file AvalancheMC.hh.

111{ m_viewer = nullptr; }

◆ DriftElectron()

bool Garfield::AvalancheMC::DriftElectron ( const double x,
const double y,
const double z,
const double t )

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

Definition at line 212 of file AvalancheMC.cc.

213 {
214
215 std::vector<std::pair<Point, Particle> > particles;
216 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
218 return TransportParticles(particles, true, false, false);
219}

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

◆ DriftHole()

bool Garfield::AvalancheMC::DriftHole ( const double x,
const double y,
const double z,
const double t )

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

Definition at line 221 of file AvalancheMC.cc.

222 {
223 std::vector<std::pair<Point, Particle> > particles;
224 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
226 return TransportParticles(particles, false, true, false);
227}

Referenced by main().

◆ DriftIon()

bool Garfield::AvalancheMC::DriftIon ( const double x,
const double y,
const double z,
const double t )

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

Definition at line 229 of file AvalancheMC.cc.

230 {
231 std::vector<std::pair<Point, Particle> > particles;
232 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
234 return TransportParticles(particles, false, true, false);
235}

◆ DriftNegativeIon()

bool Garfield::AvalancheMC::DriftNegativeIon ( const double x,
const double y,
const double z,
const double t )

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

Definition at line 237 of file AvalancheMC.cc.

238 {
239 std::vector<std::pair<Point, Particle> > particles;
240 particles.emplace_back(std::make_pair(MakePoint(x0, y0, z0, t0),
242 return TransportParticles(particles, false, true, false);
243}

◆ 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 151 of file AvalancheMC.hh.

151{ m_useAttachment = true; }

◆ EnableAttachmentMap()

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

Retrieve the attachment coefficient from the component.

Definition at line 158 of file AvalancheMC.hh.

158{ 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 164 of file AvalancheMC.hh.

164{ m_sizeCut = size; }

◆ EnableDebugging()

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

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

Definition at line 199 of file AvalancheMC.hh.

199{ m_debug = on; }

◆ EnableDiffusion()

void Garfield::AvalancheMC::EnableDiffusion ( )
inline

Switch on diffusion (default: enabled).

Definition at line 145 of file AvalancheMC.hh.

145{ m_useDiffusion = true; }

◆ EnableDriftLines()

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

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

Definition at line 114 of file AvalancheMC.hh.

114{ m_storeDriftLines = on; }

◆ EnableInducedChargeCalculation()

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

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

Definition at line 130 of file AvalancheMC.hh.

130 {
131 m_doInducedCharge = on;
132 }

◆ EnablePlotting()

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

Switch on drift line plotting.

Definition at line 77 of file AvalancheMC.cc.

77 {
78 if (!view) {
79 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
80 return;
81 }
82
83 m_viewer = view;
84}

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 140 of file AvalancheMC.hh.

140 {
141 m_doEquilibration = on;
142 }

◆ 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 136 of file AvalancheMC.hh.

136{ m_doRKF = on; }

◆ EnableSignalCalculation()

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

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

Definition at line 117 of file AvalancheMC.hh.

117{ m_doSignal = on; }

◆ EnableTownsendMap()

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

Retrieve the Townsend coefficient from the component.

Definition at line 156 of file AvalancheMC.hh.

156{ m_useTownsendMap = on; }

◆ EnableVelocityMap()

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

Retrieve the drift velocity from the component.

Definition at line 160 of file AvalancheMC.hh.

160{ 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 193 of file AvalancheMC.hh.

193 {
194 ne = m_nElectrons;
195 ni = std::max(m_nIons, m_nHoles);
196 }

◆ GetAvalancheSizeLimit()

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

Return the currently set avalanche size limit.

Definition at line 168 of file AvalancheMC.hh.

168{ return m_sizeCut; }

◆ GetElectronEndpoint()

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

195 {
196 if (i >= m_electrons.size()) {
197 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
198 return;
199 }
200
201 x0 = m_electrons[i].path.front().x;
202 y0 = m_electrons[i].path.front().y;
203 z0 = m_electrons[i].path.front().z;
204 t0 = m_electrons[i].path.front().t;
205 x1 = m_electrons[i].path.back().x;
206 y1 = m_electrons[i].path.back().y;
207 z1 = m_electrons[i].path.back().z;
208 t1 = m_electrons[i].path.back().t;
209 status = m_electrons[i].status;
210}

Referenced by GarfieldPhysics::DoIt().

◆ GetElectrons()

const std::vector< EndPoint > & Garfield::AvalancheMC::GetElectrons ( ) const
inline

Definition at line 75 of file AvalancheMC.hh.

75{ return m_electrons; }

◆ GetHoleEndpoint()

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

Definition at line 153 of file AvalancheMC.cc.

156 {
157 if (i >= m_holes.size()) {
158 std::cerr << m_className << "::GetHoleEndpoint: Index out of range.\n";
159 return;
160 }
161
162 x0 = m_holes[i].path.front().x;
163 y0 = m_holes[i].path.front().y;
164 z0 = m_holes[i].path.front().z;
165 t0 = m_holes[i].path.front().t;
166 x1 = m_holes[i].path.back().x;
167 y1 = m_holes[i].path.back().y;
168 z1 = m_holes[i].path.back().z;
169 t1 = m_holes[i].path.back().t;
170 status = m_holes[i].status;
171}

◆ GetHoles()

const std::vector< EndPoint > & Garfield::AvalancheMC::GetHoles ( ) const
inline

Definition at line 76 of file AvalancheMC.hh.

76{ return m_holes; }

◆ GetIonEndpoint()

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

Definition at line 173 of file AvalancheMC.cc.

175 {
176 if (i >= m_ions.size()) {
177 std::cerr << m_className << "::GetIonEndpoint: Index out of range.\n";
178 return;
179 }
180
181 x0 = m_ions[i].path.front().x;
182 y0 = m_ions[i].path.front().y;
183 z0 = m_ions[i].path.front().z;
184 t0 = m_ions[i].path.front().t;
185 x1 = m_ions[i].path.back().x;
186 y1 = m_ions[i].path.back().y;
187 z1 = m_ions[i].path.back().z;
188 t1 = m_ions[i].path.back().t;
189 status = m_ions[i].status;
190}

◆ GetIons()

const std::vector< EndPoint > & Garfield::AvalancheMC::GetIons ( ) const
inline

Definition at line 77 of file AvalancheMC.hh.

77{ return m_ions; }

◆ GetNegativeIons()

const std::vector< EndPoint > & Garfield::AvalancheMC::GetNegativeIons ( ) const
inline

Definition at line 78 of file AvalancheMC.hh.

78 {
79 return m_negativeIons;
80 }

◆ GetNumberOfElectronEndpoints()

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

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

Definition at line 84 of file AvalancheMC.hh.

84{ return m_electrons.size(); }

◆ GetNumberOfHoleEndpoints()

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

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

Definition at line 87 of file AvalancheMC.hh.

87{ return m_holes.size(); }

◆ GetNumberOfIonEndpoints()

size_t Garfield::AvalancheMC::GetNumberOfIonEndpoints ( ) const
inline

Return the number of ion trajectories.

Definition at line 89 of file AvalancheMC.hh.

89{ return m_ions.size(); }

◆ ResumeAvalanche()

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

Resume the simulation from the current set of charge carriers.

Definition at line 584 of file AvalancheMC.cc.

584 {
585
586 std::vector<std::pair<Point, Particle> > particles;
587 for (const auto& p: m_electrons) {
588 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
589 particles.push_back(std::make_pair(p.path.back(), Particle::Electron));
590 }
591 }
592 for (const auto& p: m_holes) {
593 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
594 particles.push_back(std::make_pair(p.path.back(), Particle::Hole));
595 }
596 }
597 for (const auto& p: m_ions) {
598 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
599 particles.push_back(std::make_pair(p.path.back(), Particle::Ion));
600 }
601 }
602 for (const auto& p: m_negativeIons) {
603 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
604 particles.push_back(std::make_pair(p.path.back(), Particle::NegativeIon));
605 }
606 }
607 return TransportParticles(particles, electrons, holes, true);
608}

◆ 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 116 of file AvalancheMC.cc.

116 {
117 m_stepModel = StepModel::CollisionTime;
118 if (n < 1) {
119 std::cerr << m_className << "::SetCollisionSteps:\n "
120 << "Number of collisions set to default value (100).\n";
121 m_nMc = 100;
122 return;
123 }
124 if (m_debug) {
125 std::cout << m_className << "::SetCollisionSteps:\n "
126 << "Number of collisions to be skipped set to " << n << ".\n";
127 }
128 m_nMc = n;
129}

◆ SetDistanceSteps()

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

Use fixed distance steps (default 10 um).

Definition at line 101 of file AvalancheMC.cc.

101 {
102 m_stepModel = StepModel::FixedDistance;
103 if (d < Small) {
104 std::cerr << m_className << "::SetDistanceSteps:\n "
105 << "Step size is too small. Using default (10 um) instead.\n";
106 m_dMc = 0.001;
107 return;
108 }
109 if (m_debug) {
110 std::cout << m_className << "::SetDistanceSteps:\n"
111 << " Step size set to " << d << " cm.\n";
112 }
113 m_dMc = d;
114}

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 186 of file AvalancheMC.hh.

186{ m_scaleE = scale; }

◆ SetHoleSignalScalingFactor()

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

Set multiplication factor for the signal induced by holes.

Definition at line 188 of file AvalancheMC.hh.

188{ m_scaleH = scale; }

◆ SetIonSignalScalingFactor()

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

Set multiplication factor for the signal induced by ions.

Definition at line 190 of file AvalancheMC.hh.

190{ m_scaleI = scale; }

◆ SetSensor()

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

Set the sensor.

Definition at line 68 of file AvalancheMC.cc.

68 {
69 if (!sensor) {
70 std::cerr << m_className << "::SetSensor: Null pointer.\n";
71 return;
72 }
73
74 m_sensor = sensor;
75}

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 122 of file AvalancheMC.hh.

122{ m_navg = navg; }

◆ SetStepDistanceFunction()

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

Retrieve the step distance from a user-supplied function.

Definition at line 131 of file AvalancheMC.cc.

132 {
133 if (!f) {
134 std::cerr << m_className << "::SetStepDistanceFunction: Null pointer.\n";
135 return;
136 }
137 m_fStep = f;
138 m_stepModel = StepModel::UserDistance;
139}

◆ SetTimeSteps()

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

Use fixed-time steps (default 20 ps).

Definition at line 86 of file AvalancheMC.cc.

86 {
87 m_stepModel = StepModel::FixedTime;
88 if (d < Small) {
89 std::cerr << m_className << "::SetTimeSteps:\n "
90 << "Step size is too small. Using default (20 ps) instead.\n";
91 m_tMc = 0.02;
92 return;
93 }
94 if (m_debug) {
95 std::cout << m_className << "::SetTimeSteps:\n"
96 << " Step size set to " << d << " ns.\n";
97 }
98 m_tMc = d;
99}

◆ 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 141 of file AvalancheMC.cc.

141 {
142 if (fabs(t1 - t0) < Small) {
143 std::cerr << m_className << "::SetTimeWindow:\n"
144 << " Time interval must be greater than zero.\n";
145 return;
146 }
147
148 m_tMin = std::min(t0, t1);
149 m_tMax = std::max(t0, t1);
150 m_hasTimeWindow = true;
151}
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 183 of file AvalancheMC.hh.

183{ 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 125 of file AvalancheMC.hh.

125 {
126 m_useWeightingPotential = on;
127 }

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