Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
TrackElectron.cc
Go to the documentation of this file.
1#include <iostream>
2#include <cmath>
3
4#include "Sensor.hh"
5#include "TrackElectron.hh"
8#include "Random.hh"
9
10namespace Garfield {
11
13 : m_ready(false),
14 m_x(0.),
15 m_y(0.),
16 m_z(0.),
17 m_t(0.),
18 m_dx(0.),
19 m_dy(0),
20 m_dz(1.),
21 m_mediumName(""),
22 m_mediumDensity(0.),
23 m_mfp(0.) {
24
25 m_className = "TrackElectron";
26
27 // Setup the particle properties.
28 m_q = -1;
29 m_spin = 1;
30 m_mass = ElectronMass;
31 m_isElectron = true;
32 SetBetaGamma(3.);
33 m_particleName = "electron";
34
35}
36
37void TrackElectron::SetParticle(const std::string& particle) {
38
39 if (particle != "electron" && particle != "e" && particle != "e-") {
40 std::cerr << m_className << "::SetParticle:\n";
41 std::cerr << " Only electrons can be transported.\n";
42 }
43}
44
45bool TrackElectron::NewTrack(const double x0, const double y0, const double z0,
46 const double t0, const double dx0,
47 const double dy0, const double dz0) {
48
49 m_ready = false;
50
51 // Make sure the sensor has been set.
52 if (!m_sensor) {
53 std::cerr << m_className << "::NewTrack:\n";
54 std::cerr << " Sensor is not defined.\n";
55 return false;
56 }
57
58 // Get the medium at this location and check if it is "ionisable".
59 Medium* medium = NULL;
60 if (!m_sensor->GetMedium(x0, y0, z0, medium)) {
61 std::cerr << m_className << "::NewTrack:\n";
62 std::cerr << " No medium at initial position.\n";
63 return false;
64 }
65 if (!medium->IsIonisable()) {
66 std::cerr << m_className << "::NewTrack:\n";
67 std::cerr << " Medium at initial position is not ionisable.\n";
68 return false;
69 }
70
71 // Check if the medium is a gas.
72 if (!medium->IsGas()) {
73 std::cerr << m_className << "::NewTrack:\n";
74 std::cerr << " Medium at initial position is not a gas.\n";
75 return false;
76 }
77
78 if (!SetupGas(medium)) {
79 std::cerr << m_className << "::NewTrack:\n";
80 std::cerr << " Properties of medium " << medium->GetName()
81 << " are not available.\n";
82 return false;
83 }
84
85 if (!UpdateCrossSection()) {
86 std::cerr << m_className << "::NewTrack:\n";
87 std::cerr << " Cross-sections could not be calculated.\n";
88 return false;
89 }
90
91 m_mediumName = medium->GetName();
92
93 m_x = x0;
94 m_y = y0;
95 m_z = z0;
96 m_t = t0;
97 const double dd = sqrt(dx0 * dx0 + dy0 * dy0 + dz0 * dz0);
98 if (dd < Small) {
99 if (m_debug) {
100 std::cout << m_className << "::NewTrack:\n";
101 std::cout << " Direction vector has zero norm.\n";
102 std::cout << " Initial direction is randomized.\n";
103 }
104 const double ctheta = 1. - 2. * RndmUniform();
105 const double stheta = sqrt(1. - ctheta * ctheta);
106 const double phi = TwoPi * RndmUniform();
107 m_dx = cos(phi) * stheta;
108 m_dy = sin(phi) * stheta;
109 m_dz = ctheta;
110 } else {
111 // Normalize the direction vector.
112 m_dx = dx0 / dd;
113 m_dy = dy0 / dd;
114 m_dz = dz0 / dd;
115 }
116
117 m_ready = true;
118 return true;
119}
120
121bool TrackElectron::GetCluster(double& xcls, double& ycls, double& zcls,
122 double& tcls, int& ncls, double& edep,
123 double& extra) {
124
125 edep = extra = 0.;
126 ncls = 0;
127
128 m_electrons.clear();
129
130 if (!m_ready) {
131 std::cerr << m_className << "::GetCluster:\n"
132 << " Track not initialized. Call NewTrack first.\n";
133 return false;
134 }
135
136 // Draw a step length and propagate the electron.
137 const double d = -m_mfp * log(RndmUniformPos());
138 m_x += d * m_dx;
139 m_y += d * m_dy;
140 m_z += d * m_dz;
141 m_t += d / (sqrt(m_beta2) * SpeedOfLight);
142
143 if (!m_sensor->IsInArea(m_x, m_y, m_z)) {
144 m_ready = false;
145 return false;
146 }
147
148 Medium* medium = NULL;
149 if (!m_sensor->GetMedium(m_x, m_y, m_z, medium)) {
150 m_ready = false;
151 return false;
152 }
153
154 if (medium->GetName() != m_mediumName ||
155 medium->GetNumberDensity() != m_mediumDensity || !medium->IsIonisable()) {
156 m_ready = false;
157 return false;
158 }
159
160 xcls = m_x;
161 ycls = m_y;
162 zcls = m_z;
163 tcls = m_t;
164 const double r = RndmUniform();
165 int iComponent = 0;
166 const int nComponents = m_components.size();
167 for (int i = 0; i < nComponents; ++i) {
168 if (r <= RndmUniform()) {
169 iComponent = i;
170 break;
171 }
172 }
173
174 // Sample secondary electron energy according to
175 // Opal-Beaty-Peterson splitting function.
176 const double e0 = ElectronMass * (sqrt(1. / (1. - m_beta2)) - 1.);
177 double esec = m_components[iComponent].wSplit *
178 tan(RndmUniform() * atan((e0 - m_components[iComponent].ethr) /
179 (2. * m_components[iComponent].wSplit)));
180 esec = m_components[iComponent].wSplit *
181 pow(esec / m_components[iComponent].wSplit, 0.9524);
182 m_electrons.resize(1);
183 m_electrons[0].energy = esec;
184 m_electrons[0].x = xcls;
185 m_electrons[0].y = ycls;
186 m_electrons[0].z = zcls;
187
188 ncls = 1;
189 edep = esec;
190
191 return true;
192}
193
195
196 if (!m_ready) {
197 std::cerr << m_className << "::GetClusterDensity:\n";
198 std::cerr << " Track has not been initialized.\n";
199 return 0.;
200 }
201
202 if (m_mfp <= 0.) {
203 std::cerr << m_className << "::GetClusterDensity:\n";
204 std::cerr << " Mean free path is not available.\n";
205 return 0.;
206 }
207
208 return 1. / m_mfp;
209}
210
212
213 if (!m_ready) {
214 std::cerr << m_className << "::GetStoppingPower:\n";
215 std::cerr << " Track has not been initialised.\n";
216 return 0.;
217 }
218
219 const double prefactor = 4 * Pi * pow(HbarC / ElectronMass, 2);
220 const double lnBg2 = log(m_beta2 / (1. - m_beta2));
221
222 double dedx = 0.;
223 // Primary energy
224 const double e0 = ElectronMass * (sqrt(1. / (1. - m_beta2)) - 1.);
225 const int nComponents = m_components.size();
226 for (int i = nComponents; i--;) {
227 // Calculate the mean number of clusters per cm.
228 const double cmean =
229 m_mediumDensity * m_components[i].fraction * (prefactor / m_beta2) *
230 (m_components[i].m2Ion * (lnBg2 - m_beta2) + m_components[i].cIon);
231 const double ew = (e0 - m_components[i].ethr) / (2 * m_components[i].wSplit);
232 // Calculate the mean secondary electron energy.
233 const double emean =
234 (m_components[i].wSplit / (2 * atan(ew))) * log(1. + ew * ew);
235 dedx += cmean * emean;
236 }
237
238 return dedx;
239}
240
241bool TrackElectron::SetupGas(Medium* gas) {
242
243 m_components.clear();
244
245 if (!gas) {
246 std::cerr << m_className << "::SetupGas:\n";
247 std::cerr << " Medium is not defined.\n";
248 return false;
249 }
250
251 m_mediumDensity = gas->GetNumberDensity();
252 const int nComponents = gas->GetNumberOfComponents();
253 if (nComponents <= 0) {
254 std::cerr << m_className << "::SetupGas:\n";
255 std::cerr << " Medium composition is not defined.\n";
256 return false;
257 }
258 m_components.resize(nComponents);
259
260 // Density correction parameters from
261 // R. M. Sternheimer, M. J. Berger, S. M. Seltzer,
262 // Atomic Data and Nuclear Data Tables 30 (1984), 261-271
263 bool ok = true;
264 for (int i = nComponents; i--;) {
265 std::string gasname = "";
266 double frac = 0.;
267 gas->GetComponent(i, gasname, frac);
268 m_components[i].fraction = frac;
269 m_components[i].p = 0.;
270 if (gasname == "CF4") {
271 m_components[i].m2Ion = 7.2;
272 m_components[i].cIon = 93.;
273 m_components[i].x0Dens = 1.;
274 m_components[i].x1Dens = 0.;
275 m_components[i].cDens = 0.;
276 m_components[i].aDens = 0.;
277 m_components[i].mDens = 0.;
278 m_components[i].ethr = 15.9;
279 m_components[i].wSplit = 19.5;
280 } else if (gasname == "Ar") {
281 m_components[i].m2Ion = 3.593;
282 m_components[i].cIon = 39.7;
283 m_components[i].x0Dens = 1.7635;
284 m_components[i].x1Dens = 4.4855;
285 m_components[i].cDens = 11.9480;
286 m_components[i].aDens = 0.19714;
287 m_components[i].mDens = 2.9618;
288 m_components[i].ethr = 15.75961;
289 m_components[i].wSplit = 15.;
290 } else if (gasname == "He") {
291 m_components[i].m2Ion = 0.489;
292 m_components[i].cIon = 5.5;
293 m_components[i].x0Dens = 2.2017;
294 m_components[i].x1Dens = 3.6122;
295 m_components[i].cDens = 11.1393;
296 m_components[i].aDens = 0.13443;
297 m_components[i].mDens = 5.8347;
298 m_components[i].ethr = 24.58739;
299 m_components[i].wSplit = 10.5;
300 } else if (gasname == "He-3") {
301 m_components[i].m2Ion = 0.489;
302 m_components[i].cIon = 5.5;
303 m_components[i].x0Dens = 2.2017;
304 m_components[i].x1Dens = 3.6122;
305 m_components[i].cDens = 11.1393;
306 m_components[i].aDens = 0.13443;
307 m_components[i].mDens = 5.8347;
308 m_components[i].ethr = 24.58739;
309 m_components[i].wSplit = 10.5;
310 } else if (gasname == "Ne") {
311 m_components[i].m2Ion = 1.69;
312 m_components[i].cIon = 17.8;
313 m_components[i].x0Dens = 2.0735;
314 m_components[i].x1Dens = 4.6421;
315 m_components[i].cDens = 11.9041;
316 m_components[i].aDens = 0.08064;
317 m_components[i].mDens = 3.5771;
318 m_components[i].ethr = 21.56454;
319 m_components[i].wSplit = 19.5;
320 } else if (gasname == "Kr") {
321 m_components[i].m2Ion = 5.5;
322 m_components[i].cIon = 56.9;
323 m_components[i].x0Dens = 1.7158;
324 m_components[i].x1Dens = 5.0748;
325 m_components[i].cDens = 12.5115;
326 m_components[i].aDens = 0.07446;
327 m_components[i].mDens = 3.4051;
328 m_components[i].ethr = 13.996;
329 m_components[i].wSplit = 21.;
330 } else if (gasname == "Xe") {
331 m_components[i].m2Ion = 8.04;
332 m_components[i].cIon = 75.25;
333 m_components[i].x0Dens = 1.5630;
334 m_components[i].x1Dens = 4.7371;
335 m_components[i].cDens = 12.7281;
336 m_components[i].aDens = 0.23314;
337 m_components[i].mDens = 2.7414;
338 m_components[i].ethr = 12.129843;
339 m_components[i].wSplit = 23.7;
340 } else if (gasname == "CH4") {
341 m_components[i].m2Ion = 3.75;
342 m_components[i].cIon = 42.5;
343 m_components[i].x0Dens = 1.6263;
344 m_components[i].x1Dens = 3.9716;
345 m_components[i].cDens = 9.5243;
346 m_components[i].aDens = 0.09253;
347 m_components[i].mDens = 3.6257;
348 m_components[i].ethr = 12.65;
349 m_components[i].wSplit = 8.;
350 } else if (gasname == "iC4H10") {
351 m_components[i].m2Ion = 15.5;
352 m_components[i].cIon = 160.;
353 m_components[i].x0Dens = 1.3788;
354 m_components[i].x1Dens = 3.7524;
355 m_components[i].cDens = 8.5633;
356 m_components[i].aDens = 0.10852;
357 m_components[i].mDens = 3.4884;
358 m_components[i].ethr = 10.67;
359 m_components[i].wSplit = 7.;
360 } else if (gasname == "CO2") {
361 m_components[i].m2Ion = 5.6;
362 m_components[i].cIon = 57.91;
363 m_components[i].x0Dens = 1.6294;
364 m_components[i].x1Dens = 4.1825;
365 m_components[i].aDens = 0.11768;
366 m_components[i].mDens = 3.3227;
367 m_components[i].ethr = 13.777;
368 m_components[i].wSplit = 13.;
369 } else if (gasname == "N2") {
370 m_components[i].m2Ion = 3.35;
371 m_components[i].cIon = 38.1;
372 m_components[i].x0Dens = 1.7378;
373 m_components[i].x1Dens = 4.1323;
374 m_components[i].cDens = 10.5400;
375 m_components[i].aDens = 0.15349;
376 m_components[i].mDens = 3.2125;
377 m_components[i].ethr = 15.581;
378 m_components[i].wSplit = 13.8;
379 } else {
380 std::cerr << m_className << "::SetupGas:\n";
381 std::cerr << " Cross-section for " << gasname
382 << " is not available.\n";
383 ok = false;
384 }
385 }
386
387 if (!ok) {
388 m_components.clear();
389 }
390
391 return true;
392}
393
394bool TrackElectron::UpdateCrossSection() {
395
396 const double prefactor = 4 * Pi * pow(HbarC / ElectronMass, 2);
397 const double lnBg2 = log(m_beta2 / (1. - m_beta2));
398 // Parameter X in the Sternheimer fit formula
399 const double eta = m_mediumDensity / LoschmidtNumber;
400 const double x = 0.5 * (lnBg2 + log(eta)) / log(10.);
401 double csSum = 0.;
402 const int nComponents = m_components.size();
403 for (int i = nComponents; i--;) {
404 double delta = 0.;
405 if (m_components[i].x0Dens < m_components[i].x1Dens &&
406 x >= m_components[i].x0Dens) {
407 delta = 2 * log(10.) * x - m_components[i].cDens;
408 if (x < m_components[i].x1Dens) {
409 delta += m_components[i].aDens *
410 pow(m_components[i].x1Dens - x, m_components[i].mDens);
411 }
412 }
413 const double cs = (m_components[i].fraction * prefactor / m_beta2) *
414 (m_components[i].m2Ion * (lnBg2 - m_beta2 - delta) +
415 m_components[i].cIon);
416 m_components[i].p = cs;
417 csSum += cs;
418 }
419
420 if (csSum <= 0.) {
421 std::cerr << m_className << "::UpdateCrossSection:\n";
422 std::cerr << " Total cross-section <= 0.\n";
423 return false;
424 }
425
426 m_mfp = 1. / (csSum * m_mediumDensity);
427
428 for (int i = 0; i < nComponents; ++i) {
429 m_components[i].p /= csSum;
430 if (i > 0) m_components[i].p += m_components[i - 1].p;
431 }
432
433 return true;
434}
435}
Abstract base class for media.
Definition: Medium.hh:11
virtual void GetComponent(const unsigned int i, std::string &label, double &f)
Definition: Medium.cc:142
virtual double GetNumberDensity() const
Definition: Medium.hh:47
unsigned int GetNumberOfComponents() const
Definition: Medium.hh:37
const std::string & GetName() const
Definition: Medium.hh:22
virtual bool IsGas() const
Definition: Medium.hh:23
bool IsIonisable() const
Definition: Medium.hh:59
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
Definition: Sensor.cc:264
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:150
virtual void SetParticle(const std::string &particle)
Set the type of particle.
virtual double GetClusterDensity()
virtual bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0)
virtual bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &ncls, double &ecls, double &extra)
virtual double GetStoppingPower()
Get the stopping power (mean energy loss [eV] per cm).
Sensor * m_sensor
Definition: Track.hh:90
void SetBetaGamma(const double bg)
Set the relative momentum of the particle.
Definition: Track.cc:116
bool m_debug
Definition: Track.hh:97
bool m_isElectron
Definition: Track.hh:87
double m_q
Definition: Track.hh:82
std::string m_particleName
Definition: Track.hh:88
double m_beta2
Definition: Track.hh:86
std::string m_className
Definition: Track.hh:80
double m_mass
Definition: Track.hh:84
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition: Random.hh:17
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337