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

#include <AvalancheGrid.hh>

Public Member Functions

 AvalancheGrid ()
 Constructor.
 
 ~AvalancheGrid ()
 Destructor.
 
void SetSensor (Sensor *sensor)
 Set the sensor.
 
void StartGridAvalanche ()
 
void SetElectronVelocity (const double vx, const double vy, const double vz)
 Set the electron drift velocity (in cm / ns).
 
void SetElectronTownsend (const double town)
 Set the electron Townsend coefficient (in 1 / cm).
 
void SetElectronAttachment (const double att)
 Set the electron attachment coefficient (in 1 / cm).
 
void SetMaxAvalancheSize (const double size)
 Set the maximum avalanche size (1e7 by default).
 
void EnableDiffusion (const double diffSigma)
 
void AvalancheElectron (const double x, const double y, const double z, const double t=0, const int n=1)
 
void ImportElectronsFromAvalancheMicroscopic (AvalancheMicroscopic *avmc)
 Import electron data from AvalancheMicroscopic class.
 
void SetGrid (const double xmin, const double xmax, const int xsteps, const double ymin, const double ymax, const int ysteps, const double zmin, const double zmax, const int zsteps)
 Import electron data from AvalancheMicroscopic class.
 
int GetAmountOfStartingElectrons ()
 Returns the initial number of electrons in the avalanche.
 
int GetAvalancheSize ()
 Returns the final number of electrons in the avalanche.
 
void AsignLayerIndex (ComponentParallelPlate *RPC)
 Asigning layer index to all Avalanche nodes.
 
void EnableDebugging ()
 
void Reset ()
 

Detailed Description

Calculate avalanches in a uniform electric field using avalanche statistics.

Definition at line 18 of file AvalancheGrid.hh.

Constructor & Destructor Documentation

◆ AvalancheGrid()

Garfield::AvalancheGrid::AvalancheGrid ( )
inline

Constructor.

Definition at line 21 of file AvalancheGrid.hh.

21{}

◆ ~AvalancheGrid()

Garfield::AvalancheGrid::~AvalancheGrid ( )
inline

Destructor.

Definition at line 23 of file AvalancheGrid.hh.

23{}

Member Function Documentation

◆ AsignLayerIndex()

void Garfield::AvalancheGrid::AsignLayerIndex ( ComponentParallelPlate * RPC)

Asigning layer index to all Avalanche nodes.

Definition at line 509 of file AvalancheGrid.cc.

509 {
510 int im = 0;
511 double epsM = 0;
512 std::vector<double> nLayer(RPC->NumberOfLayers(), 0);
513 for (AvalancheNode &node : m_activeNodes) {
514 double y = m_avgrid.ygrid[node.iy];
515 RPC->getLayer(y, im, epsM);
516 node.layer = im;
517 nLayer[im - 1] += node.n;
518 // std::cerr << m_className << "::AsignLayerIndex::im = "<<im<<".\n";
519 }
520
521 m_NLayer = nLayer;
522 m_layerIndix = true;
523}

◆ AvalancheElectron()

void Garfield::AvalancheGrid::AvalancheElectron ( const double x,
const double y,
const double z,
const double t = 0,
const int n = 1 )

Setting the starting point of an electron that .

Parameters
zz-coordinate of initial electron.
xx-coordinate of initial electron.
vspeed of initial electron.
tstarting time of avalanche.

Definition at line 360 of file AvalancheGrid.cc.

362 {
363 m_driftAvalanche = true;
364
365 if (m_avgrid.time == 0 && m_avgrid.time != t && m_debug)
366 std::cerr << m_className
367 << "::CreateAvalanche::Overwriting start time of avalanche for t "
368 "= 0 to " << t << ".\n";
369 m_avgrid.time = t;
370
371 if (SnapToGrid(m_avgrid, x, y, z, 0, n) && m_debug)
372 std::cerr << m_className
373 << "::CreateAvalanche::Electron added at (t,x,y,z) = (" << t
374 << "," << x << "," << y << "," << z << ").\n";
375}

◆ EnableDebugging()

void Garfield::AvalancheGrid::EnableDebugging ( )
inline

Definition at line 83 of file AvalancheGrid.hh.

83{ m_debug = true; }

◆ EnableDiffusion()

void Garfield::AvalancheGrid::EnableDiffusion ( const double diffSigma)
inline

Enable transverse diffusion of electrons with transverse diffusion coefficients (in √cm).

Definition at line 54 of file AvalancheGrid.hh.

54 {
55 m_diffusion = true;
56 m_DiffSigma = diffSigma;
57 }

◆ GetAmountOfStartingElectrons()

int Garfield::AvalancheGrid::GetAmountOfStartingElectrons ( )
inline

Returns the initial number of electrons in the avalanche.

Definition at line 76 of file AvalancheGrid.hh.

76{ return m_nestart; }

◆ GetAvalancheSize()

int Garfield::AvalancheGrid::GetAvalancheSize ( )
inline

Returns the final number of electrons in the avalanche.

Definition at line 78 of file AvalancheGrid.hh.

78{ return m_avgrid.N; }

◆ ImportElectronsFromAvalancheMicroscopic()

void Garfield::AvalancheGrid::ImportElectronsFromAvalancheMicroscopic ( AvalancheMicroscopic * avmc)

Import electron data from AvalancheMicroscopic class.

Definition at line 377 of file AvalancheGrid.cc.

378 {
379 // Get the information of the electrons from the AvalancheMicroscopic class.
380 if (!avmc) return;
381
382 if (!m_importAvalanche) m_importAvalanche = true;
383
384 // Get initial positions of electrons
385 for (const auto& electron : avmc->GetElectrons()) {
386 // If the electron is not stopped due to
387 // the upper bound of the time range: then skip this electron.
388 if (electron.status != -17) continue;
389 const auto& p1 = electron.path.at(0);
390 const auto& p2 = electron.path.back();
391 const double vel = (p2.z - p1.z) / (p2.t - p1.t);
392 m_avgrid.time = p2.t;
393
394 if (SnapToGrid(m_avgrid, p2.x, p2.y, p2.z, vel) && m_debug) {
395 std::cerr << m_className
396 << "::ImportElectronsFromAvalancheMicroscopic: Electron added at "
397 "(x,y,z) = (" << p2.x << "," << p2.y << "," << p2.z << ").\n";
398 }
399 }
400}

◆ Reset()

void Garfield::AvalancheGrid::Reset ( )

Definition at line 492 of file AvalancheGrid.cc.

492 {
493
494 std::cerr << m_className << "::Reset::Resetting AvalancheGrid.\n";
495 m_avgrid.time = 0;
496 m_avgrid.N = 0;
497 m_avgrid.run = true;
498
499 m_Saturated = false;
500 m_SaturationTime = -1;
501
502 m_driftAvalanche = false;
503
504 m_activeNodes.clear();
505 m_layerIndix = false;
506 m_NLayer.clear();
507}

◆ SetElectronAttachment()

void Garfield::AvalancheGrid::SetElectronAttachment ( const double att)
inline

Set the electron attachment coefficient (in 1 / cm).

Definition at line 49 of file AvalancheGrid.hh.

49{ m_Attachment = att; }

◆ SetElectronTownsend()

void Garfield::AvalancheGrid::SetElectronTownsend ( const double town)
inline

Set the electron Townsend coefficient (in 1 / cm).

Definition at line 47 of file AvalancheGrid.hh.

47{ m_Townsend = town; }

◆ SetElectronVelocity()

void Garfield::AvalancheGrid::SetElectronVelocity ( const double vx,
const double vy,
const double vz )
inline

Set the electron drift velocity (in cm / ns).

Definition at line 36 of file AvalancheGrid.hh.

36 {
37 double vel = sqrt(vx * vx + vy * vy + vz * vz);
38 if (vel != std::abs(vx) && vel != std::abs(vy) && vel != std::abs(vz))
39 return;
40 int nx = (int)vx / vel;
41 int ny = (int)vy / vel;
42 int nz = (int)vz / vel;
43 m_velNormal = {nx, ny, nz};
44 m_Velocity = -std::abs(vel);
45 }
DoubleAc sqrt(const DoubleAc &f)
Definition DoubleAc.cpp:314

◆ SetGrid()

void Garfield::AvalancheGrid::SetGrid ( const double xmin,
const double xmax,
const int xsteps,
const double ymin,
const double ymax,
const int ysteps,
const double zmin,
const double zmax,
const int zsteps )

Import electron data from AvalancheMicroscopic class.

Definition at line 51 of file AvalancheGrid.cc.

55 {
56 m_avgrid.gridset = true;
57
58 if (zmin >= zmax || zsteps <= 0 || xmin > xmax || xsteps <= 0 ||
59 ymin > ymax || ysteps <= 0) {
60 std::cerr << m_className
61 << "::SetGrid:Error: Grid is not properly defined.\n";
62 return;
63 }
64
65 // Setting grid
66
67 SetZGrid(m_avgrid, zmax, zmin, zsteps);
68 SetYGrid(m_avgrid, ymax, ymin, ysteps);
69 SetXGrid(m_avgrid, xmax, xmin, xsteps);
70
71 if (m_debug) {
72 std::cerr << m_className << "::SetGrid: Grid created:\n";
73 std::cerr << " x range = (" << xmin << "," << xmax << ").\n";
74 std::cerr << " y range = (" << ymin << "," << ymax << ").\n";
75 std::cerr << " z range = (" << zmin << "," << zmax << ").\n";
76 }
77}

◆ SetMaxAvalancheSize()

void Garfield::AvalancheGrid::SetMaxAvalancheSize ( const double size)
inline

Set the maximum avalanche size (1e7 by default).

Definition at line 51 of file AvalancheGrid.hh.

51{ m_MaxSize = size; }

◆ SetSensor()

void Garfield::AvalancheGrid::SetSensor ( Sensor * sensor)
inline

Set the sensor.

Definition at line 25 of file AvalancheGrid.hh.

25{ m_sensor = sensor; }

◆ StartGridAvalanche()

void Garfield::AvalancheGrid::StartGridAvalanche ( )

Start grid based avalanche simulation.

Parameters
zmin,zmaxz-coordinate range of grid [cm].
zstepsamount of z-coordinate points in grid.
xmin,xmaxx-coordinate range of grid [cm].
xstepsamount of x-coordinate points in grid.

Definition at line 319 of file AvalancheGrid.cc.

319 {
320 // Start the AvalancheGrid algorithm.
321 if ((!m_importAvalanche && !m_driftAvalanche) || !m_sensor) return;
322
323 std::cerr << m_className
324 << "::StartGridAvalanche::Starting grid based simulation with "
325 << m_avgrid.N << " initial electrons.\n";
326 if (m_avgrid.N <= 0) {
327 std::cerr << m_className << "::StartGridAvalanche::Cancelled.\n";
328 return;
329 }
330
331 m_nestart = m_avgrid.N;
332
333 // Main loop.
334 while (m_avgrid.run == true) {
335 if (m_debug)
336 std::cerr
337 << "============ \n" << m_className
338 << "::StartGridAvalanche: Looping over nodes.\n ============ \n";
339 NextAvalancheGridPoint(m_avgrid);
340 }
341
342 std::vector<double> tlist = {};
343 for (AvalancheNode &node : m_activeNodes) {
344 tlist.push_back(node.time);
345 }
346 double maxTime = *max_element(std::begin(tlist), std::end(tlist));
347
348 if (m_Saturated)
349 std::cerr << m_className
350 << "::StartGridAvalanche::Avalanche maximum size of " << m_MaxSize
351 << " electrons reached at " << m_SaturationTime << " ns.\n";
352
353 std::cerr << m_className
354 << "::StartGridAvalanche::Final avalanche size = " << m_avgrid.N
355 << " ended at t = " << maxTime << " ns.\n";
356
357 return;
358}

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