Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
AvalancheGridSpaceCharge.hh
Go to the documentation of this file.
1#ifndef GARFIELD_AVALANCHEGRIDSPACECHARGE_HH
2#define GARFIELD_AVALANCHEGRIDSPACECHARGE_HH
3
4#include <algorithm>
5#include <string>
6#include <utility>
7#include <vector>
8
9namespace Garfield {
10class Sensor;
13
18 public:
22 AvalancheGridSpaceCharge(Sensor *sensor);
25
27 void Reset();
28
30 void EnableDebugging(const bool option = true) { m_bDebug = option; }
31
33 void EnableStickyAnode(const bool option = true) { m_bStick = option; }
34
36 void EnableTOF(const bool option = true) { m_bUseTOF = option; }
37
39 void EnableDiffusion(const bool option = true) { m_bDiffusion = option; }
40
42 void EnableSpaceChargeEffect(const bool option = true) {
43 m_bSpaceCharge = option;
44 }
45
47 void EnableAdaptiveTimeStepping(const bool option = true) {
48 m_bAdaptiveTime = option;
49 }
50
55 void EnableMC(const bool option = true) { m_bMC = option; }
56
59 void SetNCrit(const long NCrit = 1e8) { m_lNCrit = NCrit; }
60
66 void SetFieldCalculation(const std::string &option = "coulomb",
67 const int nof_approx = 1.) {
68 m_sFieldOption = std::move(option);
69 m_iFieldApprox = std::move(nof_approx);
70 }
71
74 void SetK(float option = 0.95) { m_fStreamerK = option; }
75
77 void SetStopAtK(bool option = true) { m_bStopAtK = option; }
78
80 void SetSensor(Sensor *sensor);
81
90 void Set2dGrid(double zmin, double zmax, int zsteps, double rmax, int rsteps);
91
95
97 void AddElectron(double x, double y, double z, double t = 0, int n = 1);
98
101 void AddExtraElectron(double y, int n = 1);
102
105 void StartGridAvalanche(double dtime = -1);
106
108 long GetAvalancheSize() const { return m_nTotPosIons; }
109
112
114 [[nodiscard]] long ReachedKPercent() const {
115 if (m_bFieldK)
116 return m_lElectronsK;
117 else
118 return -1;
119 }
120
122 [[nodiscard]] const std::vector<std::pair<double, long>>
125 }
126
131 void ExportGrid(const std::string &filename);
132
133 private:
134 struct GridNode {
135 long nElectron = 0;
136 double nPosIon = 0;
137 double nNegIon = 0;
138 // holder memories for stepping in time:
140 double nPosIonHolder = 0;
141 double nNegIonHolder = 0;
142
143 double townsend = 0;
144 double attachment = 0;
146 double velocity = 0.;
148 double dSigmaL = 0;
150 double dSigmaT = 0;
151
152 double Wv = 0;
153 double Wr = 0;
155 double townsendPT = 0;
157 double attachmentPT = 0;
159 double eFieldR = 0;
161 double eFieldZ = 0;
162
163 double time = 0.;
164
165 bool anode = false;
167 int layerIndex = 0;
169 int gasGapIndex = 0;
170 bool isGasGap = true;
171 };
172
173 struct Point {
174 double x, y, z;
175 double t;
176
178 };
179
180 // Prepare grid and place stored electrons from AvalancheMicroscopic import
182
183 // Assign electron to the closest grid point
184 bool SnapTo2dGrid(double x, double y, double z, long n = 1, int gasLayer = 0);
185
186 // Prepare the mesh with the ComponentParallelPlate
188
189 // Transports the electrons/nodes a timestep
191
192 // Diffuses the electrons/nodes a timestep
193 void DiffuseTimeStep(double dx, long nElectron, double nPosIon,
194 double nNegIon, int iz, int ir, int gasGap);
195
196 // Redistributes the charges
197 void DistributeCharges(long nElectron, double nPosIon, double nNegIon, int iz,
198 int ir, double stepZ, double stepR, int gasGap);
199
200 // Calculate the field from all the contributions to the bin of interest. May
201 // need much more functionalities/tables.
202 void GetLocalField(int iz, int ir, double &eFieldZ, double &eFieldR,
203 const std::string &fieldOption, int gasGap);
204
205 // Calculate the field of charged ring in vacuum using coulomb potentials and
206 // indices
207 void GetFreeChargedRing(int iz, int ir, int fz, int fr, double &eFieldZ,
208 double &eFieldR);
209
210 // Calculate the field of charged ring in vacuum using coulomb potentials and
211 // coordinates
212 void GetFreeChargedRing(double zi, double ri, double zf, double rf,
213 double &eFieldZ, double &eFieldR);
214
215 // Get field at (zi, ri) from N charges at (zf, rf) either as a ring or a
216 // coulomb ball (rf = 0) if i and f are too close it is considered as self
217 // interaction and not included
218 bool AddFieldFromChargeAt(int iz, int ir, int fz, int fr, double N,
219 double &eFieldZ, double &eFieldR);
220
221 // Get field at (zi, ri) from N charges at (zf, rf) either as a ring or a
222 // coulomb ball (rf = 0) if i and f are too close it is considered as self
223 // interaction and not included
224 bool AddFieldFromChargeAt(int iz, int ir, double zf, double rf, double N,
225 double &eFieldZ, double &eFieldR);
226
227 // Get swarm parameters at electric field magnitude
228 void GetSwarmParameters(double MagEField, double &alpha, double &eta,
229 double &drift, double &dSigmaL, double &dSigmaT,
230 double &wv, double &wr, double &alphaPT,
231 double &etaPT, int gasGap);
232
233 // Change from 2dGrid to Global coordinates
234 void GetGlobalCoordinates(double r, double z, double phi, double &xg,
235 double &yg, double &zg, int gasGap);
236
237 // Import elliptic integral values
238 void ImportEllipticIntegralValues(const std::string &filename);
239
240 // Gets elliptic integrals via list
241 void GetEllipticIntegrals(double x, double &K, double &E);
242
243 // Get from index the gas gap number, else -1
244 int GetGasGapNumber(int layerIndex) {
245 auto it =
246 std::find(m_vIndexGasGaps.begin(), m_vIndexGasGaps.end(), layerIndex);
247 return (it != m_vIndexGasGaps.end())
248 ? std::distance(m_vIndexGasGaps.begin(), it)
249 : -1;
250 }
251
252 private:
253 std::string m_className = "AvalancheGridSpaceCharge";
254
255 bool m_bDebug = false;
256 bool m_bDiffusion = false;
257 bool m_bStick = true;
258 // boolean for AvalancheElectron
259 bool m_bDriftAvalanche = false;
260 // boolean for ImportElectronsFromAvalancheMicroscopic
261 bool m_bImportAvalanche = false;
263 long m_lNCrit = 1e8;
264 bool m_bSpaceCharge = true;
265
266 float m_fStreamerK = 0.95;
267 bool m_bStopAtK = false;
268 bool m_bFieldK = false;
270
271 bool m_bAdaptiveTime = true;
272 bool m_bImportElliptic = false;
275 bool m_bUseTOF = true;
277 bool m_bWrAvailable = true;
279 bool m_bRatesAvailable = true;
280
281 bool m_bMC = true;
282
283 int m_iFieldApprox = 1; //< order of approximation in Set(1,2,3,...)
284 double m_dMinGroups = 50; //< same values as lippmann
285
287 Sensor *m_sensor = nullptr;
288
289 std::vector<double> m_zGrid;
290 int m_zSteps = 0;
291 double m_zStepSize = 0.;
292
293 std::vector<double> m_rGrid;
294 int m_rSteps = 0.;
295 double m_rStepSize = 0.;
296
297 bool m_isgridset = false;
298 long m_nTotElectron = 0;
299 long m_nTotPosIons = 0;
300
301 double m_time = 0.;
302 double m_time0 = 0.;
303 double m_dt = 0.;
305 bool m_run = true;
306
307 std::vector<std::vector<int>>
309
310 std::vector<std::vector<GridNode>> m_grid;
312 std::vector<std::vector<Point>> m_vElectrons;
313
314 std::vector<std::pair<double, long>> m_vNElectronEvolution;
315 std::vector<long> m_vGroupSizes = {
316 1500, 800, 400, 200, 100,
317 50, 20, 10, 5, 2};
318
319 std::vector<int> m_vIndexGasGaps = {0};
320
323 std::vector<std::vector<double>> m_vCoNGasLayer{};
325 std::vector<double> m_vYPointInGasGap{};
326
328 std::vector<double> m_ezBkg = {0};
330 std::vector<int> m_vSaturatedGaps{};
331
332 std::string m_sFieldOption = "coulomb";
333 std::vector<double> m_vXElliptic;
334 std::vector<double> m_vKElliptic;
335 std::vector<double> m_vEElliptic;
336};
337
338} // namespace Garfield
339
340#endif // GARFIELD_AVALANCHEGRIDSPACECHARGE_HH
std::vector< std::vector< Point > > m_vElectrons
Electrons to transfer onto grid.
std::vector< double > m_zGrid
Grid points of z-coordinate.
void EnableStickyAnode(const bool option=true)
Enable sticky anode (default on)
bool m_isgridset
Distance between the grid points.
bool AddFieldFromChargeAt(int iz, int ir, int fz, int fr, double N, double &eFieldZ, double &eFieldR)
void Set2dGrid(double zmin, double zmax, int zsteps, double rmax, int rsteps)
void EnableMC(const bool option=true)
Enable Monte Carlo (gain/diffusion) up to a certain total number of electrons e.g.
void ExportGrid(const std::string &filename)
Export the current grid to a txt file (electron, ions numbers and field magnitude) Take care: only wh...
void AddExtraElectron(double y, int n=1)
After calling AddElectron, add more electrons on the same transversal line (y-freedom).
std::vector< std::pair< double, long > > m_vNElectronEvolution
std::vector< double > m_vYPointInGasGap
Example point (y-coord) in each gas gap.
std::vector< std::vector< int > > m_zGasGapBoundaries
[k] -> {izLeft, ..., izRight}
std::vector< double > m_ezBkg
Uniform background field in z direction, can be negative.
std::vector< int > m_vSaturatedGaps
Which gas gaps are saturated if saturation is on.
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 (rel...
std::vector< double > m_rGrid
Distance between the grid points.
bool m_bRatesAvailable
Flag if temporal rates are available to the simulation.
void EnableDiffusion(const bool option=true)
Enable diffusion (default off)
std::vector< std::vector< GridNode > > m_grid
grid with nodes on it
void GetEllipticIntegrals(double x, double &K, double &E)
bool m_bWrAvailable
Flag if bulk drift velocity is available to the simulation.
void SetStopAtK(bool option=true)
Stop the avalanche if K % field is reached.
std::vector< std::vector< double > > m_vCoNGasLayer
Coordinates of center of electron number.
void SetSensor(Sensor *sensor)
Set the sensor (and determine if it includes a parallel-plate component).
long m_nTotElectron
Total amount of electrons at time step.
void EnableSpaceChargeEffect(const bool option=true)
Enable space charge calculations (default on)
long m_nTotPosIons
total amount of charge created
std::vector< int > m_vIndexGasGaps
Which layer indices are gas layers.
void DiffuseTimeStep(double dx, long nElectron, double nPosIon, double nNegIon, int iz, int ir, int gasGap)
bool m_run
Tracking if the charges are still in the drift gap.
long ReachedKPercent() const
Returns if 100 * K % background charge has been reached.
bool SnapTo2dGrid(double x, double y, double z, long n=1, int gasLayer=0)
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 EnableAdaptiveTimeStepping(const bool option=true)
Enable adaptive time stepping (default on)
void GetFreeChargedRing(double zi, double ri, double zf, double rf, double &eFieldZ, double &eFieldR)
~AvalancheGridSpaceCharge()=default
Destructor.
void ImportEllipticIntegralValues(const std::string &filename)
void AddElectrons(AvalancheMicroscopic *avmc)
Import electron (no ions) data from AvalancheMicroscopic class to axi-symmetric grid.
void DistributeCharges(long nElectron, double nPosIon, double nNegIon, int iz, int ir, double stepZ, double stepR, int gasGap)
std::vector< long > m_vGroupSizes
same values as Lippmann & Riegler
bool AddFieldFromChargeAt(int iz, int ir, double zf, double rf, double N, double &eFieldZ, double &eFieldR)
double GetMeanDistance()
Return current mean distance of the electrons on the grid.
void SetNCrit(const long NCrit=1e8)
Set electron saturation value applied for each gas gap if space charge effect is turned off (default ...
long GetAvalancheSize() const
Returns the total positive charge in the gap's.
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 EnableTOF(const bool option=true)
Enable to use TOF swarm parameters (default on)
AvalancheGridSpaceCharge(Sensor *sensor)
Constructor.
void Reset()
Reset the charges.
void SetK(float option=0.95)
Set the streamer-inception criterion constant K in the interval (0, ∞) s.t.
bool m_bUseTOF
Flag if TOF parameters should be used, else Magboltz drift and SST spatial coefficients.
void AddElectron(double x, double y, double z, double t=0, int n=1)
Set n electrons onto the grid.
void StartGridAvalanche(double dtime=-1)
Starts the simulation with the imported electrons for a time step dt (dt = -1 until there are no elec...
void GetGlobalCoordinates(double r, double z, double phi, double &xg, double &yg, double &zg, int gasGap)
const std::vector< std::pair< double, long > > GetElectronEvolution() const
Returns the total electron number evolution.
void EnableDebugging(const bool option=true)
Enable/disable debugging (log messages) (default off)
Calculate electron drift lines and avalanches using microscopic tracking.
Component for parallel-plate geometries.
double attachmentPT
Attachment rate from TOF experiment 1/ns -> 1/cm.
int gasGapIndex
Gas gap index: -1 if not gas gap; starts with 0, 1, ...
double velocity
Magnitude of velocity of the node (not negative) cm/ns.
double eFieldR
Space-charge electric field in R direction (can be negative)
double townsendPT
Ionization rate from TOF experiment 1/ns -> 1/cm.
double nPosIon
pos ion on node (smeared values allowed)
double dSigmaT
Diffusion transverse to E (radial, phi dir is net 0).
double eFieldZ
Space-charge electric field in Z direction (can be negative)
int layerIndex
LayerIndex in ParallelPlate convention != gas gap index.
double nNegIon
neg ion on node (smeared values allowed)