Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::AvalancheMC Class Reference

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

#include <AvalancheMC.hh>

Classes

struct  EndPoint
struct  Point
struct  Seed

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, const size_t w=1)
 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, const size_t w=1)
 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, const size_t w=1)
 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, const size_t w=1)
 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, const size_t w=1)
 Simulate an avalanche initiated by an electron at a given starting point.
bool AvalancheHole (const double x, const double y, const double z, const double t, const bool electron=false, const size_t w=1)
 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, const size_t w=1)
 Simulate an avalanche initiated by an electron-hole pair.
void AddElectron (const double x, const double y, const double z, const double t, const size_t w=1)
 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, const size_t w=1)
 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, const size_t w=1)
 Add an ion to the list of particles to be transported.
void AddNegativeIon (const double x, const double y, const double z, const double t, const size_t w=1)
 Add an negative 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
 Return the number of electron trajectories in the last simulated avalanche (including captured electrons).
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
 Return the coordinates and time of start and end point of a given electron drift line.
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 GetNegativeIonEndpoint (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)
 Set the number of points to be used when averaging the signal vector over a time bin in the Sensor class.
void UseWeightingPotential (const bool on=true)
 Use the weighting potential (as opposed to the weighting field) for calculating the induced signal.
void EnableInducedChargeCalculation (const bool on=true)
 Switch on calculation of induced charge (default: disabled).
void EnableRKFSteps (const bool on=true)
 Switch on Runge-Kutta-Fehlberg stepping (as opposed to simple straight-line steps.
void EnableProjectedPathIntegration (const bool on=true)
 Switch on equilibration of multiplication and attachment over the drift line (default: enabled).
void EnableDiffusion ()
 Switch on diffusion (default: enabled).
void DisableDiffusion ()
 Switch off diffusion.
void EnableAttachment ()
 Switch on attachment (default: enabled).
void DisableAttachment ()
 Switch off attachment and multiplication.
void EnableMultiplication (const bool on)
 Switch multiplication on/off.
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 EnableMobilityMap (const bool on=true)
 Retrieve the (low-field) mobility from the component.
void EnableVelocityMap (const bool on=true)
 Retrieve the drift velocity from the component.
void EnableAvalancheSizeLimit (const unsigned int size)
 Set a maximum avalanche size (ignore further multiplication once this size has been reached).
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)
 Use exponentially distributed time steps with mean equal to the specified multiple of the collision time (default model).
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.
std::pair< unsigned int, unsigned int > GetAvalancheSize () 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).

Private Types

enum class  StepModel { FixedTime , FixedDistance , CollisionTime , UserDistance }

Private Member Functions

int DriftLine (const Seed &seed, std::vector< Point > &path, std::vector< Seed > &secondaries, const bool aval, const bool signal) const
 Compute a single drift line.
bool TransportParticles (std::vector< Seed > &stack, const bool withElectrons, const bool withHoles, const bool aval)
 Compute an avalanche.
int GetField (const std::array< double, 3 > &x, std::array< double, 3 > &e, std::array< double, 3 > &b, Medium *&medium) const
 Compute electric and magnetic field at a given position.
double GetMobility (const Particle particle, Medium *medium, const std::array< double, 3 > &x) const
 Retrieve the low-field mobility.
bool GetVelocity (const Particle particle, Medium *medium, const std::array< double, 3 > &x, const std::array< double, 3 > &e, const std::array< double, 3 > &b, std::array< double, 3 > &v) const
 Compute the drift velocity.
bool GetDiffusion (const Particle particle, Medium *medium, const std::array< double, 3 > &e, const std::array< double, 3 > &b, double &dl, double &dt) const
 Compute the diffusion coefficients.
double GetAttachment (const Particle particle, Medium *medium, const std::array< double, 3 > &x, const std::array< double, 3 > &e, const std::array< double, 3 > &b) const
 Compute the attachment coefficient.
double GetTownsend (const Particle particle, Medium *medium, const std::array< double, 3 > &x, const std::array< double, 3 > &e, const std::array< double, 3 > &b) const
 Compute the Townsend coefficient.
void StepRKF (const Particle particle, const std::array< double, 3 > &x0, const std::array< double, 3 > &v0, const double dt, std::array< double, 3 > &xf, std::array< double, 3 > &vf, int &status) const
 Compute end point and effective velocity for a step.
void AddDiffusion (const double step, const double dl, const double dt, std::array< double, 3 > &x, const std::array< double, 3 > &v) const
 Add a diffusion step.
void Terminate (const std::array< double, 3 > &x0, const double t0, std::array< double, 3 > &x, double &t) const
 Terminate a drift line close to the boundary.
bool ComputeGainLoss (const Particle ptype, const size_t w, std::vector< Point > &path, int &status, std::vector< Seed > &secondaries, const bool semiconductor=false) const
 Compute multiplication and losses along the current drift line.
bool ComputeAlphaEta (const Particle ptype, std::vector< Point > &path, std::vector< double > &alphas, std::vector< double > &etas) const
 Compute Townsend and attachment coefficients along the current drift line.
bool Equilibrate (std::vector< double > &alphas) const
void ComputeSignal (const double q, const std::vector< Point > &path) const
 Compute the induced signal for the current drift line.
void ComputeInducedCharge (const double q, const std::vector< Point > &path) const
 Compute the induced charge for the current drift line.
void PrintError (const std::string &fcn, const std::string &par, const Particle particle, const std::array< double, 3 > &x) const

Private Attributes

std::string m_className = "AvalancheMC"
Sensor * m_sensor = nullptr
StepModel m_stepModel = StepModel::CollisionTime
 Step size model.
double m_tMc = 0.02
 Fixed time step.
double m_dMc = 0.001
 Fixed distance step.
int m_nMc = 100
 Sample step size according to collision time.
double(* m_fStep )(double x, double y, double z) = nullptr
 User function returning the step size.
bool m_hasTimeWindow = false
 Flag whether a time window should be used.
double m_tMin = 0.
 Lower limit of the time window.
double m_tMax = 0.
 Upper limit of the time window.
unsigned int m_sizeCut = 0
 Max. avalanche size.
unsigned int m_nElectrons = 0
 Number of electrons produced.
unsigned int m_nHoles = 0
 Number of holes produced.
unsigned int m_nIons = 0
 Number of ions produced.
unsigned int m_nNegativeIons = 0
 Number of negative ions produced.
std::vector< EndPointm_electrons
 Start/end points of all electrons in the avalanche (including captured ones).
std::vector< EndPointm_holes
 Start/end points of all holes in the avalanche (including captured ones).
std::vector< EndPointm_ions
 Start/end points of all positive ions in the avalanche.
std::vector< EndPointm_negativeIons
 Start/end points of all negative ions in the avalanche.
ViewDriftm_viewer = nullptr
bool m_storeDriftLines = false
bool m_doSignal = true
unsigned int m_navg = 1
bool m_useWeightingPotential = true
bool m_doInducedCharge = false
bool m_doEquilibration = true
bool m_doRKF = false
bool m_useDiffusion = true
bool m_useAttachment = true
bool m_useMultiplication = true
double m_scaleE = 1.
 Scaling factor for electron signals.
double m_scaleH = 1.
 Scaling factor for hole signals.
double m_scaleI = 1.
 Scaling factor for ion signals.
bool m_useTownsendMap = false
 Take Townsend coefficients from the component.
bool m_useAttachmentMap = false
 Take attachment coefficients from the component.
bool m_useMobilityMap = false
 Take mobility coefficients from the component.
bool m_useVelocityMap = false
 Take the drift velocities from the component.
bool m_debug = false

Detailed Description

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

Definition at line 18 of file AvalancheMC.hh.

Member Enumeration Documentation

◆ StepModel

enum class Garfield::AvalancheMC::StepModel
strongprivate
Enumerator
FixedTime 
FixedDistance 
CollisionTime 
UserDistance 

Definition at line 224 of file AvalancheMC.hh.

224 {
225 FixedTime,
226 FixedDistance,
227 CollisionTime,
228 UserDistance
229 };

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.

◆ AvalancheMC() [2/2]

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

Constructor.

◆ ~AvalancheMC()

Garfield::AvalancheMC::~AvalancheMC ( )
inline

Destructor.

Definition at line 25 of file AvalancheMC.hh.

25{}

Member Function Documentation

◆ AddDiffusion()

void Garfield::AvalancheMC::AddDiffusion ( const double step,
const double dl,
const double dt,
std::array< double, 3 > & x,
const std::array< double, 3 > & v ) const
private

Add a diffusion step.

◆ AddElectron()

void Garfield::AvalancheMC::AddElectron ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ AddHole()

void Garfield::AvalancheMC::AddHole ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ AddIon()

void Garfield::AvalancheMC::AddIon ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ AddNegativeIon()

void Garfield::AvalancheMC::AddNegativeIon ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ AvalancheElectron()

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

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
wmultiplicity of the initial electron

◆ AvalancheElectronHole()

bool Garfield::AvalancheMC::AvalancheElectronHole ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

Simulate an avalanche initiated by an electron-hole pair.

◆ AvalancheHole()

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

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

◆ ComputeAlphaEta()

bool Garfield::AvalancheMC::ComputeAlphaEta ( const Particle ptype,
std::vector< Point > & path,
std::vector< double > & alphas,
std::vector< double > & etas ) const
private

Compute Townsend and attachment coefficients along the current drift line.

◆ ComputeGainLoss()

bool Garfield::AvalancheMC::ComputeGainLoss ( const Particle ptype,
const size_t w,
std::vector< Point > & path,
int & status,
std::vector< Seed > & secondaries,
const bool semiconductor = false ) const
private

Compute multiplication and losses along the current drift line.

◆ ComputeInducedCharge()

void Garfield::AvalancheMC::ComputeInducedCharge ( const double q,
const std::vector< Point > & path ) const
private

Compute the induced charge for the current drift line.

◆ ComputeSignal()

void Garfield::AvalancheMC::ComputeSignal ( const double q,
const std::vector< Point > & path ) const
private

Compute the induced signal for the current drift line.

◆ DisableAttachment()

void Garfield::AvalancheMC::DisableAttachment ( )
inline

Switch off attachment and multiplication.

Definition at line 164 of file AvalancheMC.hh.

164{ m_useAttachment = false; }

◆ DisableAvalancheSizeLimit()

void Garfield::AvalancheMC::DisableAvalancheSizeLimit ( )
inline

Do not limit the maximum avalanche size.

Definition at line 181 of file AvalancheMC.hh.

181{ m_sizeCut = 0; }
unsigned int m_sizeCut
Max. avalanche size.

◆ DisableDiffusion()

void Garfield::AvalancheMC::DisableDiffusion ( )
inline

Switch off diffusion.

Definition at line 159 of file AvalancheMC.hh.

159{ m_useDiffusion = false; }

◆ DisablePlotting()

void Garfield::AvalancheMC::DisablePlotting ( )
inline

Switch off drift line plotting.

Definition at line 123 of file AvalancheMC.hh.

123{ m_viewer = nullptr; }

◆ DriftElectron()

bool Garfield::AvalancheMC::DriftElectron ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ DriftHole()

bool Garfield::AvalancheMC::DriftHole ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ DriftIon()

bool Garfield::AvalancheMC::DriftIon ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ DriftLine()

int Garfield::AvalancheMC::DriftLine ( const Seed & seed,
std::vector< Point > & path,
std::vector< Seed > & secondaries,
const bool aval,
const bool signal ) const
private

Compute a single drift line.

◆ DriftNegativeIon()

bool Garfield::AvalancheMC::DriftNegativeIon ( const double x,
const double y,
const double z,
const double t,
const size_t w = 1 )

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

◆ EnableAttachment()

void Garfield::AvalancheMC::EnableAttachment ( )
inline

Switch on attachment (default: enabled).

Definition at line 162 of file AvalancheMC.hh.

162{ m_useAttachment = true; }

◆ EnableAttachmentMap()

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

Retrieve the attachment coefficient from the component.

Definition at line 171 of file AvalancheMC.hh.

171{ m_useAttachmentMap = on; }
bool m_useAttachmentMap
Take attachment coefficients from the component.

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

179{ m_sizeCut = size; }

◆ EnableDebugging()

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

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

Definition at line 217 of file AvalancheMC.hh.

217{ m_debug = on; }

◆ EnableDiffusion()

void Garfield::AvalancheMC::EnableDiffusion ( )
inline

Switch on diffusion (default: enabled).

Definition at line 157 of file AvalancheMC.hh.

157{ m_useDiffusion = true; }

◆ EnableDriftLines()

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

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

Definition at line 126 of file AvalancheMC.hh.

◆ EnableInducedChargeCalculation()

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

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

Definition at line 142 of file AvalancheMC.hh.

142 {
144 }

◆ EnableMobilityMap()

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

Retrieve the (low-field) mobility from the component.

Definition at line 173 of file AvalancheMC.hh.

173{ m_useMobilityMap = on; }
bool m_useMobilityMap
Take mobility coefficients from the component.

◆ EnableMultiplication()

void Garfield::AvalancheMC::EnableMultiplication ( const bool on)
inline

Switch multiplication on/off.

Definition at line 166 of file AvalancheMC.hh.

◆ EnablePlotting()

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

Switch on drift line plotting.

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

152 {
154 }

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

148{ 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 129 of file AvalancheMC.hh.

129{ m_doSignal = on; }

◆ EnableTownsendMap()

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

Retrieve the Townsend coefficient from the component.

Definition at line 169 of file AvalancheMC.hh.

169{ m_useTownsendMap = on; }
bool m_useTownsendMap
Take Townsend coefficients from the component.

◆ EnableVelocityMap()

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

Retrieve the drift velocity from the component.

Definition at line 175 of file AvalancheMC.hh.

175{ m_useVelocityMap = on; }
bool m_useVelocityMap
Take the drift velocities from the component.

◆ Equilibrate()

bool Garfield::AvalancheMC::Equilibrate ( std::vector< double > & alphas) const
private

◆ GetAttachment()

double Garfield::AvalancheMC::GetAttachment ( const Particle particle,
Medium * medium,
const std::array< double, 3 > & x,
const std::array< double, 3 > & e,
const std::array< double, 3 > & b ) const
private

Compute the attachment coefficient.

◆ GetAvalancheSize() [1/2]

std::pair< unsigned int, unsigned int > Garfield::AvalancheMC::GetAvalancheSize ( ) const
inline

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

Definition at line 213 of file AvalancheMC.hh.

213 {
214 return std::make_pair(m_nElectrons, std::max(m_nIons, m_nHoles));
215 }
unsigned int m_nElectrons
Number of electrons produced.
unsigned int m_nIons
Number of ions produced.
unsigned int m_nHoles
Number of holes produced.

◆ GetAvalancheSize() [2/2]

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

208 {
209 ne = m_nElectrons;
210 ni = std::max(m_nIons, m_nHoles);
211 }

◆ GetAvalancheSizeLimit()

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

Return the currently set avalanche size limit.

Definition at line 183 of file AvalancheMC.hh.

183{ return m_sizeCut; }

◆ GetDiffusion()

bool Garfield::AvalancheMC::GetDiffusion ( const Particle particle,
Medium * medium,
const std::array< double, 3 > & e,
const std::array< double, 3 > & b,
double & dl,
double & dt ) const
private

Compute the diffusion coefficients.

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

◆ GetElectrons()

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

Definition at line 90 of file AvalancheMC.hh.

90{ return m_electrons; }
std::vector< EndPoint > m_electrons
Start/end points of all electrons in the avalanche (including captured ones).

◆ GetField()

int Garfield::AvalancheMC::GetField ( const std::array< double, 3 > & x,
std::array< double, 3 > & e,
std::array< double, 3 > & b,
Medium *& medium ) const
private

Compute electric and magnetic field at a given position.

◆ GetHoles()

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

Definition at line 91 of file AvalancheMC.hh.

91{ return m_holes; }
std::vector< EndPoint > m_holes
Start/end points of all holes in the avalanche (including captured ones).

◆ 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

◆ GetIons()

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

Definition at line 92 of file AvalancheMC.hh.

92{ return m_ions; }
std::vector< EndPoint > m_ions
Start/end points of all positive ions in the avalanche.

◆ GetMobility()

double Garfield::AvalancheMC::GetMobility ( const Particle particle,
Medium * medium,
const std::array< double, 3 > & x ) const
private

Retrieve the low-field mobility.

◆ GetNegativeIonEndpoint()

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

◆ GetNegativeIons()

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

Definition at line 93 of file AvalancheMC.hh.

93 {
94 return m_negativeIons;
95 }
std::vector< EndPoint > m_negativeIons
Start/end points of all negative ions in the avalanche.

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

99{ return m_electrons.size(); }

◆ GetNumberOfIonEndpoints()

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

Return the number of ion trajectories.

Definition at line 101 of file AvalancheMC.hh.

101{ return m_ions.size(); }

◆ GetTownsend()

double Garfield::AvalancheMC::GetTownsend ( const Particle particle,
Medium * medium,
const std::array< double, 3 > & x,
const std::array< double, 3 > & e,
const std::array< double, 3 > & b ) const
private

Compute the Townsend coefficient.

◆ GetVelocity()

bool Garfield::AvalancheMC::GetVelocity ( const Particle particle,
Medium * medium,
const std::array< double, 3 > & x,
const std::array< double, 3 > & e,
const std::array< double, 3 > & b,
std::array< double, 3 > & v ) const
private

Compute the drift velocity.

◆ PrintError()

void Garfield::AvalancheMC::PrintError ( const std::string & fcn,
const std::string & par,
const Particle particle,
const std::array< double, 3 > & x ) const
private

◆ ResumeAvalanche()

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

Resume the simulation from the current set of charge carriers.

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

◆ SetDistanceSteps()

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

Use fixed distance steps (default 10 um).

◆ SetElectronSignalScalingFactor()

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

Set multiplication factor for the signal induced by electrons.

Definition at line 201 of file AvalancheMC.hh.

201{ m_scaleE = scale; }
double m_scaleE
Scaling factor for electron signals.

◆ SetHoleSignalScalingFactor()

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

Set multiplication factor for the signal induced by holes.

Definition at line 203 of file AvalancheMC.hh.

203{ m_scaleH = scale; }
double m_scaleH
Scaling factor for hole signals.

◆ SetIonSignalScalingFactor()

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

Set multiplication factor for the signal induced by ions.

Definition at line 205 of file AvalancheMC.hh.

205{ m_scaleI = scale; }
double m_scaleI
Scaling factor for ion signals.

◆ SetSensor()

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

Set the sensor.

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

134{ m_navg = navg; }

◆ SetStepDistanceFunction()

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

Retrieve the step distance from a user-supplied function.

◆ SetTimeSteps()

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

Use fixed-time steps (default 20 ps).

◆ SetTimeWindow()

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

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

◆ StepRKF()

void Garfield::AvalancheMC::StepRKF ( const Particle particle,
const std::array< double, 3 > & x0,
const std::array< double, 3 > & v0,
const double dt,
std::array< double, 3 > & xf,
std::array< double, 3 > & vf,
int & status ) const
private

Compute end point and effective velocity for a step.

◆ Terminate()

void Garfield::AvalancheMC::Terminate ( const std::array< double, 3 > & x0,
const double t0,
std::array< double, 3 > & x,
double & t ) const
private

Terminate a drift line close to the boundary.

◆ TransportParticles()

bool Garfield::AvalancheMC::TransportParticles ( std::vector< Seed > & stack,
const bool withElectrons,
const bool withHoles,
const bool aval )
private

Compute an avalanche.

◆ UnsetTimeWindow()

void Garfield::AvalancheMC::UnsetTimeWindow ( )
inline

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

Definition at line 198 of file AvalancheMC.hh.

198{ m_hasTimeWindow = false; }
bool m_hasTimeWindow
Flag whether a time window should be used.

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

137 {
139 }

Member Data Documentation

◆ m_className

std::string Garfield::AvalancheMC::m_className = "AvalancheMC"
private

Definition at line 220 of file AvalancheMC.hh.

◆ m_debug

bool Garfield::AvalancheMC::m_debug = false
private

Definition at line 300 of file AvalancheMC.hh.

◆ m_dMc

double Garfield::AvalancheMC::m_dMc = 0.001
private

Fixed distance step.

Definition at line 236 of file AvalancheMC.hh.

◆ m_doEquilibration

bool Garfield::AvalancheMC::m_doEquilibration = true
private

Definition at line 279 of file AvalancheMC.hh.

◆ m_doInducedCharge

bool Garfield::AvalancheMC::m_doInducedCharge = false
private

Definition at line 278 of file AvalancheMC.hh.

◆ m_doRKF

bool Garfield::AvalancheMC::m_doRKF = false
private

Definition at line 280 of file AvalancheMC.hh.

◆ m_doSignal

bool Garfield::AvalancheMC::m_doSignal = true
private

Definition at line 275 of file AvalancheMC.hh.

◆ m_electrons

std::vector<EndPoint> Garfield::AvalancheMC::m_electrons
private

Start/end points of all electrons in the avalanche (including captured ones).

Definition at line 263 of file AvalancheMC.hh.

◆ m_fStep

double(* Garfield::AvalancheMC::m_fStep) (double x, double y, double z) = nullptr
private

User function returning the step size.

Definition at line 240 of file AvalancheMC.hh.

◆ m_hasTimeWindow

bool Garfield::AvalancheMC::m_hasTimeWindow = false
private

Flag whether a time window should be used.

Definition at line 243 of file AvalancheMC.hh.

◆ m_holes

std::vector<EndPoint> Garfield::AvalancheMC::m_holes
private

Start/end points of all holes in the avalanche (including captured ones).

Definition at line 266 of file AvalancheMC.hh.

◆ m_ions

std::vector<EndPoint> Garfield::AvalancheMC::m_ions
private

Start/end points of all positive ions in the avalanche.

Definition at line 268 of file AvalancheMC.hh.

◆ m_navg

unsigned int Garfield::AvalancheMC::m_navg = 1
private

Definition at line 276 of file AvalancheMC.hh.

◆ m_negativeIons

std::vector<EndPoint> Garfield::AvalancheMC::m_negativeIons
private

Start/end points of all negative ions in the avalanche.

Definition at line 270 of file AvalancheMC.hh.

◆ m_nElectrons

unsigned int Garfield::AvalancheMC::m_nElectrons = 0
private

Number of electrons produced.

Definition at line 253 of file AvalancheMC.hh.

◆ m_nHoles

unsigned int Garfield::AvalancheMC::m_nHoles = 0
private

Number of holes produced.

Definition at line 255 of file AvalancheMC.hh.

◆ m_nIons

unsigned int Garfield::AvalancheMC::m_nIons = 0
private

Number of ions produced.

Definition at line 257 of file AvalancheMC.hh.

◆ m_nMc

int Garfield::AvalancheMC::m_nMc = 100
private

Sample step size according to collision time.

Definition at line 238 of file AvalancheMC.hh.

◆ m_nNegativeIons

unsigned int Garfield::AvalancheMC::m_nNegativeIons = 0
private

Number of negative ions produced.

Definition at line 259 of file AvalancheMC.hh.

◆ m_scaleE

double Garfield::AvalancheMC::m_scaleE = 1.
private

Scaling factor for electron signals.

Definition at line 285 of file AvalancheMC.hh.

◆ m_scaleH

double Garfield::AvalancheMC::m_scaleH = 1.
private

Scaling factor for hole signals.

Definition at line 287 of file AvalancheMC.hh.

◆ m_scaleI

double Garfield::AvalancheMC::m_scaleI = 1.
private

Scaling factor for ion signals.

Definition at line 289 of file AvalancheMC.hh.

◆ m_sensor

Sensor* Garfield::AvalancheMC::m_sensor = nullptr
private

Definition at line 222 of file AvalancheMC.hh.

◆ m_sizeCut

unsigned int Garfield::AvalancheMC::m_sizeCut = 0
private

Max. avalanche size.

Definition at line 250 of file AvalancheMC.hh.

◆ m_stepModel

StepModel Garfield::AvalancheMC::m_stepModel = StepModel::CollisionTime
private

Step size model.

Definition at line 231 of file AvalancheMC.hh.

◆ m_storeDriftLines

bool Garfield::AvalancheMC::m_storeDriftLines = false
private

Definition at line 274 of file AvalancheMC.hh.

◆ m_tMax

double Garfield::AvalancheMC::m_tMax = 0.
private

Upper limit of the time window.

Definition at line 247 of file AvalancheMC.hh.

◆ m_tMc

double Garfield::AvalancheMC::m_tMc = 0.02
private

Fixed time step.

Definition at line 234 of file AvalancheMC.hh.

◆ m_tMin

double Garfield::AvalancheMC::m_tMin = 0.
private

Lower limit of the time window.

Definition at line 245 of file AvalancheMC.hh.

◆ m_useAttachment

bool Garfield::AvalancheMC::m_useAttachment = true
private

Definition at line 282 of file AvalancheMC.hh.

◆ m_useAttachmentMap

bool Garfield::AvalancheMC::m_useAttachmentMap = false
private

Take attachment coefficients from the component.

Definition at line 294 of file AvalancheMC.hh.

◆ m_useDiffusion

bool Garfield::AvalancheMC::m_useDiffusion = true
private

Definition at line 281 of file AvalancheMC.hh.

◆ m_useMobilityMap

bool Garfield::AvalancheMC::m_useMobilityMap = false
private

Take mobility coefficients from the component.

Definition at line 296 of file AvalancheMC.hh.

◆ m_useMultiplication

bool Garfield::AvalancheMC::m_useMultiplication = true
private

Definition at line 283 of file AvalancheMC.hh.

◆ m_useTownsendMap

bool Garfield::AvalancheMC::m_useTownsendMap = false
private

Take Townsend coefficients from the component.

Definition at line 292 of file AvalancheMC.hh.

◆ m_useVelocityMap

bool Garfield::AvalancheMC::m_useVelocityMap = false
private

Take the drift velocities from the component.

Definition at line 298 of file AvalancheMC.hh.

◆ m_useWeightingPotential

bool Garfield::AvalancheMC::m_useWeightingPotential = true
private

Definition at line 277 of file AvalancheMC.hh.

◆ m_viewer

ViewDrift* Garfield::AvalancheMC::m_viewer = nullptr
private

Definition at line 272 of file AvalancheMC.hh.


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