Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
AvalancheGrid.hh
Go to the documentation of this file.
1#ifndef G_AVALANCHE_GRID_H
2#define G_AVALANCHE_GRID_H
3
4#include <iostream>
5#include <string>
6#include <vector>
7
10#include "GarfieldConstants.hh"
11#include "Sensor.hh"
12
13namespace Garfield {
14
15/// Calculate avalanches in a uniform electric field using avalanche
16/// statistics.
17
19 public:
20 /// Constructor
22 /// Destructor
24 /// Set the sensor.
25 void SetSensor(Sensor *sensor) { m_sensor = sensor; }
26
27 /** Start grid based avalanche simulation.
28 *
29 * \param zmin,zmax z-coordinate range of grid [cm].
30 * \param zsteps amount of z-coordinate points in grid.
31 * \param xmin,xmax x-coordinate range of grid [cm].
32 * \param xsteps amount of x-coordinate points in grid.
33 */
34 void StartGridAvalanche();
35 /// Set the electron drift velocity (in cm / ns).
36 void SetElectronVelocity(const double vx, const double vy, const double vz) {
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 }
46 /// Set the electron Townsend coefficient (in 1 / cm).
47 void SetElectronTownsend(const double town) { m_Townsend = town; }
48 /// Set the electron attachment coefficient (in 1 / cm).
49 void SetElectronAttachment(const double att) { m_Attachment = att; }
50 /// Set the maximum avalanche size (1e7 by default).
51 void SetMaxAvalancheSize(const double size) { m_MaxSize = size; }
52 /// Enable transverse diffusion of electrons with transverse diffusion
53 /// coefficients (in √cm).
54 void EnableDiffusion(const double diffSigma) {
55 m_diffusion = true;
56 m_DiffSigma = diffSigma;
57 }
58 /** Setting the starting point of an electron that .
59 *
60 * \param z z-coordinate of initial electron.
61 * \param x x-coordinate of initial electron.
62 * \param v speed of initial electron.
63 * \param t starting time of avalanche.
64 */
65 void AvalancheElectron(const double x, const double y, const double z,
66 const double t = 0, const int n = 1);
67 /// Import electron data from AvalancheMicroscopic class
69
70 /// Import electron data from AvalancheMicroscopic class
71 void SetGrid(const double xmin, const double xmax, const int xsteps,
72 const double ymin, const double ymax, const int ysteps,
73 const double zmin, const double zmax, const int zsteps);
74
75 /// Returns the initial number of electrons in the avalanche.
76 int GetAmountOfStartingElectrons() { return m_nestart; }
77 /// Returns the final number of electrons in the avalanche.
78 int GetAvalancheSize() { return m_avgrid.N; }
79
80 /// Asigning layer index to all Avalanche nodes.
82
83 void EnableDebugging() { m_debug = true; }
84
85 void Reset();
86
87 private:
88 bool m_debug = false;
89
90 double m_Townsend = -1; // [1/cm]
91
92 double m_Attachment = -1; // [1/cm]
93
94 double m_Velocity = 0.; // [cm/ns]
95
96 std::vector<int> m_velNormal = {0, 0, 0};
97
98 double m_MaxSize = 1.6e7; // Saturations size
99
100 bool m_Saturated = false; // Check if avalanche has reached maximum size
101
102 double m_SaturationTime =
103 -1.; // Time when the avalanche has reached maximum size
104
105 bool m_diffusion = false; // Check if transverse diffusion is enabled.
106
107 double m_DiffSigma = 0.; // Transverse diffusion coefficients (in √cm).
108
109 int m_nestart = 0.;
110
111 bool m_driftAvalanche = false;
112 bool m_importAvalanche = false;
113
114 bool m_layerIndix = false;
115 std::vector<double> m_NLayer;
116
117 std::string m_className = "AvalancheGrid";
118
119 Sensor *m_sensor = nullptr;
120
121 bool m_printPar = false;
122
123 struct Grid {
124 std::vector<double> zgrid; ///< Grid points of z-coordinate.
125 int zsteps = 0; ///< Amount of grid points.
126 double zStepSize =
127 0.; ///< Distance between the grid points of z-coordinate.
128
129 std::vector<double> ygrid; ///< Grid points of y-coordinate.
130 double yStepSize = 0.; ///< Amount of grid points.
131 int ysteps = 0.; ///< Distance between the grid points of y-coordinate.
132
133 std::vector<double> xgrid; ///< Grid points of x-coordinate.
134 double xStepSize = 0.; ///< Amount of grid points.
135 int xsteps = 0.; ///< Distance between the grid points of x-coordinate.
136
137 bool gridset = false; ///< Keeps track if the grid has been defined.
138 int N = 0; ///< Total amount of charge.
139
140 double time = 0; ///< Clock.
141
142 bool run = true; ///< Tracking if the charges are still in the drift gap.
143 };
144
145 struct AvalancheNode {
146 double ix = 0;
147 double iy = 0;
148 double iz = 0;
149
150 int n = 1;
151
152 int layer = 0;
153
154 double townsend = 0;
155 double attachment = 0;
156 double velocity = 0;
157
158 double stepSize = 0;
159 std::vector<int> velNormal = {0, 0, 0};
160
161 double time = 0.; ///< Clock.
162 double dt = -1.; ///< time step.
163
164 bool active = true;
165 double dSigmaL = 0;
166 double dSigmaT = 0;
167 };
168
169 std::vector<AvalancheNode> m_activeNodes = {};
170
171 Grid m_avgrid;
172 // Setting z-coordinate grid.
173 void SetZGrid(Grid &av, const double top, const double bottom,
174 const int steps);
175 // Setting y-coordinate grid.
176 void SetYGrid(Grid &av, const double top, const double bottom,
177 const int steps);
178 // Setting x-coordinate grid.
179 void SetXGrid(Grid &av, const double top, const double bottom,
180 const int steps);
181 // Get size of avalanche when going from z to z-dz.
182 int GetAvalancheSize(double dz, const int nsize, const double alpha,
183 const double eta);
184 // Assign electron to the closest grid point.
185 bool SnapToGrid(Grid &av, const double x, const double y, const double z,
186 const double v, const int n = 1);
187 // Go to next time step.
188 void NextAvalancheGridPoint(Grid &av);
189 // Obtain the Townsend coef., Attachment coef. and velocity vector from
190 // sensor class.
191 bool GetParameters(AvalancheNode &newNode);
192
193 void DeactivateNode(AvalancheNode &node);
194};
195} // namespace Garfield
196
197#endif
void SetMaxAvalancheSize(const double size)
Set the maximum avalanche size (1e7 by default).
void AvalancheElectron(const double x, const double y, const double z, const double t=0, const int n=1)
void SetSensor(Sensor *sensor)
Set the sensor.
void SetElectronTownsend(const double town)
Set the electron Townsend coefficient (in 1 / cm).
void ImportElectronsFromAvalancheMicroscopic(AvalancheMicroscopic *avmc)
Import electron data from AvalancheMicroscopic class.
int GetAmountOfStartingElectrons()
Returns the initial number of electrons in the avalanche.
void SetElectronAttachment(const double att)
Set the electron attachment coefficient (in 1 / cm).
AvalancheGrid()
Constructor.
void SetElectronVelocity(const double vx, const double vy, const double vz)
Set the electron drift velocity (in cm / ns).
void AsignLayerIndex(ComponentParallelPlate *RPC)
Asigning layer index to all Avalanche nodes.
void EnableDiffusion(const double diffSigma)
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 GetAvalancheSize()
Returns the final number of electrons in the avalanche.
Calculate electron drift lines and avalanches using microscopic tracking.
Component for parallel-plate geometries.