Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::AvalancheGrid Class Reference

#include <AvalancheGrid.hh>

Public Member Functions

 AvalancheGrid ()
 Constructor.
 
 ~AvalancheGrid ()
 Destructor.
 
void SetSensor (Sensor *sensor)
 Set the sensor.
 
void SetAvalancheMicroscopic (AvalancheMicroscopic *avmc)
 Set the AvalancheMicroscopic.
 
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 CreateAvalanche (const double x, const double y, const double z, const double t=0, const int n=1)
 
void GetElectronsFromAvalancheMicroscopic ()
 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 ()
 
void EnableDebugging ()
 
void Reset ()
 

Detailed Description

Calculate avalanches in a uniform electric field using avalanche statistics.

Definition at line 17 of file AvalancheGrid.hh.

Constructor & Destructor Documentation

◆ AvalancheGrid()

Garfield::AvalancheGrid::AvalancheGrid ( )
inline

Constructor.

Definition at line 20 of file AvalancheGrid.hh.

20{}

◆ ~AvalancheGrid()

Garfield::AvalancheGrid::~AvalancheGrid ( )
inline

Destructor.

Definition at line 22 of file AvalancheGrid.hh.

22{}

Member Function Documentation

◆ CreateAvalanche()

void Garfield::AvalancheGrid::CreateAvalanche ( 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 483 of file AvalancheGrid.cc.

485 {
486 m_driftAvalanche = true;
487
488 if (m_avgrid.time == 0 && m_avgrid.time != t && m_debug)
489 std::cerr << m_className
490 << "::CreateAvalanche::Overwriting start time of avalanche for t "
491 "= 0 to "
492 << t << ".\n";
493 m_avgrid.time = t;
494
495 if (SnapToGrid(m_avgrid, x, y, z, 0, n) && m_debug)
496 std::cerr << m_className
497 << "::CreateAvalanche::Electron added at (t,x,y,z) = (" << t
498 << "," << x << "," << y << "," << z << ").\n";
499
500 // std::cerr<< m_className<< "::CreateAvalanche::expected contribution is "<<
501 // exp((m_Townsend-m_Attachment)*z) << ".\n";
502}

◆ EnableDebugging()

void Garfield::AvalancheGrid::EnableDebugging ( )
inline

Definition at line 76 of file AvalancheGrid.hh.

76{ 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

Definition at line 74 of file AvalancheGrid.hh.

74{ return m_nestart; }

◆ GetElectronsFromAvalancheMicroscopic()

void Garfield::AvalancheGrid::GetElectronsFromAvalancheMicroscopic ( )

Import electron data from AvalancheMicroscopic class.

Definition at line 503 of file AvalancheGrid.cc.

503 {
504 // Get the information of the electrons from the AvalancheMicroscopic class.
505 if (!m_avmc) return;
506
507 if (!m_importAvalanche) m_importAvalanche = true;
508
509 int np = m_avmc->GetNumberOfElectronEndpoints();
510
511 if (np == 0) return;
512
513 // Get initial positions of electrons
514 double x1, y1, z1, t1;
515 double x2, y2, z2, t2;
516 double e1, e2;
517 int status;
518
519 double vel = 0.;
520
521 for (int i = 0; i < np; ++i) {
522 m_avmc->GetElectronEndpoint(i, x1, y1, z1, t1, e1, x2, y2, z2, t2, e2,
523 status);
524
525 vel = (z2 - z1) / (t2 - t1);
526
527 m_avgrid.time = t2;
528
529 if (SnapToGrid(m_avgrid, x2, y2, z2, vel) && m_debug)
530 std::cerr << m_className
531 << "::GetElectronsFromAvalancheMicroscopic::Electron added at "
532 "(x,y,z) = ("
533 << x2 << "," << y2 << "," << z2 << ").\n";
534 }
535}
void GetElectronEndpoint(const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, int &status) const
unsigned int GetNumberOfElectronEndpoints() const

Referenced by StartGridAvalanche().

◆ Reset()

void Garfield::AvalancheGrid::Reset ( )

Definition at line 588 of file AvalancheGrid.cc.

588 {
589 std::cerr << m_className << "::Reset::Resetting AvalancheGrid.\n";
590
591 m_avgrid.n.clear();
592 m_avgrid.transverseDiffusion.clear();
593 m_avgrid.time = 0;
594 m_avgrid.N = 0;
595 m_avgrid.run = true;
596
597 m_Saturated = false;
598 m_SaturationTime = -1;
599
600 m_driftAvalanche = false;
601
602 std::vector<int> nhx(m_avgrid.xsteps, 0);
603 std::vector<std::vector<int>> nhy(m_avgrid.ysteps, nhx);
604 std::vector<std::vector<std::vector<int>>> nhz(m_avgrid.zsteps, nhy);
605 m_avgrid.n = nhz;
606
607 for (int i = 0; i < 3; i++) {
608 m_avgrid.gridPosition[i].clear();
609 }
610
611 // Get the diffusion factors for neighboring points on the grid.
612 DiffusionFactors(m_avgrid);
613}

◆ SetAvalancheMicroscopic()

void Garfield::AvalancheGrid::SetAvalancheMicroscopic ( AvalancheMicroscopic avmc)
inline

Set the AvalancheMicroscopic.

Definition at line 26 of file AvalancheGrid.hh.

26{ m_avmc = avmc; }

◆ 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 37 of file AvalancheGrid.hh.

37 {
38 double vel = sqrt(vx * vx + vy * vy + vz * vz);
39 if (vel != std::abs(vx) && vel != std::abs(vy) && vel != std::abs(vz)) 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 58 of file AvalancheGrid.cc.

62 {
63 m_avgrid.gridset = true;
64
65 if (zmin >= zmax || zsteps <= 0 || xmin > xmax || xsteps <= 0 ||
66 ymin > ymax || ysteps <= 0) {
67 std::cerr << m_className
68 << "::SetGrid:Error: Grid is not properly defined.\n";
69 return;
70 }
71
72 // Setting grid
73
74 SetZGrid(m_avgrid, zmax, zmin, zsteps);
75 SetYGrid(m_avgrid, ymax, ymin, ysteps);
76 SetXGrid(m_avgrid, xmax, xmin, xsteps);
77
78 if (m_sensor) GetParametersFromSensor();
79}

◆ 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 24 of file AvalancheGrid.hh.

24{ 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 381 of file AvalancheGrid.cc.

381 {
382 // Start the AvalancheGrid algorithm.
383 if ((!m_avmc && !m_driftAvalanche) || !m_sensor) return;
384
385 if (!m_importAvalanche && m_avmc) GetElectronsFromAvalancheMicroscopic();
386
387 std::cerr << m_className
388 << "::StartGridAvalanche::Starting grid based simulation with "
389 << m_avgrid.N << " initial electrons.\n";
390 if (m_avgrid.N <= 0) {
391 std::cerr << m_className << "::StartGridAvalanche::Cancelled.\n";
392 return;
393 }
394
395 m_nestart = m_avgrid.N;
396
397 // The vector containing the indexes of the z-coordinate grid points of the
398 // initial electrons needs to be ordered from small to large values. All
399 // duplicate values need to be removed.
400
401 GetParametersFromSensor();
402
403 SortPositionVector();
404
405 if (m_debug) {
406 std::cerr << m_className
407 << "::StartGridAvalanche::m_avgrid.gridPosition at iz = ";
408 for (size_t i = 0; i < m_avgrid.gridPosition[0].size(); i++) {
409 std::cerr << m_avgrid.gridPosition[0][i] << ",";
410 }
411 std::cerr << ".\n";
412 }
413
414 // Set velocity if given.
415 m_avgrid.velocity = m_Velocity;
416 // Main loop.
417 while (m_avgrid.run == true) {
418 NextAvalancheGridPoint(m_avgrid);
419 }
420
421 if (m_Saturated)
422 std::cerr << m_className
423 << "::StartGridAvalanche::Avalanche maximum size of " << m_MaxSize
424 << " electrons reached at " << m_SaturationTime << " ns.\n";
425
426 std::cerr << m_className
427 << "::StartGridAvalanche::Final avalanche size = " << m_avgrid.N
428 << " at t = " << m_avgrid.time << " ns.\n";
429
430 return;
431}
void GetElectronsFromAvalancheMicroscopic()
Import electron data from AvalancheMicroscopic class.

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