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

Propagates avalanches with the 2d (axi-symmetric) space-charge routine from Lippmann, Riegler (2004) in uniform background fields. More...

#include <AvalancheGridSpaceCharge.hh>

Classes

struct  GridNode
struct  Point

Public Member Functions

 AvalancheGridSpaceCharge ()
 Default constructor.
 AvalancheGridSpaceCharge (Sensor *sensor)
 Constructor.
 ~AvalancheGridSpaceCharge ()=default
 Destructor.
void Reset ()
 Reset the charges.
void EnableDebugging (const bool option=true)
 Enable/disable debugging (log messages) (default off)
void EnableStickyAnode (const bool option=true)
 Enable sticky anode (default on)
void EnableTOF (const bool option=true)
 Enable to use TOF swarm parameters (default on)
void EnableDiffusion (const bool option=true)
 Enable diffusion (default off)
void EnableSpaceChargeEffect (const bool option=true)
 Enable space charge calculations (default on)
void EnableAdaptiveTimeStepping (const bool option=true)
 Enable adaptive time stepping (default on)
void EnableMC (const bool option=true)
 Enable Monte Carlo (gain/diffusion) up to a certain total number of electrons e.g.
void SetNCrit (const long NCrit=1e8)
 Set electron saturation value applied for each gas gap if space charge effect is turned off (default 1e8)
void SetFieldCalculation (const std::string &option="coulomb", const int nof_approx=1.)
 Sets the different options to calculate the space charge field coulomb: free field approximation (relaxation: fdm-relaxation methode (with a condition to stop)) mirror: symmetric 3 layer single gap rpc with metal - resistive layer - gas gap - r.
void SetK (float option=0.95)
 Set the streamer-inception criterion constant K in the interval (0, ∞) s.t.
void SetStopAtK (bool option=true)
 Stop the avalanche if K % field is reached.
void SetSensor (Sensor *sensor)
 Set the sensor (and determine if it includes a parallel-plate component).
void Set2dGrid (double zmin, double zmax, int zsteps, double rmax, int rsteps)
void AddElectrons (AvalancheMicroscopic *avmc)
 Import electron (no ions) data from AvalancheMicroscopic class to axi-symmetric grid.
void AddElectron (double x, double y, double z, double t=0, int n=1)
 Set n electrons onto the grid.
void AddExtraElectron (double y, int n=1)
 After calling AddElectron, add more electrons on the same transversal line (y-freedom).
void StartGridAvalanche (double dtime=-1)
 Starts the simulation with the imported electrons for a time step dt (dt = -1 until there are no electrons left in the gap).
long GetAvalancheSize () const
 Returns the total positive charge in the gap's.
double GetMeanDistance ()
 Return current mean distance of the electrons on the grid.
long ReachedKPercent () const
 Returns if 100 * K % background charge has been reached.
const std::vector< std::pair< double, long > > GetElectronEvolution () const
 Returns the total electron number evolution.
void ExportGrid (const std::string &filename)
 Export the current grid to a txt file (electron, ions numbers and field magnitude) Take care: only where electrons are located is the space-charge field evaluated!

Private Types

enum class  Elliptic : std::size_t { X , K , E }

Private Member Functions

void PrepareElectronsFromMicroscopicAvalanche ()
bool SnapTo2dGrid (double x, double y, double z, long n=1, int gasLayer=0)
void Prepare2dMesh ()
bool TransportTimeStep ()
void DiffuseTimeStep (double dx, long nElectron, double nPosIon, double nNegIon, int iz, int ir, int gasGap)
void DistributeCharges (long nElectron, double nPosIon, double nNegIon, int iz, int ir, double stepZ, double stepR, int gasGap)
void GetLocalField (int iz, int ir, double &eFieldZ, double &eFieldR, const std::string &fieldOption, int gasGap)
void GetFreeChargedRing (int iz, int ir, int fz, int fr, double &eFieldZ, double &eFieldR)
void GetFreeChargedRing (double zi, double ri, double zf, double rf, double &eFieldZ, double &eFieldR)
bool AddFieldFromChargeAt (int iz, int ir, int fz, int fr, double N, double &eFieldZ, double &eFieldR)
bool AddFieldFromChargeAt (int iz, int ir, double zf, double rf, double N, double &eFieldZ, double &eFieldR)
void GetSwarmParameters (double MagEField, double &alpha, double &eta, double &drift, double &dSigmaL, double &dSigmaT, double &wv, double &wr, double &alphaPT, double &etaPT, int gasGap)
void GetGlobalCoordinates (double r, double z, double phi, double &xg, double &yg, double &zg, int gasGap)
void ImportEllipticIntegralValues (const std::string &filename)
void GetEllipticIntegrals (double x, double &K, double &E)
int GetGasGapNumber (int layerIndex)

Private Attributes

std::string m_className = "AvalancheGridSpaceCharge"
bool m_bDebug = false
bool m_bDiffusion = false
bool m_bStick = true
bool m_bDriftAvalanche = false
bool m_bImportAvalanche = false
bool m_bPreparedImportAvalanche = false
long m_lNCrit = 1e8
bool m_bSpaceCharge = true
float m_fStreamerK = 0.95
bool m_bStopAtK = false
bool m_bFieldK = false
long m_lElectronsK
bool m_bAdaptiveTime = true
bool m_bImportElliptic = false
bool m_bUseTOF = true
 Flag if TOF parameters should be used, else Magboltz drift and SST spatial coefficients.
bool m_bWrAvailable = true
 Flag if bulk drift velocity is available to the simulation.
bool m_bRatesAvailable = true
 Flag if temporal rates are available to the simulation.
bool m_bMC = true
int m_iFieldApprox = 1
double m_dMinGroups = 50
ComponentParallelPlatem_pp = nullptr
Sensorm_sensor = nullptr
std::vector< double > m_zGrid
 Grid points of z-coordinate.
int m_zSteps = 0
 Number of grid points.
double m_zStepSize = 0.
std::vector< double > m_rGrid
 Distance between the grid points.
int m_rSteps = 0.
 Number of grid points.
double m_rStepSize = 0.
bool m_isgridset = false
 Distance between the grid points.
long m_nTotElectron = 0
 Total amount of electrons at time step.
long m_nTotPosIons = 0
 total amount of charge created
double m_time = 0.
 Clock.
double m_time0 = 0.
 Initial time.
double m_dt = 0.
 Time step.
bool m_run = true
 Tracking if the charges are still in the drift gap.
std::vector< std::vector< int > > m_zGasGapBoundaries
 [k] -> {izLeft, ..., izRight}
std::vector< std::vector< GridNode > > m_grid
 grid with nodes on it
std::vector< std::vector< Point > > m_vElectrons
 Electrons to transfer onto grid.
std::vector< std::pair< double, long > > m_vNElectronEvolution
std::vector< long > m_vGroupSizes
 same values as Lippmann & Riegler
std::vector< int > m_vIndexGasGaps = {0}
 Which layer indices are gas layers.
std::vector< std::vector< double > > m_vCoNGasLayer {}
 Coordinates of center of electron number.
std::vector< double > m_vYPointInGasGap {}
 Example point (y-coord) in each gas gap.
std::vector< double > m_ezBkg = {0}
 Uniform background field in z direction, can be negative.
std::vector< int > m_vSaturatedGaps {}
 Which gas gaps are saturated if saturation is on.
std::string m_sFieldOption = "coulomb"

Static Private Attributes

static const constexpr std::size_t elliptic_size {29981}
static const std::array< std::array< double, 3 >, elliptic_sizem_elliptic

Detailed Description

Propagates avalanches with the 2d (axi-symmetric) space-charge routine from Lippmann, Riegler (2004) in uniform background fields.

Different options to calculate space-charge-fields can be chosen.

Definition at line 18 of file AvalancheGridSpaceCharge.hh.

Member Enumeration Documentation

◆ Elliptic

enum class Garfield::AvalancheGridSpaceCharge::Elliptic : std::size_t
strongprivate
Enumerator

Definition at line 334 of file AvalancheGridSpaceCharge.hh.

335 {
336 X,
337 K,
338 E
339 };

Constructor & Destructor Documentation

◆ AvalancheGridSpaceCharge() [1/2]

Garfield::AvalancheGridSpaceCharge::AvalancheGridSpaceCharge ( )
inline

Default constructor.

Definition at line 21 of file AvalancheGridSpaceCharge.hh.

21: AvalancheGridSpaceCharge(nullptr) {}

◆ AvalancheGridSpaceCharge() [2/2]

Garfield::AvalancheGridSpaceCharge::AvalancheGridSpaceCharge ( Sensor * sensor)

Constructor.

◆ ~AvalancheGridSpaceCharge()

Garfield::AvalancheGridSpaceCharge::~AvalancheGridSpaceCharge ( )
default

Destructor.

Member Function Documentation

◆ AddElectron()

void Garfield::AvalancheGridSpaceCharge::AddElectron ( double x,
double y,
double z,
double t = 0,
int n = 1 )

Set n electrons onto the grid.

◆ AddElectrons()

void Garfield::AvalancheGridSpaceCharge::AddElectrons ( AvalancheMicroscopic * avmc)

Import electron (no ions) data from AvalancheMicroscopic class to axi-symmetric grid.

◆ AddExtraElectron()

void Garfield::AvalancheGridSpaceCharge::AddExtraElectron ( double y,
int n = 1 )

After calling AddElectron, add more electrons on the same transversal line (y-freedom).

◆ AddFieldFromChargeAt() [1/2]

bool Garfield::AvalancheGridSpaceCharge::AddFieldFromChargeAt ( int iz,
int ir,
double zf,
double rf,
double N,
double & eFieldZ,
double & eFieldR )
private

◆ AddFieldFromChargeAt() [2/2]

bool Garfield::AvalancheGridSpaceCharge::AddFieldFromChargeAt ( int iz,
int ir,
int fz,
int fr,
double N,
double & eFieldZ,
double & eFieldR )
private

◆ DiffuseTimeStep()

void Garfield::AvalancheGridSpaceCharge::DiffuseTimeStep ( double dx,
long nElectron,
double nPosIon,
double nNegIon,
int iz,
int ir,
int gasGap )
private

◆ DistributeCharges()

void Garfield::AvalancheGridSpaceCharge::DistributeCharges ( long nElectron,
double nPosIon,
double nNegIon,
int iz,
int ir,
double stepZ,
double stepR,
int gasGap )
private

◆ EnableAdaptiveTimeStepping()

void Garfield::AvalancheGridSpaceCharge::EnableAdaptiveTimeStepping ( const bool option = true)
inline

Enable adaptive time stepping (default on)

Definition at line 48 of file AvalancheGridSpaceCharge.hh.

◆ EnableDebugging()

void Garfield::AvalancheGridSpaceCharge::EnableDebugging ( const bool option = true)
inline

Enable/disable debugging (log messages) (default off)

Definition at line 31 of file AvalancheGridSpaceCharge.hh.

◆ EnableDiffusion()

void Garfield::AvalancheGridSpaceCharge::EnableDiffusion ( const bool option = true)
inline

Enable diffusion (default off)

Definition at line 40 of file AvalancheGridSpaceCharge.hh.

◆ EnableMC()

void Garfield::AvalancheGridSpaceCharge::EnableMC ( const bool option = true)
inline

Enable Monte Carlo (gain/diffusion) up to a certain total number of electrons e.g.

(1e5) (default on) if disabled: Mean values are considered instead of the statistical processes

Definition at line 56 of file AvalancheGridSpaceCharge.hh.

◆ EnableSpaceChargeEffect()

void Garfield::AvalancheGridSpaceCharge::EnableSpaceChargeEffect ( const bool option = true)
inline

Enable space charge calculations (default on)

Definition at line 43 of file AvalancheGridSpaceCharge.hh.

◆ EnableStickyAnode()

void Garfield::AvalancheGridSpaceCharge::EnableStickyAnode ( const bool option = true)
inline

Enable sticky anode (default on)

Definition at line 34 of file AvalancheGridSpaceCharge.hh.

◆ EnableTOF()

void Garfield::AvalancheGridSpaceCharge::EnableTOF ( const bool option = true)
inline

Enable to use TOF swarm parameters (default on)

Definition at line 37 of file AvalancheGridSpaceCharge.hh.

37{ m_bUseTOF = option; }
bool m_bUseTOF
Flag if TOF parameters should be used, else Magboltz drift and SST spatial coefficients.

◆ ExportGrid()

void Garfield::AvalancheGridSpaceCharge::ExportGrid ( const std::string & filename)

Export the current grid to a txt file (electron, ions numbers and field magnitude) Take care: only where electrons are located is the space-charge field evaluated!

◆ GetAvalancheSize()

long Garfield::AvalancheGridSpaceCharge::GetAvalancheSize ( ) const
inline

Returns the total positive charge in the gap's.

Definition at line 109 of file AvalancheGridSpaceCharge.hh.

109{ return m_nTotPosIons; }
long m_nTotPosIons
total amount of charge created

◆ GetElectronEvolution()

const std::vector< std::pair< double, long > > Garfield::AvalancheGridSpaceCharge::GetElectronEvolution ( ) const
inlinenodiscard

Returns the total electron number evolution.

Definition at line 124 of file AvalancheGridSpaceCharge.hh.

124 {
126 }
std::vector< std::pair< double, long > > m_vNElectronEvolution

◆ GetEllipticIntegrals()

void Garfield::AvalancheGridSpaceCharge::GetEllipticIntegrals ( double x,
double & K,
double & E )
private

◆ GetFreeChargedRing() [1/2]

void Garfield::AvalancheGridSpaceCharge::GetFreeChargedRing ( double zi,
double ri,
double zf,
double rf,
double & eFieldZ,
double & eFieldR )
private

◆ GetFreeChargedRing() [2/2]

void Garfield::AvalancheGridSpaceCharge::GetFreeChargedRing ( int iz,
int ir,
int fz,
int fr,
double & eFieldZ,
double & eFieldR )
private

◆ GetGasGapNumber()

int Garfield::AvalancheGridSpaceCharge::GetGasGapNumber ( int layerIndex)
inlineprivate

Definition at line 245 of file AvalancheGridSpaceCharge.hh.

245 {
246 auto it =
247 std::find(m_vIndexGasGaps.begin(), m_vIndexGasGaps.end(), layerIndex);
248 return (it != m_vIndexGasGaps.end())
249 ? std::distance(m_vIndexGasGaps.begin(), it)
250 : -1;
251 }
std::vector< int > m_vIndexGasGaps
Which layer indices are gas layers.

◆ GetGlobalCoordinates()

void Garfield::AvalancheGridSpaceCharge::GetGlobalCoordinates ( double r,
double z,
double phi,
double & xg,
double & yg,
double & zg,
int gasGap )
private

◆ GetLocalField()

void Garfield::AvalancheGridSpaceCharge::GetLocalField ( int iz,
int ir,
double & eFieldZ,
double & eFieldR,
const std::string & fieldOption,
int gasGap )
private

◆ GetMeanDistance()

double Garfield::AvalancheGridSpaceCharge::GetMeanDistance ( )

Return current mean distance of the electrons on the grid.

◆ GetSwarmParameters()

void Garfield::AvalancheGridSpaceCharge::GetSwarmParameters ( double MagEField,
double & alpha,
double & eta,
double & drift,
double & dSigmaL,
double & dSigmaT,
double & wv,
double & wr,
double & alphaPT,
double & etaPT,
int gasGap )
private

◆ ImportEllipticIntegralValues()

void Garfield::AvalancheGridSpaceCharge::ImportEllipticIntegralValues ( const std::string & filename)
private

◆ Prepare2dMesh()

void Garfield::AvalancheGridSpaceCharge::Prepare2dMesh ( )
private

◆ PrepareElectronsFromMicroscopicAvalanche()

void Garfield::AvalancheGridSpaceCharge::PrepareElectronsFromMicroscopicAvalanche ( )
private

◆ ReachedKPercent()

long Garfield::AvalancheGridSpaceCharge::ReachedKPercent ( ) const
inlinenodiscard

Returns if 100 * K % background charge has been reached.

Definition at line 115 of file AvalancheGridSpaceCharge.hh.

◆ Reset()

void Garfield::AvalancheGridSpaceCharge::Reset ( )

Reset the charges.

◆ Set2dGrid()

void Garfield::AvalancheGridSpaceCharge::Set2dGrid ( double zmin,
double zmax,
int zsteps,
double rmax,
int rsteps )
Parameters
zmincoordinate along direction of background field
zmaxshould more or less be at the end of the gap
zstepsnumber of steps
rmaxcoordinate perpendicular to background field
rstepsnumber of steps

◆ SetFieldCalculation()

void Garfield::AvalancheGridSpaceCharge::SetFieldCalculation ( const std::string & option = "coulomb",
const int nof_approx = 1. )
inline

Sets the different options to calculate the space charge field coulomb: free field approximation (relaxation: fdm-relaxation methode (with a condition to stop)) mirror: symmetric 3 layer single gap rpc with metal - resistive layer - gas gap - r.

l. - m.

Definition at line 67 of file AvalancheGridSpaceCharge.hh.

68 {
69 m_sFieldOption = std::move(option);
70 m_iFieldApprox = std::move(nof_approx);
71 }

◆ SetK()

void Garfield::AvalancheGridSpaceCharge::SetK ( float option = 0.95)
inline

Set the streamer-inception criterion constant K in the interval (0, ∞) s.t.

1 = 100%

Definition at line 75 of file AvalancheGridSpaceCharge.hh.

◆ SetNCrit()

void Garfield::AvalancheGridSpaceCharge::SetNCrit ( const long NCrit = 1e8)
inline

Set electron saturation value applied for each gas gap if space charge effect is turned off (default 1e8)

Definition at line 60 of file AvalancheGridSpaceCharge.hh.

◆ SetSensor()

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

Set the sensor (and determine if it includes a parallel-plate component).

◆ SetStopAtK()

void Garfield::AvalancheGridSpaceCharge::SetStopAtK ( bool option = true)
inline

Stop the avalanche if K % field is reached.

Definition at line 78 of file AvalancheGridSpaceCharge.hh.

◆ SnapTo2dGrid()

bool Garfield::AvalancheGridSpaceCharge::SnapTo2dGrid ( double x,
double y,
double z,
long n = 1,
int gasLayer = 0 )
private

◆ StartGridAvalanche()

void Garfield::AvalancheGridSpaceCharge::StartGridAvalanche ( double dtime = -1)

Starts the simulation with the imported electrons for a time step dt (dt = -1 until there are no electrons left in the gap).

◆ TransportTimeStep()

bool Garfield::AvalancheGridSpaceCharge::TransportTimeStep ( )
private

Member Data Documentation

◆ elliptic_size

const constexpr std::size_t Garfield::AvalancheGridSpaceCharge::elliptic_size {29981}
staticconstexprprivate

Definition at line 340 of file AvalancheGridSpaceCharge.hh.

340{29981};

◆ m_bAdaptiveTime

bool Garfield::AvalancheGridSpaceCharge::m_bAdaptiveTime = true
private

Definition at line 272 of file AvalancheGridSpaceCharge.hh.

◆ m_bDebug

bool Garfield::AvalancheGridSpaceCharge::m_bDebug = false
private

Definition at line 256 of file AvalancheGridSpaceCharge.hh.

◆ m_bDiffusion

bool Garfield::AvalancheGridSpaceCharge::m_bDiffusion = false
private

Definition at line 257 of file AvalancheGridSpaceCharge.hh.

◆ m_bDriftAvalanche

bool Garfield::AvalancheGridSpaceCharge::m_bDriftAvalanche = false
private

Definition at line 260 of file AvalancheGridSpaceCharge.hh.

◆ m_bFieldK

bool Garfield::AvalancheGridSpaceCharge::m_bFieldK = false
private

Definition at line 269 of file AvalancheGridSpaceCharge.hh.

◆ m_bImportAvalanche

bool Garfield::AvalancheGridSpaceCharge::m_bImportAvalanche = false
private

Definition at line 262 of file AvalancheGridSpaceCharge.hh.

◆ m_bImportElliptic

bool Garfield::AvalancheGridSpaceCharge::m_bImportElliptic = false
private

Definition at line 273 of file AvalancheGridSpaceCharge.hh.

◆ m_bMC

bool Garfield::AvalancheGridSpaceCharge::m_bMC = true
private

Definition at line 282 of file AvalancheGridSpaceCharge.hh.

◆ m_bPreparedImportAvalanche

bool Garfield::AvalancheGridSpaceCharge::m_bPreparedImportAvalanche = false
private

Definition at line 263 of file AvalancheGridSpaceCharge.hh.

◆ m_bRatesAvailable

bool Garfield::AvalancheGridSpaceCharge::m_bRatesAvailable = true
private

Flag if temporal rates are available to the simulation.

Definition at line 280 of file AvalancheGridSpaceCharge.hh.

◆ m_bSpaceCharge

bool Garfield::AvalancheGridSpaceCharge::m_bSpaceCharge = true
private

Definition at line 265 of file AvalancheGridSpaceCharge.hh.

◆ m_bStick

bool Garfield::AvalancheGridSpaceCharge::m_bStick = true
private

Definition at line 258 of file AvalancheGridSpaceCharge.hh.

◆ m_bStopAtK

bool Garfield::AvalancheGridSpaceCharge::m_bStopAtK = false
private

Definition at line 268 of file AvalancheGridSpaceCharge.hh.

◆ m_bUseTOF

bool Garfield::AvalancheGridSpaceCharge::m_bUseTOF = true
private

Flag if TOF parameters should be used, else Magboltz drift and SST spatial coefficients.

Definition at line 276 of file AvalancheGridSpaceCharge.hh.

◆ m_bWrAvailable

bool Garfield::AvalancheGridSpaceCharge::m_bWrAvailable = true
private

Flag if bulk drift velocity is available to the simulation.

Definition at line 278 of file AvalancheGridSpaceCharge.hh.

◆ m_className

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

Definition at line 254 of file AvalancheGridSpaceCharge.hh.

◆ m_dMinGroups

double Garfield::AvalancheGridSpaceCharge::m_dMinGroups = 50
private

Definition at line 285 of file AvalancheGridSpaceCharge.hh.

◆ m_dt

double Garfield::AvalancheGridSpaceCharge::m_dt = 0.
private

Time step.

Definition at line 304 of file AvalancheGridSpaceCharge.hh.

◆ m_elliptic

const std::array< std::array< double, 3 >, AvalancheGridSpaceCharge::elliptic_size > Garfield::AvalancheGridSpaceCharge::m_elliptic
staticprivate

Definition at line 341 of file AvalancheGridSpaceCharge.hh.

◆ m_ezBkg

std::vector<double> Garfield::AvalancheGridSpaceCharge::m_ezBkg = {0}
private

Uniform background field in z direction, can be negative.

Definition at line 329 of file AvalancheGridSpaceCharge.hh.

329{0};

◆ m_fStreamerK

float Garfield::AvalancheGridSpaceCharge::m_fStreamerK = 0.95
private

Definition at line 267 of file AvalancheGridSpaceCharge.hh.

◆ m_grid

std::vector<std::vector<GridNode> > Garfield::AvalancheGridSpaceCharge::m_grid
private

grid with nodes on it

Definition at line 311 of file AvalancheGridSpaceCharge.hh.

◆ m_iFieldApprox

int Garfield::AvalancheGridSpaceCharge::m_iFieldApprox = 1
private

Definition at line 284 of file AvalancheGridSpaceCharge.hh.

◆ m_isgridset

bool Garfield::AvalancheGridSpaceCharge::m_isgridset = false
private

Distance between the grid points.

Keeps track if the grid has been defined.

Definition at line 298 of file AvalancheGridSpaceCharge.hh.

◆ m_lElectronsK

long Garfield::AvalancheGridSpaceCharge::m_lElectronsK
private

Definition at line 270 of file AvalancheGridSpaceCharge.hh.

◆ m_lNCrit

long Garfield::AvalancheGridSpaceCharge::m_lNCrit = 1e8
private

Definition at line 264 of file AvalancheGridSpaceCharge.hh.

◆ m_nTotElectron

long Garfield::AvalancheGridSpaceCharge::m_nTotElectron = 0
private

Total amount of electrons at time step.

Definition at line 299 of file AvalancheGridSpaceCharge.hh.

◆ m_nTotPosIons

long Garfield::AvalancheGridSpaceCharge::m_nTotPosIons = 0
private

total amount of charge created

Definition at line 300 of file AvalancheGridSpaceCharge.hh.

◆ m_pp

ComponentParallelPlate* Garfield::AvalancheGridSpaceCharge::m_pp = nullptr
private

Definition at line 287 of file AvalancheGridSpaceCharge.hh.

◆ m_rGrid

std::vector<double> Garfield::AvalancheGridSpaceCharge::m_rGrid
private

Distance between the grid points.

Grid points of r-coordinate.

Definition at line 294 of file AvalancheGridSpaceCharge.hh.

◆ m_rSteps

int Garfield::AvalancheGridSpaceCharge::m_rSteps = 0.
private

Number of grid points.

Definition at line 295 of file AvalancheGridSpaceCharge.hh.

◆ m_rStepSize

double Garfield::AvalancheGridSpaceCharge::m_rStepSize = 0.
private

Definition at line 296 of file AvalancheGridSpaceCharge.hh.

◆ m_run

bool Garfield::AvalancheGridSpaceCharge::m_run = true
private

Tracking if the charges are still in the drift gap.

Definition at line 306 of file AvalancheGridSpaceCharge.hh.

◆ m_sensor

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

Definition at line 288 of file AvalancheGridSpaceCharge.hh.

◆ m_sFieldOption

std::string Garfield::AvalancheGridSpaceCharge::m_sFieldOption = "coulomb"
private

Definition at line 333 of file AvalancheGridSpaceCharge.hh.

◆ m_time

double Garfield::AvalancheGridSpaceCharge::m_time = 0.
private

Clock.

Definition at line 302 of file AvalancheGridSpaceCharge.hh.

◆ m_time0

double Garfield::AvalancheGridSpaceCharge::m_time0 = 0.
private

Initial time.

Definition at line 303 of file AvalancheGridSpaceCharge.hh.

◆ m_vCoNGasLayer

std::vector<std::vector<double> > Garfield::AvalancheGridSpaceCharge::m_vCoNGasLayer {}
private

Coordinates of center of electron number.

Required: y in [zmin, zmax]

Definition at line 324 of file AvalancheGridSpaceCharge.hh.

324{};

◆ m_vElectrons

std::vector<std::vector<Point> > Garfield::AvalancheGridSpaceCharge::m_vElectrons
private

Electrons to transfer onto grid.

Definition at line 313 of file AvalancheGridSpaceCharge.hh.

◆ m_vGroupSizes

std::vector<long> Garfield::AvalancheGridSpaceCharge::m_vGroupSizes
private
Initial value:
= {
1500, 800, 400, 200, 100,
50, 20, 10, 5, 2}

same values as Lippmann & Riegler

Definition at line 316 of file AvalancheGridSpaceCharge.hh.

316 {
317 1500, 800, 400, 200, 100,
318 50, 20, 10, 5, 2};

◆ m_vIndexGasGaps

std::vector<int> Garfield::AvalancheGridSpaceCharge::m_vIndexGasGaps = {0}
private

Which layer indices are gas layers.

Definition at line 320 of file AvalancheGridSpaceCharge.hh.

320{0};

◆ m_vNElectronEvolution

std::vector<std::pair<double, long> > Garfield::AvalancheGridSpaceCharge::m_vNElectronEvolution
private

Definition at line 315 of file AvalancheGridSpaceCharge.hh.

◆ m_vSaturatedGaps

std::vector<int> Garfield::AvalancheGridSpaceCharge::m_vSaturatedGaps {}
private

Which gas gaps are saturated if saturation is on.

Definition at line 331 of file AvalancheGridSpaceCharge.hh.

331{};

◆ m_vYPointInGasGap

std::vector<double> Garfield::AvalancheGridSpaceCharge::m_vYPointInGasGap {}
private

Example point (y-coord) in each gas gap.

Definition at line 326 of file AvalancheGridSpaceCharge.hh.

326{};

◆ m_zGasGapBoundaries

std::vector<std::vector<int> > Garfield::AvalancheGridSpaceCharge::m_zGasGapBoundaries
private

[k] -> {izLeft, ..., izRight}

Definition at line 309 of file AvalancheGridSpaceCharge.hh.

◆ m_zGrid

std::vector<double> Garfield::AvalancheGridSpaceCharge::m_zGrid
private

Grid points of z-coordinate.

Definition at line 290 of file AvalancheGridSpaceCharge.hh.

◆ m_zSteps

int Garfield::AvalancheGridSpaceCharge::m_zSteps = 0
private

Number of grid points.

Definition at line 291 of file AvalancheGridSpaceCharge.hh.

◆ m_zStepSize

double Garfield::AvalancheGridSpaceCharge::m_zStepSize = 0.
private

Definition at line 292 of file AvalancheGridSpaceCharge.hh.


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