Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
HeedDeltaElectron.cpp
Go to the documentation of this file.
7
8// 2003, I. Smirnov
9
10#define USE_ADJUSTED_W
11#define RANDOM_POIS
12
13namespace {
14
15long findInterval(Heed::EnergyMesh* emesh, const double energy) {
16
17 const long n = emesh->get_interval_number_between_centers(energy);
18 return std::min(std::max(n, 0L), emesh->get_q() - 2);
19}
20
21double interpolate(Heed::EnergyMesh* emesh, const double x,
22 const std::vector<double>& y) {
23
24 const long n = findInterval(emesh, x);
25 const double x1 = emesh->get_ec(n);
26 const double x2 = emesh->get_ec(n + 1);
27 const double y1 = y[n];
28 const double y2 = y[n + 1];
29 return y1 + (x - x1) * (y2 - y1) / (x2 - x1);
30}
31
32double sample_ctheta(const double sigma) {
33 double ctheta = 0.;
34 do {
35 double y = Heed::SRANLUX();
36 while (y == 1.) {
37 y = Heed::SRANLUX();
38 }
39 const double x = sigma * (-log(1.0 - y));
40 ctheta = 1. - x;
41 } while (ctheta <= -1.0);
42 return ctheta;
43}
44
45}
46
47namespace Heed {
48
49using CLHEP::degree;
50using CLHEP::cm;
51using CLHEP::eV;
52using CLHEP::keV;
53using CLHEP::MeV;
54using CLHEP::c_light;
55using CLHEP::c_squared;
56
57bool HeedDeltaElectron::s_low_mult_scattering = true;
58bool HeedDeltaElectron::s_high_mult_scattering = true;
59bool HeedDeltaElectron::s_direct_low_if_little = true;
60
62 const vec& vel, vfloat ftime,
63 long fparent_particle_number,
64 HeedFieldMap* fieldmap,
65 bool fprint_listing)
66 : eparticle(primvol, pt, vel, ftime, &electron_def, fieldmap),
67 parent_particle_number(fparent_particle_number),
68 m_particle_number(s_counter++),
69 m_print_listing(fprint_listing) {
70 mfunname("HeedDeltaElectron::HeedDeltaElectron(...)");
71}
72
73void HeedDeltaElectron::physics_mrange(double& fmrange) {
74 mfunname("void HeedDeltaElectron::physics_mrange(double& fmrange)");
75 if (m_print_listing) mcout << "HeedDeltaElectron::physics_mrange\n";
76
77 m_mult_low_path_length = false;
78 m_q_low_path_length = 0.0;
79 m_path_length = false;
80 if (fmrange <= 0.0) return;
81 if (m_curr_ekin <= 0.0) {
82 fmrange = 0.0;
83 return;
84 }
85 // Get local volume and convert it to a cross-section object.
86 const absvol* av = m_currpos.volume();
87 auto hdecs = dynamic_cast<const HeedDeltaElectronCS*>(av);
88 if (!hdecs) return;
89 if (m_print_listing) Iprintnf(mcout, fmrange);
90 const double ek = m_curr_ekin / MeV;
91 // Get the dE/dx at this kinetic energy.
92 EnergyMesh* emesh = hdecs->hmd->energy_mesh;
93 const double dedx = interpolate(emesh, ek, hdecs->eLoss);
94 // Min. loss 50 eV.
95 double eloss = std::max(0.1 * ek, 0.00005);
96 if (eloss > ek) {
97 eloss = ek;
98 m_stop_eloss = true;
99 } else {
100 m_stop_eloss = false;
101 }
102 fmrange = std::min(fmrange, (eloss / dedx) * cm);
103 if (m_print_listing) Iprint2nf(mcout, fmrange, ek);
104 const double ek_restr = std::max(ek, 0.0005);
105 if (m_print_listing) Iprintnf(mcout, ek_restr / keV);
106
107 double low_path_length = 0.; // in internal units
108 if (s_low_mult_scattering) {
109 low_path_length = interpolate(emesh, ek_restr, hdecs->low_lambda) * cm;
110 if (m_print_listing) Iprintnf(mcout, low_path_length / cm);
111 long qscat = hdecs->eesls->get_qscat();
112 const double sigma_ctheta = hdecs->get_sigma(ek_restr, qscat);
113 // Reduce the number of scatterings, if the angle is too large.
114 if (sigma_ctheta > 0.3) qscat = long(qscat * 0.3 / sigma_ctheta);
115 const double mult_low_path_length = qscat * low_path_length;
116 if (m_print_listing) Iprintnf(mcout, mult_low_path_length);
117 if (fmrange > mult_low_path_length) {
118 fmrange = mult_low_path_length;
119 m_mult_low_path_length = true;
120 m_q_low_path_length = hdecs->eesls->get_qscat();
121 m_stop_eloss = false;
122 } else {
123 m_mult_low_path_length = false;
124 m_q_low_path_length = fmrange / low_path_length;
125 }
126 if (m_print_listing) Iprint2nf(mcout, fmrange, m_q_low_path_length);
127 }
128
129 if (s_high_mult_scattering) {
130 const double mean_path = interpolate(emesh, ek_restr, hdecs->lambda);
131 if (m_print_listing) Iprintnf(mcout, mean_path);
132 const double path_length = -mean_path * cm * log(1.0 - SRANLUX());
133 if (m_print_listing) Iprintnf(mcout, path_length);
134 if (fmrange > path_length) {
135 fmrange = path_length;
136 m_path_length = true;
137 m_mult_low_path_length = true;
138 if (s_low_mult_scattering) {
139 m_q_low_path_length = fmrange / low_path_length;
140 if (m_print_listing) Iprintnf(mcout, m_q_low_path_length);
141 }
142 m_stop_eloss = false;
143 } else {
144 m_path_length = false;
145 }
146 if (m_print_listing) Iprintnf(mcout, fmrange);
147 }
148 m_phys_mrange = fmrange;
149}
150
152 std::vector<gparticle*>& /*secondaries*/) {
153 mfunname("void HeedDeltaElectron::physics_after_new_speed()");
154 if (m_print_listing) {
155 mcout << "HeedDeltaElectron::physics_after_new_speed\n";
157 }
159 if (m_currpos.prange <= 0.0) {
160 if (m_curr_ekin <= 0.0) {
161 // Get local volume.
162 absvol* av = m_currpos.volume();
163 if (av && av->s_sensitive && m_fieldMap->inside(m_currpos.ptloc)) {
164 if (m_print_listing) mcout << "Convert to conduction electron.\n";
165 conduction_electrons.emplace_back(
167 }
168 m_alive = false;
169 }
170 if (m_print_listing) mcout << "exit due to currpos.prange <= 0.0\n";
171 return;
172 }
173 // Get local volume and convert it to a cross-section object.
174 const absvol* av = m_currpos.volume();
175 auto hdecs = dynamic_cast<const HeedDeltaElectronCS*>(av);
176 if (!hdecs) return;
177 double ek = m_curr_ekin / MeV;
178 if (m_print_listing) {
179 Iprintnf(mcout, ek);
180 Iprint3n(mcout, m_stop_eloss, m_phys_mrange, m_currpos.prange);
181 }
182 // Calculate dE/dx and energy loss. Update the kinetic energy.
183 double dedx;
184 double Eloss = 0.;
185 if (m_stop_eloss && m_phys_mrange == m_currpos.prange) {
186 Eloss = m_curr_ekin;
187 m_curr_ekin = 0.0;
188 dedx = Eloss / m_currpos.prange / (MeV / cm);
189 } else {
190 EnergyMesh* emesh = hdecs->hmd->energy_mesh;
191 dedx = interpolate(emesh, ek, hdecs->eLoss);
192 Eloss = std::min(m_currpos.prange * dedx * MeV / cm, m_curr_ekin);
193 m_total_eloss += Eloss;
194 m_curr_ekin -= Eloss;
195 }
196 if (m_print_listing)
197 Iprint3nf(mcout, m_prev_ekin / eV, m_curr_ekin / eV, Eloss / eV);
198 if (m_curr_ekin <= 0.0) {
199 if (m_print_listing) mcout << "m_curr_ekin <= 0.0\n";
200 m_curr_ekin = 0.0;
201 m_curr_gamma_1 = 0.0;
202 m_currpos.speed = 0.0;
203 m_alive = false;
204 } else {
205 const double resten = m_mass * c_squared;
206 m_curr_gamma_1 = m_curr_ekin / resten;
207 m_currpos.speed = c_light * lorbeta(m_curr_gamma_1);
208 }
209 absvol* vav = m_currpos.volume();
210 if (vav && vav->s_sensitive) {
211 if (m_print_listing) {
212 mcout << "volume is sensitive\n";
213 Iprint2nf(mcout, Eloss / eV, m_necessary_energy / eV);
214 }
215 if (Eloss > 0.0) ionisation(Eloss, dedx, hdecs->pairprod);
216 }
217 if (m_print_listing) {
218 mcout << '\n';
220 }
221 if (!m_alive) {
222 // Done tracing the delta electron. Create the last conduction electron.
223 vav = m_currpos.volume();
224 if (vav && vav->s_sensitive && m_fieldMap->inside(m_currpos.ptloc)) {
225 if (m_print_listing) mcout << "Last conduction electron\n";
226 conduction_electrons.emplace_back(
228 }
229 return;
230 }
231
232 if (m_print_listing) mcout << "\nstart to rotate by low angle\n";
233 double ek_restr = std::max(ek, 0.0005);
234 if (m_print_listing) Iprint2nf(mcout, m_currpos.prange, m_phys_mrange);
235 if (m_currpos.prange < m_phys_mrange) {
236 // recalculate scatterings
237 m_path_length = false;
238 if (s_low_mult_scattering) {
239 EnergyMesh* emesh = hdecs->hmd->energy_mesh;
240 const double low_path_length = interpolate(emesh, ek_restr, hdecs->low_lambda) * cm;
241 if (m_print_listing) Iprintnf(mcout, low_path_length / cm);
242 m_mult_low_path_length = false;
243 m_q_low_path_length = m_currpos.prange / low_path_length;
244 if (m_print_listing) Iprintnf(mcout, m_q_low_path_length);
245 }
246 }
247 if (m_print_listing) Iprintnf(mcout, m_q_low_path_length);
248#ifdef RANDOM_POIS
249 if (m_q_low_path_length > 0.0) {
250 m_q_low_path_length = pois(m_q_low_path_length);
251 if (m_print_listing) {
252 mcout << "After pois:\n";
253 Iprintnf(mcout, m_q_low_path_length);
254 }
255 }
256#endif
257 if (m_q_low_path_length > 0) {
258 if (s_direct_low_if_little && m_q_low_path_length < 5) {
259 // direct modeling
260 if (m_print_listing) {
261 mcout << "direct modeling of low scatterings\n";
262 Iprint(mcout, m_currpos.dir);
263 }
264 EnergyMesh* emesh = hdecs->hmd->energy_mesh;
265 const long n1r = findInterval(emesh, ek_restr);
266 for (long nscat = 0; nscat < m_q_low_path_length; ++nscat) {
267 if (m_print_listing) Iprintn(mcout, nscat);
268 const double theta =
269 hdecs->low_angular_points_ran[n1r].ran(SRANLUX()) * degree;
270 if (m_print_listing) Iprintnf(mcout, theta);
271 turn(cos(theta), sin(theta));
272 }
273 } else {
274 const double sigma = hdecs->get_sigma(ek_restr, m_q_low_path_length);
275 // actually it is mean(1-cos(theta)) or
276 // sqrt(mean(square(1-cos(theta)))) depending on USE_MEAN_COEF
277 if (m_print_listing) Iprintnf(mcout, sigma);
278 // Gauss:
279 // double ctheta = 1.0 - fabs(rnorm_improved() * sigma);
280 // Exponential distribution fits better:
281#ifdef USE_MEAN_COEF
282 const double ctheta = sample_ctheta(sigma);
283#else
284 const double ctheta = sample_ctheta(sigma / sqrt(2.));
285#endif
286 if (m_print_listing) Iprintnf(mcout, ctheta);
287 const double theta = acos(ctheta);
288 if (m_print_listing) Iprint2nf(mcout, theta, theta / degree);
289 turn(ctheta, sin(theta));
290 }
291 }
292 if (m_path_length) {
293 if (m_print_listing) {
294 mcout << "\nstarting to rotate by large angle" << std::endl;
295 Iprintnf(mcout, m_path_length);
296 }
297 EnergyMesh* emesh = hdecs->hmd->energy_mesh;
298 const long n1r = findInterval(emesh, ek_restr);
299 const double theta = hdecs->angular_points_ran[n1r].ran(SRANLUX()) * degree;
300 if (m_print_listing) Iprintnf(mcout, theta);
301 turn(cos(theta), sin(theta));
302 }
303 if (m_print_listing) Iprint2nf(mcout, m_currpos.dir, m_currpos.dirloc);
304}
305
306void HeedDeltaElectron::ionisation(const double eloss, const double dedx,
307 PairProd* pairprod) {
308
309 if (eloss < m_necessary_energy) {
310 m_necessary_energy -= eloss;
311 return;
312 }
313
314 if (m_print_listing) mcout << "\nstart to leave conduction electrons\n";
315 if (m_necessary_energy <= 0.0) {
316#ifdef USE_ADJUSTED_W
317 m_necessary_energy = pairprod->get_eloss(m_prev_ekin / eV) * eV;
318#else
319 m_necessary_energy = pairprod->get_eloss() * eV;
320#endif
321 }
322 if (m_print_listing) Iprintnf(mcout, m_necessary_energy / eV);
323 double eloss_left = eloss;
324 point curpt = m_prevpos.pt;
325 vec dir = m_prevpos.dir; // this approximation ignores curvature
326 double ekin = m_prev_ekin;
327 if (m_print_listing) Iprintnf(mcout, curpt);
328 while (eloss_left >= m_necessary_energy) {
329 const double step_length = m_necessary_energy / (dedx * MeV / cm);
330 if (m_print_listing) Iprintnf(mcout, step_length);
331 curpt = curpt + dir * step_length;
332 if (m_print_listing) Iprintf(mcout, curpt);
333 point ptloc = curpt;
334 m_prevpos.tid.up_absref(&ptloc);
335 if (m_print_listing) mcout << "New conduction electron\n";
336 if (m_fieldMap->inside(ptloc)) {
337 conduction_electrons.emplace_back(HeedCondElectron(ptloc, m_currpos.time));
338 conduction_ions.emplace_back(HeedCondElectron(ptloc, m_currpos.time));
339 }
340 eloss_left -= m_necessary_energy;
341 ekin -= m_necessary_energy;
342 if (ekin < 0.) break;
343 // Generate next random energy
344#ifdef USE_ADJUSTED_W
345 m_necessary_energy = eV * pairprod->get_eloss(ekin / eV);
346#else
347 m_necessary_energy = pairprod->get_eloss() * eV;
348#endif
349 if (m_print_listing) {
350 Iprintnf(mcout, eloss_left / eV);
351 Iprint2nf(mcout, ekin / eV, m_necessary_energy / eV);
352 }
353 }
354 m_necessary_energy -= eloss_left;
355 if (m_print_listing) Iprintnf(mcout, m_necessary_energy / eV);
356}
357
358void HeedDeltaElectron::print(std::ostream& file, int l) const {
359 if (l < 0) return;
360 file << "HeedDeltaElectron: particle_number=" << m_particle_number << "\n";
361 if (l <= 1) return;
362 file << " s_low_mult_scattering=" << s_low_mult_scattering
363 << " s_high_mult_scattering=" << s_high_mult_scattering << '\n'
364 << " phys_mrange=" << m_phys_mrange
365 << " stop_eloss=" << m_stop_eloss
366 << " mult_low_path_length=" << m_mult_low_path_length << '\n'
367 << " q_low_path_length=" << m_q_low_path_length
368 << " path_length=" << m_path_length
369 << " necessary_energy/eV=" << m_necessary_energy / eV << '\n'
370 << " parent_particle_number=" << parent_particle_number << '\n';
371 mparticle::print(file, l - 1);
372}
373}
#define check_econd11(a, signb, stream)
#define mfunname(string)
long get_q() const
Return number of bins.
Definition EnergyMesh.h:42
double get_ec(long n) const
Return center of a given bin.
Definition EnergyMesh.h:50
long get_interval_number_between_centers(const double ener) const
std::vector< HeedCondElectron > conduction_electrons
std::vector< HeedCondElectron > conduction_ions
void physics_mrange(double &fmrange) override
void print(std::ostream &file, int l) const override
Print-out.
HeedDeltaElectron()=default
Default constructor.
void physics_after_new_speed(std::vector< gparticle * > &secondaries) override
Apply any other processes (turn the trajectory, kill the particle, ...).
Retrieve electric and magnetic field from Sensor.
double get_eloss() const
Calculate energy loss (in eV).
Definition PairProd.cpp:43
bool s_sensitive
Definition volume.h:73
eparticle()=default
Default constructor.
HeedFieldMap * m_fieldMap
Pointer to field map.
Definition eparticle.h:34
stvpoint m_prevpos
Previous point.
Definition gparticle.h:274
void turn(const double ctheta, const double stheta)
Change the direction of the particle.
stvpoint m_currpos
Current point.
Definition gparticle.h:276
bool m_alive
Status flag whether the particle is active.
Definition gparticle.h:259
static std::atomic< long > s_counter
Instance counter.
Definition gparticle.h:256
Abstract base classs for volume "manipulators".
Definition volume.h:128
double m_prev_ekin
Previous kinetic energy.
Definition mparticle.h:77
double m_curr_ekin
Current kinetic energy.
Definition mparticle.h:73
void print(std::ostream &file, int l) const override
Print-out.
double m_curr_gamma_1
Current .
Definition mparticle.h:80
double m_mass
Mass (not mass * speed_of_light^2)
Definition mparticle.h:70
Point.
Definition vec.h:368
Definition BGMesh.cpp:6
DoubleAc cos(const DoubleAc &f)
Definition DoubleAc.cpp:432
double lorbeta(const double gamma_1)
as function of .
Definition lorgamma.cpp:23
particle_def electron_def("electron", "e-", electron_mass_c2/c_squared, electron_charge, 0.5)
int vecerror
Definition vec.cpp:29
DoubleAc acos(const DoubleAc &f)
Definition DoubleAc.cpp:490
long pois(const double amu)
Definition pois.cpp:9
DoubleAc sin(const DoubleAc &f)
Definition DoubleAc.cpp:384
DoubleAc sqrt(const DoubleAc &f)
Definition DoubleAc.cpp:314
double vfloat
Definition vfloat.h:16
#define Iprint3nf(file, name1, name2, name3)
Definition prstream.h:236
#define mcout
Definition prstream.h:126
#define Iprint3n(file, name1, name2, name3)
Definition prstream.h:232
#define Iprint(file, name)
Definition prstream.h:197
#define mcerr
Definition prstream.h:128
#define Iprintn(file, name)
Definition prstream.h:204
#define Iprintf(file, name)
Definition prstream.h:199
#define Iprint2nf(file, name1, name2)
Definition prstream.h:222
#define Iprint2n(file, name1, name2)
Definition prstream.h:219
#define Iprintnf(file, name)
Definition prstream.h:206