Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
MediumSilicon.cc
Go to the documentation of this file.
1#include <algorithm>
2#include <cmath>
3#include <fstream>
4#include <iostream>
5#include <sstream>
6
10#include "Garfield/Random.hh"
11#include "Garfield/Utilities.hh"
12
13namespace Garfield {
14
16 : Medium(),
17 m_eStepXL(m_eFinalXL / nEnergyStepsXL),
18 m_eStepG(m_eFinalG / nEnergyStepsG),
19 m_eStepV(m_eFinalV / nEnergyStepsV) {
20 m_className = "MediumSilicon";
21 m_name = "Si";
22
23 SetTemperature(300.);
25 SetAtomicNumber(14.);
26 SetAtomicWeight(28.0855);
27 SetMassDensity(2.329);
28
31 m_microscopic = true;
32
33 m_w = 3.6;
34 m_fano = 0.11;
35
36 m_ieMinL = int(m_eMinL / m_eStepXL) + 1;
37 m_ieMinG = int(m_eMinG / m_eStepG) + 1;
38
39 // Load the density of states table.
40 InitialiseDensityOfStates();
41}
42
43void MediumSilicon::SetDoping(const char type, const double c) {
44 if (toupper(type) == 'N') {
45 m_dopingType = 'n';
46 if (c > Small) {
47 m_dopingConcentration = c;
48 } else {
49 std::cerr << m_className << "::SetDoping:\n"
50 << " Doping concentration must be greater than zero.\n"
51 << " Using default value for n-type silicon "
52 << "(10^12 cm-3) instead.\n";
53 m_dopingConcentration = 1.e12;
54 }
55 } else if (toupper(type) == 'P') {
56 m_dopingType = 'p';
57 if (c > Small) {
58 m_dopingConcentration = c;
59 } else {
60 std::cerr << m_className << "::SetDoping:\n"
61 << " Doping concentration must be greater than zero.\n"
62 << " Using default value for p-type silicon "
63 << "(10^18 cm-3) instead.\n";
64 m_dopingConcentration = 1.e18;
65 }
66 } else if (toupper(type) == 'I') {
67 m_dopingType = 'i';
68 m_dopingConcentration = 0.;
69 } else {
70 std::cerr << m_className << "::SetDoping:\n"
71 << " Unknown dopant type (" << type << ").\n"
72 << " Available types are n, p and i (intrinsic).\n";
73 return;
74 }
75
76 m_isChanged = true;
77}
78
79void MediumSilicon::GetDoping(char& type, double& c) const {
80 type = m_dopingType;
81 c = m_dopingConcentration;
82}
83
84void MediumSilicon::SetTrapCrossSection(const double ecs, const double hcs) {
85 if (ecs < 0.) {
86 std::cerr << m_className << "::SetTrapCrossSection:\n"
87 << " Capture cross-section [cm2] must non-negative.\n";
88 } else {
89 m_eTrapCs = ecs;
90 }
91
92 if (hcs < 0.) {
93 std::cerr << m_className << "::SetTrapCrossSection:\n"
94 << " Capture cross-section [cm2] must be non-negative.n";
95 } else {
96 m_hTrapCs = hcs;
97 }
98
99 m_trappingModel = 0;
100 m_isChanged = true;
101}
102
103void MediumSilicon::SetTrapDensity(const double n) {
104 if (n < 0.) {
105 std::cerr << m_className << "::SetTrapDensity:\n"
106 << " Trap density [cm-3] must be non-negative.\n";
107 } else {
108 m_eTrapDensity = n;
109 m_hTrapDensity = n;
110 }
111
112 m_trappingModel = 0;
113 m_isChanged = true;
114}
115
116void MediumSilicon::SetTrappingTime(const double etau, const double htau) {
117 if (etau <= 0.) {
118 std::cerr << m_className << "::SetTrappingTime:\n"
119 << " Trapping time [ns-1] must be positive.\n";
120 } else {
121 m_eTrapTime = etau;
122 }
123
124 if (htau <= 0.) {
125 std::cerr << m_className << "::SetTrappingTime:\n"
126 << " Trapping time [ns-1] must be positive.\n";
127 } else {
128 m_hTrapTime = htau;
129 }
130
131 m_trappingModel = 1;
132 m_isChanged = true;
133}
134
135bool MediumSilicon::ElectronVelocity(const double ex, const double ey,
136 const double ez, const double bx,
137 const double by, const double bz,
138 double& vx, double& vy, double& vz) {
139 vx = vy = vz = 0.;
140 if (m_isChanged) {
141 if (!UpdateTransportParameters()) {
142 std::cerr << m_className << "::ElectronVelocity:\n"
143 << " Error calculating the transport parameters.\n";
144 return false;
145 }
146 m_isChanged = false;
147 }
148
149 if (!m_eVelE.empty()) {
150 // Interpolation in user table.
151 return Medium::ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
152 }
153
154 const double e = sqrt(ex * ex + ey * ey + ez * ez);
155
156 // Calculate the mobility
157 double mu;
158 switch (m_highFieldMobilityModel) {
159 case HighFieldMobility::Minimos:
160 ElectronMobilityMinimos(e, mu);
161 break;
162 case HighFieldMobility::Canali:
163 ElectronMobilityCanali(e, mu);
164 break;
165 case HighFieldMobility::Reggiani:
166 ElectronMobilityReggiani(e, mu);
167 break;
168 default:
169 mu = m_eMobility;
170 break;
171 }
172 mu = -mu;
173
174 const double b = sqrt(bx * bx + by * by + bz * bz);
175 if (b < Small) {
176 vx = mu * ex;
177 vy = mu * ey;
178 vz = mu * ez;
179 } else {
180 // Hall mobility
181 const double muH = m_eHallFactor * mu;
182 const double eb = bx * ex + by * ey + bz * ez;
183 const double mub = muH * b;
184 const double f = mu / (1. + mub * mub);
185 const double muH2 = muH * muH;
186 // Compute the drift velocity using the Langevin equation.
187 vx = f * (ex + muH * (ey * bz - ez * by) + muH2 * bx * eb);
188 vy = f * (ey + muH * (ez * bx - ex * bz) + muH2 * by * eb);
189 vz = f * (ez + muH * (ex * by - ey * bx) + muH2 * bz * eb);
190 }
191 return true;
192}
193
194bool MediumSilicon::ElectronTownsend(const double ex, const double ey,
195 const double ez, const double bx,
196 const double by, const double bz,
197 double& alpha) {
198 alpha = 0.;
199 if (m_isChanged) {
200 if (!UpdateTransportParameters()) {
201 std::cerr << m_className << "::ElectronTownsend:\n"
202 << " Error calculating the transport parameters.\n";
203 return false;
204 }
205 m_isChanged = false;
206 }
207
208 if (!m_eAlp.empty()) {
209 // Interpolation in user table.
210 return Medium::ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
211 }
212
213 const double e = sqrt(ex * ex + ey * ey + ez * ez);
214
215 switch (m_impactIonisationModel) {
216 case ImpactIonisation::VanOverstraeten:
217 return ElectronImpactIonisationVanOverstraetenDeMan(e, alpha);
218 break;
219 case ImpactIonisation::Grant:
220 return ElectronImpactIonisationGrant(e, alpha);
221 break;
222 case ImpactIonisation::Massey:
223 return ElectronImpactIonisationMassey(e, alpha);
224 break;
225 default:
226 std::cerr << m_className << "::ElectronTownsend: Unknown model. Bug!\n";
227 break;
228 }
229 return false;
230}
231
232bool MediumSilicon::ElectronAttachment(const double ex, const double ey,
233 const double ez, const double bx,
234 const double by, const double bz,
235 double& eta) {
236 eta = 0.;
237 if (m_isChanged) {
238 if (!UpdateTransportParameters()) {
239 std::cerr << m_className << "::ElectronAttachment:\n"
240 << " Error calculating the transport parameters.\n";
241 return false;
242 }
243 m_isChanged = false;
244 }
245
246 if (!m_eAtt.empty()) {
247 // Interpolation in user table.
248 return Medium::ElectronAttachment(ex, ey, ez, bx, by, bz, eta);
249 }
250
251 switch (m_trappingModel) {
252 case 0:
253 eta = m_eTrapCs * m_eTrapDensity;
254 break;
255 case 1:
256 double vx, vy, vz;
257 ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
258 eta = m_eTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
259 if (eta > 0.) eta = 1. / eta;
260 break;
261 default:
262 std::cerr << m_className << "::ElectronAttachment: Unknown model. Bug!\n";
263 return false;
264 break;
265 }
266
267 return true;
268}
269
270bool MediumSilicon::HoleVelocity(const double ex, const double ey,
271 const double ez, const double bx,
272 const double by, const double bz, double& vx,
273 double& vy, double& vz) {
274 vx = vy = vz = 0.;
275 if (m_isChanged) {
276 if (!UpdateTransportParameters()) {
277 std::cerr << m_className << "::HoleVelocity:\n"
278 << " Error calculating the transport parameters.\n";
279 return false;
280 }
281 m_isChanged = false;
282 }
283
284 if (!m_hVelE.empty()) {
285 // Interpolation in user table.
286 return Medium::HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
287 }
288
289 const double e = sqrt(ex * ex + ey * ey + ez * ez);
290
291 // Calculate the mobility
292 double mu;
293 switch (m_highFieldMobilityModel) {
294 case HighFieldMobility::Minimos:
295 HoleMobilityMinimos(e, mu);
296 break;
297 case HighFieldMobility::Canali:
298 HoleMobilityCanali(e, mu);
299 break;
300 case HighFieldMobility::Reggiani:
301 HoleMobilityReggiani(e, mu);
302 break;
303 default:
304 mu = m_hMobility;
305 }
306
307 const double b = sqrt(bx * bx + by * by + bz * bz);
308 if (b < Small) {
309 vx = mu * ex;
310 vy = mu * ey;
311 vz = mu * ez;
312 } else {
313 // Hall mobility
314 const double muH = m_hHallFactor * mu;
315 const double eb = bx * ex + by * ey + bz * ez;
316 const double mub = muH * b;
317 const double f = mu / (1. + mub * mub);
318 const double muH2 = muH * muH;
319 // Compute the drift velocity using the Langevin equation.
320 vx = f * (ex + muH * (ey * bz - ez * by) + muH2 * bx * eb);
321 vy = f * (ey + muH * (ez * bx - ex * bz) + muH2 * by * eb);
322 vz = f * (ez + muH * (ex * by - ey * bx) + muH2 * bz * eb);
323 }
324 return true;
325}
326
327bool MediumSilicon::HoleTownsend(const double ex, const double ey,
328 const double ez, const double bx,
329 const double by, const double bz,
330 double& alpha) {
331 alpha = 0.;
332 if (m_isChanged) {
333 if (!UpdateTransportParameters()) {
334 std::cerr << m_className << "::HoleTownsend:\n"
335 << " Error calculating the transport parameters.\n";
336 return false;
337 }
338 m_isChanged = false;
339 }
340
341 if (!m_hAlp.empty()) {
342 // Interpolation in user table.
343 return Medium::HoleTownsend(ex, ey, ez, bx, by, bz, alpha);
344 }
345
346 const double e = sqrt(ex * ex + ey * ey + ez * ez);
347
348 switch (m_impactIonisationModel) {
349 case ImpactIonisation::VanOverstraeten:
350 return HoleImpactIonisationVanOverstraetenDeMan(e, alpha);
351 break;
352 case ImpactIonisation::Grant:
353 return HoleImpactIonisationGrant(e, alpha);
354 break;
355 case ImpactIonisation::Massey:
356 return HoleImpactIonisationMassey(e, alpha);
357 break;
358 default:
359 std::cerr << m_className << "::HoleTownsend: Unknown model. Bug!\n";
360 break;
361 }
362 return false;
363}
364
365bool MediumSilicon::HoleAttachment(const double ex, const double ey,
366 const double ez, const double bx,
367 const double by, const double bz,
368 double& eta) {
369 eta = 0.;
370 if (m_isChanged) {
371 if (!UpdateTransportParameters()) {
372 std::cerr << m_className << "::HoleAttachment:\n"
373 << " Error calculating the transport parameters.\n";
374 return false;
375 }
376 m_isChanged = false;
377 }
378
379 if (!m_hAtt.empty()) {
380 // Interpolation in user table.
381 return Medium::HoleAttachment(ex, ey, ez, bx, by, bz, eta);
382 }
383
384 switch (m_trappingModel) {
385 case 0:
386 eta = m_hTrapCs * m_hTrapDensity;
387 break;
388 case 1:
389 double vx, vy, vz;
390 HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
391 eta = m_hTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
392 if (eta > 0.) eta = 1. / eta;
393 break;
394 default:
395 std::cerr << m_className << "::HoleAttachment: Unknown model. Bug!\n";
396 return false;
397 break;
398 }
399
400 return true;
401}
402
403void MediumSilicon::SetLowFieldMobility(const double mue, const double muh) {
404 if (mue <= 0. || muh <= 0.) {
405 std::cerr << m_className << "::SetLowFieldMobility:\n"
406 << " Mobility must be greater than zero.\n";
407 return;
408 }
409
410 m_eMobility = mue;
411 m_hMobility = muh;
412 m_hasUserMobility = true;
413 m_isChanged = true;
414}
415
417 m_latticeMobilityModel = LatticeMobility::Minimos;
418 m_hasUserMobility = false;
419 m_isChanged = true;
420}
421
423 m_latticeMobilityModel = LatticeMobility::Sentaurus;
424 m_hasUserMobility = false;
425 m_isChanged = true;
426}
427
429 m_latticeMobilityModel = LatticeMobility::Reggiani;
430 m_hasUserMobility = false;
431 m_isChanged = true;
432}
433
435 m_dopingMobilityModel = DopingMobility::Minimos;
436 m_hasUserMobility = false;
437 m_isChanged = true;
438}
439
441 m_dopingMobilityModel = DopingMobility::Masetti;
442 m_hasUserMobility = false;
443 m_isChanged = true;
444}
445
447 const double vsath) {
448 if (vsate <= 0. || vsath <= 0.) {
449 std::cout << m_className << "::SetSaturationVelocity:\n"
450 << " Restoring default values.\n";
451 m_hasUserSaturationVelocity = false;
452 } else {
453 m_eSatVel = vsate;
454 m_hSatVel = vsath;
455 m_hasUserSaturationVelocity = true;
456 }
457
458 m_isChanged = true;
459}
460
462 m_saturationVelocityModel = SaturationVelocity::Minimos;
463 m_hasUserSaturationVelocity = false;
464 m_isChanged = true;
465}
466
468 m_saturationVelocityModel = SaturationVelocity::Canali;
469 m_hasUserSaturationVelocity = false;
470 m_isChanged = true;
471}
472
474 m_saturationVelocityModel = SaturationVelocity::Reggiani;
475 m_hasUserSaturationVelocity = false;
476 m_isChanged = true;
477}
478
480 m_highFieldMobilityModel = HighFieldMobility::Minimos;
481 m_isChanged = true;
482}
483
485 m_highFieldMobilityModel = HighFieldMobility::Canali;
486 m_isChanged = true;
487}
488
490 m_highFieldMobilityModel = HighFieldMobility::Reggiani;
491 m_isChanged = true;
492}
493
495 m_highFieldMobilityModel = HighFieldMobility::Constant;
496}
497
499 m_impactIonisationModel = ImpactIonisation::VanOverstraeten;
500 m_isChanged = true;
501}
502
504 m_impactIonisationModel = ImpactIonisation::Grant;
505 m_isChanged = true;
506}
507
509 m_impactIonisationModel = ImpactIonisation::Massey;
510 m_isChanged = true;
511}
512
514 if (e <= m_eMinG + Small) {
515 std::cerr << m_className << "::SetMaxElectronEnergy:\n"
516 << " Requested upper electron energy limit (" << e
517 << " eV) is too small.\n";
518 return false;
519 }
520
521 m_eFinalG = e;
522 // Determine the energy interval size.
523 m_eStepG = m_eFinalG / nEnergyStepsG;
524
525 m_isChanged = true;
526
527 return true;
528}
529
530double MediumSilicon::GetElectronEnergy(const double px, const double py,
531 const double pz, double& vx, double& vy,
532 double& vz, const int band) {
533 // Effective masses
534 double mx = ElectronMass, my = ElectronMass, mz = ElectronMass;
535 // Energy offset
536 double e0 = 0.;
537 if (band >= 0 && band < m_nValleysX) {
538 // X valley
539 if (m_anisotropic) {
540 switch (band) {
541 case 0:
542 case 1:
543 // X 100, -100
544 mx *= m_mLongX;
545 my *= m_mTransX;
546 mz *= m_mTransX;
547 break;
548 case 2:
549 case 3:
550 // X 010, 0-10
551 mx *= m_mTransX;
552 my *= m_mLongX;
553 mz *= m_mTransX;
554 break;
555 case 4:
556 case 5:
557 // X 001, 00-1
558 mx *= m_mTransX;
559 my *= m_mTransX;
560 mz *= m_mLongX;
561 break;
562 default:
563 std::cerr << m_className << "::GetElectronEnergy:\n"
564 << " Unexpected band index " << band << "!\n";
565 break;
566 }
567 } else {
568 // Conduction effective mass
569 const double mc = 3. / (1. / m_mLongX + 2. / m_mTransX);
570 mx *= mc;
571 my *= mc;
572 mz *= mc;
573 }
574 } else if (band < m_nValleysX + m_nValleysL) {
575 // L valley, isotropic approximation
576 e0 = m_eMinL;
577 // Effective mass
578 const double mc = 3. / (1. / m_mLongL + 2. / m_mTransL);
579 mx *= mc;
580 my *= mc;
581 mz *= mc;
582 } else if (band == m_nValleysX + m_nValleysL) {
583 // Higher band(s)
584 }
585
586 if (m_nonParabolic) {
587 // Non-parabolicity parameter
588 double alpha = 0.;
589 if (band < m_nValleysX) {
590 // X valley
591 alpha = m_alphaX;
592 } else if (band < m_nValleysX + m_nValleysL) {
593 // L valley
594 alpha = m_alphaL;
595 }
596
597 const double p2 = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
598 if (alpha > 0.) {
599 const double e = 0.5 * (sqrt(1. + 4 * alpha * p2) - 1.) / alpha;
600 const double a = SpeedOfLight / (1. + 2 * alpha * e);
601 vx = a * px / mx;
602 vy = a * py / my;
603 vz = a * pz / mz;
604 return e0 + e;
605 }
606 }
607
608 const double e = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
609 vx = SpeedOfLight * px / mx;
610 vy = SpeedOfLight * py / my;
611 vz = SpeedOfLight * pz / mz;
612 return e0 + e;
613}
614
615void MediumSilicon::GetElectronMomentum(const double e, double& px, double& py,
616 double& pz, int& band) {
617 int nBands = m_nValleysX;
618 if (e > m_eMinL) nBands += m_nValleysL;
619 if (e > m_eMinG) ++nBands;
620
621 // If the band index is out of range, choose one at random.
622 if (band < 0 || band > m_nValleysX + m_nValleysL ||
623 (e < m_eMinL || band >= m_nValleysX) ||
624 (e < m_eMinG || band == m_nValleysX + m_nValleysL)) {
625 if (e < m_eMinL) {
626 band = int(m_nValleysX * RndmUniform());
627 if (band >= m_nValleysX) band = m_nValleysX - 1;
628 } else {
629 const double dosX = GetConductionBandDensityOfStates(e, 0);
630 const double dosL = GetConductionBandDensityOfStates(e, m_nValleysX);
631 const double dosG =
632 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
633 const double dosSum = m_nValleysX * dosX + m_nValleysL * dosL + dosG;
634 if (dosSum < Small) {
635 band = m_nValleysX + m_nValleysL;
636 } else {
637 const double r = RndmUniform() * dosSum;
638 if (r < dosX) {
639 band = int(m_nValleysX * RndmUniform());
640 if (band >= m_nValleysX) band = m_nValleysX - 1;
641 } else if (r < dosX + dosL) {
642 band = m_nValleysX + int(m_nValleysL * RndmUniform());
643 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysL - 1;
644 } else {
645 band = m_nValleysX + m_nValleysL;
646 }
647 }
648 }
649 if (m_debug) {
650 std::cout << m_className << "::GetElectronMomentum:\n"
651 << " Randomised band index: " << band << "\n";
652 }
653 }
654 if (band < m_nValleysX) {
655 // X valleys
656 double pstar = sqrt(2. * ElectronMass * e);
657 if (m_nonParabolic) {
658 const double alpha = m_alphaX;
659 pstar *= sqrt(1. + alpha * e);
660 }
661
662 const double ctheta = 1. - 2. * RndmUniform();
663 const double stheta = sqrt(1. - ctheta * ctheta);
664 const double phi = TwoPi * RndmUniform();
665
666 if (m_anisotropic) {
667 const double pl = pstar * sqrt(m_mLongX);
668 const double pt = pstar * sqrt(m_mTransX);
669 switch (band) {
670 case 0:
671 case 1:
672 // 100
673 px = pl * ctheta;
674 py = pt * cos(phi) * stheta;
675 pz = pt * sin(phi) * stheta;
676 break;
677 case 2:
678 case 3:
679 // 010
680 px = pt * sin(phi) * stheta;
681 py = pl * ctheta;
682 pz = pt * cos(phi) * stheta;
683 break;
684 case 4:
685 case 5:
686 // 001
687 px = pt * cos(phi) * stheta;
688 py = pt * sin(phi) * stheta;
689 pz = pl * ctheta;
690 break;
691 default:
692 // Other band; should not occur.
693 std::cerr << m_className << "::GetElectronMomentum:\n"
694 << " Unexpected band index (" << band << ").\n";
695 px = pstar * stheta * cos(phi);
696 py = pstar * stheta * sin(phi);
697 pz = pstar * ctheta;
698 break;
699 }
700 } else {
701 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
702 px = pstar * cos(phi) * stheta;
703 py = pstar * sin(phi) * stheta;
704 pz = pstar * ctheta;
705 }
706 } else if (band < m_nValleysX + m_nValleysL) {
707 // L valleys
708 double pstar = sqrt(2. * ElectronMass * (e - m_eMinL));
709 if (m_nonParabolic) {
710 const double alpha = m_alphaL;
711 pstar *= sqrt(1. + alpha * (e - m_eMinL));
712 }
713 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
714 RndmDirection(px, py, pz, pstar);
715 } else if (band == m_nValleysX + m_nValleysL) {
716 // Higher band
717 const double pstar = sqrt(2. * ElectronMass * e);
718 RndmDirection(px, py, pz, pstar);
719 }
720}
721
723 if (m_isChanged) {
724 if (!UpdateTransportParameters()) {
725 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
726 << " Error calculating the collision rates table.\n";
727 return 0.;
728 }
729 m_isChanged = false;
730 }
731
732 if (band >= 0 && band < m_nValleysX) {
733 return m_cfNullElectronsX;
734 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
735 return m_cfNullElectronsL;
736 } else if (band == m_nValleysX + m_nValleysL) {
737 return m_cfNullElectronsG;
738 }
739 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
740 << " Band index (" << band << ") out of range.\n";
741 return 0.;
742}
743
744double MediumSilicon::GetElectronCollisionRate(const double e, const int band) {
745 if (e <= 0.) {
746 std::cerr << m_className << "::GetElectronCollisionRate:\n"
747 << " Electron energy must be positive.\n";
748 return 0.;
749 }
750
751 if (e > m_eFinalG) {
752 std::cerr << m_className << "::GetElectronCollisionRate:\n"
753 << " Collision rate at " << e << " eV (band " << band
754 << ") is not included in the current table.\n"
755 << " Increasing energy range to " << 1.05 * e << " eV.\n";
756 SetMaxElectronEnergy(1.05 * e);
757 }
758
759 if (m_isChanged) {
760 if (!UpdateTransportParameters()) {
761 std::cerr << m_className << "::GetElectronCollisionRate:\n"
762 << " Error calculating the collision rates table.\n";
763 return 0.;
764 }
765 m_isChanged = false;
766 }
767
768 if (band >= 0 && band < m_nValleysX) {
769 int iE = int(e / m_eStepXL);
770 if (iE >= nEnergyStepsXL)
771 iE = nEnergyStepsXL - 1;
772 else if (iE < 0)
773 iE = 0;
774 return m_cfTotElectronsX[iE];
775 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
776 int iE = int(e / m_eStepXL);
777 if (iE >= nEnergyStepsXL)
778 iE = nEnergyStepsXL - 1;
779 else if (iE < m_ieMinL)
780 iE = m_ieMinL;
781 return m_cfTotElectronsL[iE];
782 } else if (band == m_nValleysX + m_nValleysL) {
783 int iE = int(e / m_eStepG);
784 if (iE >= nEnergyStepsG)
785 iE = nEnergyStepsG - 1;
786 else if (iE < m_ieMinG)
787 iE = m_ieMinG;
788 return m_cfTotElectronsG[iE];
789 }
790
791 std::cerr << m_className << "::GetElectronCollisionRate:\n"
792 << " Band index (" << band << ") out of range.\n";
793 return 0.;
794}
795
797 const double e, int& type, int& level, double& e1, double& px, double& py,
798 double& pz, std::vector<std::pair<int, double> >& secondaries, int& ndxc,
799 int& band) {
800 if (e > m_eFinalG) {
801 std::cerr << m_className << "::GetElectronCollision:\n"
802 << " Requested electron energy (" << e << " eV) exceeds the "
803 << "current energy range (" << m_eFinalG << " eV).\n"
804 << " Increasing energy range to " << 1.05 * e << " eV.\n";
805 SetMaxElectronEnergy(1.05 * e);
806 } else if (e <= 0.) {
807 std::cerr << m_className << "::GetElectronCollision:\n"
808 << " Electron energy must be greater than zero.\n";
809 return false;
810 }
811
812 if (m_isChanged) {
813 if (!UpdateTransportParameters()) {
814 std::cerr << m_className << "::GetElectronCollision:\n"
815 << " Error calculating the collision rates table.\n";
816 return false;
817 }
818 m_isChanged = false;
819 }
820
821 // Energy loss
822 double loss = 0.;
823 // Sample the scattering process.
824 if (band >= 0 && band < m_nValleysX) {
825 // X valley
826 // Get the energy interval.
827 int iE = int(e / m_eStepXL);
828 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
829 if (iE < 0) iE = 0;
830 // Select the scattering process.
831 const double r = RndmUniform();
832 if (r <= m_cfElectronsX[iE][0]) {
833 level = 0;
834 } else if (r >= m_cfElectronsX[iE][m_nLevelsX - 1]) {
835 level = m_nLevelsX - 1;
836 } else {
837 const auto begin = m_cfElectronsX[iE].cbegin();
838 level = std::lower_bound(begin, begin + m_nLevelsX, r) - begin;
839 }
840
841 // Get the collision type.
842 type = m_scatTypeElectronsX[level];
843 // Fill the collision counters.
844 ++m_nCollElectronDetailed[level];
845 ++m_nCollElectronBand[band];
846 if (type == ElectronCollisionTypeAcousticPhonon) {
847 ++m_nCollElectronAcoustic;
848 } else if (type == ElectronCollisionTypeIntervalleyG) {
849 // Intervalley scattering (g type)
850 ++m_nCollElectronIntervalley;
851 // Final valley is in opposite direction.
852 switch (band) {
853 case 0:
854 band = 1;
855 break;
856 case 1:
857 band = 0;
858 break;
859 case 2:
860 band = 3;
861 break;
862 case 3:
863 band = 2;
864 break;
865 case 4:
866 band = 5;
867 break;
868 case 5:
869 band = 4;
870 break;
871 default:
872 break;
873 }
874 } else if (type == ElectronCollisionTypeIntervalleyF) {
875 // Intervalley scattering (f type)
876 ++m_nCollElectronIntervalley;
877 // Final valley is perpendicular.
878 switch (band) {
879 case 0:
880 case 1:
881 band = int(RndmUniform() * 4) + 2;
882 break;
883 case 2:
884 case 3:
885 band = int(RndmUniform() * 4);
886 if (band > 1) band += 2;
887 break;
888 case 4:
889 case 5:
890 band = int(RndmUniform() * 4);
891 break;
892 }
893 } else if (type == ElectronCollisionTypeInterbandXL) {
894 // XL scattering
895 ++m_nCollElectronIntervalley;
896 // Final valley is in L band.
897 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
898 if (band >= m_nValleysX + m_nValleysL)
899 band = m_nValleysX + m_nValleysL - 1;
900 } else if (type == ElectronCollisionTypeInterbandXG) {
901 ++m_nCollElectronIntervalley;
902 band = m_nValleysX + m_nValleysL;
903 } else if (type == ElectronCollisionTypeImpurity) {
904 ++m_nCollElectronImpurity;
905 } else if (type == ElectronCollisionTypeIonisation) {
906 ++m_nCollElectronIonisation;
907 } else {
908 std::cerr << m_className << "::GetElectronCollision:\n";
909 std::cerr << " Unexpected collision type (" << type << ").\n";
910 }
911
912 // Get the energy loss.
913 loss = m_energyLossElectronsX[level];
914
915 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
916 // L valley
917 // Get the energy interval.
918 int iE = int(e / m_eStepXL);
919 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
920 if (iE < m_ieMinL) iE = m_ieMinL;
921 // Select the scattering process.
922 const double r = RndmUniform();
923 if (r <= m_cfElectronsL[iE][0]) {
924 level = 0;
925 } else if (r >= m_cfElectronsL[iE][m_nLevelsL - 1]) {
926 level = m_nLevelsL - 1;
927 } else {
928 const auto begin = m_cfElectronsL[iE].cbegin();
929 level = std::lower_bound(begin, begin + m_nLevelsL, r) - begin;
930 }
931
932 // Get the collision type.
933 type = m_scatTypeElectronsL[level];
934 // Fill the collision counters.
935 ++m_nCollElectronDetailed[m_nLevelsX + level];
936 ++m_nCollElectronBand[band];
937 if (type == ElectronCollisionTypeAcousticPhonon) {
938 ++m_nCollElectronAcoustic;
939 } else if (type == ElectronCollisionTypeOpticalPhonon) {
940 ++m_nCollElectronOptical;
941 } else if (type == ElectronCollisionTypeIntervalleyG ||
942 type == ElectronCollisionTypeIntervalleyF) {
943 // Equivalent intervalley scattering
944 ++m_nCollElectronIntervalley;
945 // Randomise the final valley.
946 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
947 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysX + m_nValleysL;
948 } else if (type == ElectronCollisionTypeInterbandXL) {
949 // LX scattering
950 ++m_nCollElectronIntervalley;
951 // Randomise the final valley.
952 band = int(RndmUniform() * m_nValleysX);
953 if (band >= m_nValleysX) band = m_nValleysX - 1;
954 } else if (type == ElectronCollisionTypeInterbandLG) {
955 // LG scattering
956 ++m_nCollElectronIntervalley;
957 band = m_nValleysX + m_nValleysL;
958 } else if (type == ElectronCollisionTypeImpurity) {
959 ++m_nCollElectronImpurity;
960 } else if (type == ElectronCollisionTypeIonisation) {
961 ++m_nCollElectronIonisation;
962 } else {
963 std::cerr << m_className << "::GetElectronCollision:\n";
964 std::cerr << " Unexpected collision type (" << type << ").\n";
965 }
966
967 // Get the energy loss.
968 loss = m_energyLossElectronsL[level];
969 } else if (band == m_nValleysX + m_nValleysL) {
970 // Higher bands
971 // Get the energy interval.
972 int iE = int(e / m_eStepG);
973 if (iE >= nEnergyStepsG) iE = nEnergyStepsG - 1;
974 if (iE < m_ieMinG) iE = m_ieMinG;
975 // Select the scattering process.
976 const double r = RndmUniform();
977 if (r <= m_cfElectronsG[iE][0]) {
978 level = 0;
979 } else if (r >= m_cfElectronsG[iE][m_nLevelsG - 1]) {
980 level = m_nLevelsG - 1;
981 } else {
982 const auto begin = m_cfElectronsG[iE].cbegin();
983 level = std::lower_bound(begin, begin + m_nLevelsG, r) - begin;
984 }
985
986 // Get the collision type.
987 type = m_scatTypeElectronsG[level];
988 // Fill the collision counters.
989 ++m_nCollElectronDetailed[m_nLevelsX + m_nLevelsL + level];
990 ++m_nCollElectronBand[band];
991 if (type == ElectronCollisionTypeAcousticPhonon) {
992 ++m_nCollElectronAcoustic;
993 } else if (type == ElectronCollisionTypeOpticalPhonon) {
994 ++m_nCollElectronOptical;
995 } else if (type == ElectronCollisionTypeIntervalleyG ||
996 type == ElectronCollisionTypeIntervalleyF) {
997 // Equivalent intervalley scattering
998 ++m_nCollElectronIntervalley;
999 } else if (type == ElectronCollisionTypeInterbandXG) {
1000 // GX scattering
1001 ++m_nCollElectronIntervalley;
1002 // Randomise the final valley.
1003 band = int(RndmUniform() * m_nValleysX);
1004 if (band >= m_nValleysX) band = m_nValleysX - 1;
1005 } else if (type == ElectronCollisionTypeInterbandLG) {
1006 // GL scattering
1007 ++m_nCollElectronIntervalley;
1008 // Randomise the final valley.
1009 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
1010 if (band >= m_nValleysX + m_nValleysL)
1011 band = m_nValleysX + m_nValleysL - 1;
1012 } else if (type == ElectronCollisionTypeIonisation) {
1013 ++m_nCollElectronIonisation;
1014 } else {
1015 std::cerr << m_className << "::GetElectronCollision:\n";
1016 std::cerr << " Unexpected collision type (" << type << ").\n";
1017 }
1018
1019 // Get the energy loss.
1020 loss = m_energyLossElectronsG[level];
1021 } else {
1022 std::cerr << m_className << "::GetElectronCollision:\n";
1023 std::cerr << " Band index (" << band << ") out of range.\n";
1024 return false;
1025 }
1026
1027 // Secondaries
1028 ndxc = 0;
1029 // Ionising collision
1030 if (type == ElectronCollisionTypeIonisation) {
1031 double ee = 0., eh = 0.;
1032 ComputeSecondaries(e, ee, eh);
1033 loss = ee + eh + m_bandGap;
1034 // Add the secondary electron.
1035 secondaries.emplace_back(std::make_pair(IonProdTypeElectron, ee));
1036 // Add the hole.
1037 secondaries.emplace_back(std::make_pair(IonProdTypeHole, eh));
1038 }
1039
1040 if (e < loss) loss = e - 0.00001;
1041 // Update the energy.
1042 e1 = e - loss;
1043 if (e1 < Small) e1 = Small;
1044
1045 // Update the momentum.
1046 if (band >= 0 && band < m_nValleysX) {
1047 // X valleys
1048 double pstar = sqrt(2. * ElectronMass * e1);
1049 if (m_nonParabolic) {
1050 const double alpha = m_alphaX;
1051 pstar *= sqrt(1. + alpha * e1);
1052 }
1053
1054 const double ctheta = 1. - 2. * RndmUniform();
1055 const double stheta = sqrt(1. - ctheta * ctheta);
1056 const double phi = TwoPi * RndmUniform();
1057
1058 if (m_anisotropic) {
1059 const double pl = pstar * sqrt(m_mLongX);
1060 const double pt = pstar * sqrt(m_mTransX);
1061 switch (band) {
1062 case 0:
1063 case 1:
1064 // 100
1065 px = pl * ctheta;
1066 py = pt * cos(phi) * stheta;
1067 pz = pt * sin(phi) * stheta;
1068 break;
1069 case 2:
1070 case 3:
1071 // 010
1072 px = pt * sin(phi) * stheta;
1073 py = pl * ctheta;
1074 pz = pt * cos(phi) * stheta;
1075 break;
1076 case 4:
1077 case 5:
1078 // 001
1079 px = pt * cos(phi) * stheta;
1080 py = pt * sin(phi) * stheta;
1081 pz = pl * ctheta;
1082 break;
1083 default:
1084 return false;
1085 break;
1086 }
1087 } else {
1088 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
1089 px = pstar * cos(phi) * stheta;
1090 py = pstar * sin(phi) * stheta;
1091 pz = pstar * ctheta;
1092 }
1093 return true;
1094
1095 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
1096 // L valleys
1097 double pstar = sqrt(2. * ElectronMass * (e1 - m_eMinL));
1098 if (m_nonParabolic) {
1099 const double alpha = m_alphaL;
1100 pstar *= sqrt(1. + alpha * (e1 - m_eMinL));
1101 }
1102 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
1103 RndmDirection(px, py, pz, pstar);
1104 return true;
1105 } else {
1106 const double pstar = sqrt(2. * ElectronMass * e1);
1107 RndmDirection(px, py, pz, pstar);
1108 return true;
1109 }
1110
1111 std::cerr << m_className << "::GetElectronCollision:"
1112 << " Band index (" << band << ") out of range.\n";
1113 e1 = e;
1114 type = 0;
1115 return false;
1116}
1117
1119 m_nCollElectronAcoustic = m_nCollElectronOptical = 0;
1120 m_nCollElectronIntervalley = 0;
1121 m_nCollElectronImpurity = 0;
1122 m_nCollElectronIonisation = 0;
1123 const auto nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1124 m_nCollElectronDetailed.assign(nLevels, 0);
1125 const auto nBands = m_nValleysX + m_nValleysL + 1;
1126 m_nCollElectronBand.assign(nBands, 0);
1127}
1128
1130 return m_nCollElectronAcoustic + m_nCollElectronOptical +
1131 m_nCollElectronIntervalley + m_nCollElectronImpurity +
1132 m_nCollElectronIonisation;
1133}
1134
1136 return m_nLevelsX + m_nLevelsL + m_nLevelsG;
1137}
1138
1140 const unsigned int level) const {
1141 const unsigned int nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1142 if (level >= nLevels) {
1143 std::cerr << m_className << "::GetNumberOfElectronCollisions:\n"
1144 << " Scattering rate term (" << level << ") does not exist.\n";
1145 return 0;
1146 }
1147 return m_nCollElectronDetailed[level];
1148}
1149
1151 return m_nValleysX + m_nValleysL + 1;
1152}
1153
1155 const int nBands = m_nValleysX + m_nValleysL + 1;
1156 if (band < 0 || band >= nBands) {
1157 std::cerr << m_className << "::GetElectronBandPopulation:\n";
1158 std::cerr << " Band index (" << band << ") out of range.\n";
1159 return 0;
1160 }
1161 return m_nCollElectronBand[band];
1162}
1163
1164bool MediumSilicon::GetOpticalDataRange(double& emin, double& emax,
1165 const unsigned int i) {
1166 if (i != 0) {
1167 std::cerr << m_className << "::GetOpticalDataRange: Index out of range.\n";
1168 }
1169
1170 // Make sure the optical data table has been loaded.
1171 if (m_opticalDataEnergies.empty()) {
1172 if (!LoadOpticalData(m_opticalDataFile)) {
1173 std::cerr << m_className << "::GetOpticalDataRange:\n"
1174 << " Optical data table could not be loaded.\n";
1175 return false;
1176 }
1177 }
1178
1179 emin = m_opticalDataEnergies.front();
1180 emax = m_opticalDataEnergies.back();
1181 if (m_debug) {
1182 std::cout << m_className << "::GetOpticalDataRange:\n"
1183 << " " << emin << " < E [eV] < " << emax << "\n";
1184 }
1185 return true;
1186}
1187
1188bool MediumSilicon::GetDielectricFunction(const double e, double& eps1,
1189 double& eps2, const unsigned int i) {
1190 if (i != 0) {
1191 std::cerr << m_className + "::GetDielectricFunction: Index out of range.\n";
1192 return false;
1193 }
1194
1195 // Make sure the optical data table has been loaded.
1196 if (m_opticalDataEnergies.empty()) {
1197 if (!LoadOpticalData(m_opticalDataFile)) {
1198 std::cerr << m_className << "::GetDielectricFunction:\n";
1199 std::cerr << " Optical data table could not be loaded.\n";
1200 return false;
1201 }
1202 }
1203
1204 // Make sure the requested energy is within the range of the table.
1205 const double emin = m_opticalDataEnergies.front();
1206 const double emax = m_opticalDataEnergies.back();
1207 if (e < emin || e > emax) {
1208 std::cerr << m_className << "::GetDielectricFunction:\n"
1209 << " Requested energy (" << e << " eV) "
1210 << " is outside the range of the optical data table.\n"
1211 << " " << emin << " < E [eV] < " << emax << "\n";
1212 eps1 = eps2 = 0.;
1213 return false;
1214 }
1215
1216 // Locate the requested energy in the table.
1217 const auto begin = m_opticalDataEnergies.cbegin();
1218 const auto it1 = std::upper_bound(begin, m_opticalDataEnergies.cend(), e);
1219 if (it1 == begin) {
1220 eps1 = m_opticalDataEpsilon.front().first;
1221 eps2 = m_opticalDataEpsilon.front().second;
1222 return true;
1223 }
1224 const auto it0 = std::prev(it1);
1225
1226 // Interpolate the real part of dielectric function.
1227 const double x0 = *it0;
1228 const double x1 = *it1;
1229 const double lnx0 = log(*it0);
1230 const double lnx1 = log(*it1);
1231 const double lnx = log(e);
1232 const double y0 = m_opticalDataEpsilon[it0 - begin].first;
1233 const double y1 = m_opticalDataEpsilon[it1 - begin].first;
1234 if (y0 <= 0. || y1 <= 0.) {
1235 // Use linear interpolation if one of the values is negative.
1236 eps1 = y0 + (e - x0) * (y1 - y0) / (x1 - x0);
1237 } else {
1238 // Otherwise use log-log interpolation.
1239 const double lny0 = log(y0);
1240 const double lny1 = log(y1);
1241 eps1 = lny0 + (lnx - lnx0) * (lny1 - lny0) / (lnx1 - lnx0);
1242 eps1 = exp(eps1);
1243 }
1244
1245 // Interpolate the imaginary part of dielectric function,
1246 // using log-log interpolation.
1247 const double lnz0 = log(m_opticalDataEpsilon[it0 - begin].second);
1248 const double lnz1 = log(m_opticalDataEpsilon[it1 - begin].second);
1249 eps2 = lnz0 + (lnx - lnx0) * (lnz1 - lnz0) / (lnx1 - lnx0);
1250 eps2 = exp(eps2);
1251 return true;
1252}
1253
1255 if (!m_isChanged) {
1256 if (m_debug) {
1257 std::cerr << m_className << "::Initialise: Nothing changed.\n";
1258 }
1259 return true;
1260 }
1261 if (!UpdateTransportParameters()) {
1262 std::cerr << m_className << "::Initialise:\n Error preparing "
1263 << "transport parameters/calculating collision rates.\n";
1264 return false;
1265 }
1266 return true;
1267}
1268
1269bool MediumSilicon::UpdateTransportParameters() {
1270 // Calculate impact ionisation coefficients
1271 switch (m_impactIonisationModel) {
1272 case ImpactIonisation::VanOverstraeten:
1273 UpdateImpactIonisationVanOverstraetenDeMan();
1274 break;
1275 case ImpactIonisation::Grant:
1276 UpdateImpactIonisationGrant();
1277 break;
1278 case ImpactIonisation::Massey:
1279 break;
1280 default:
1281 std::cerr << m_className << "::UpdateTransportParameters:\n "
1282 << "Unknown impact ionisation model. Program bug!\n";
1283 break;
1284 }
1285
1286 if (!m_hasUserMobility) {
1287 // Calculate lattice mobility
1288 switch (m_latticeMobilityModel) {
1289 case LatticeMobility::Minimos:
1290 UpdateLatticeMobilityMinimos();
1291 break;
1292 case LatticeMobility::Sentaurus:
1293 UpdateLatticeMobilitySentaurus();
1294 break;
1295 case LatticeMobility::Reggiani:
1296 UpdateLatticeMobilityReggiani();
1297 break;
1298 default:
1299 std::cerr << m_className << "::UpdateTransportParameters:\n "
1300 << "Unknown lattice mobility model. Program bug!\n";
1301 break;
1302 }
1303
1304 // Calculate doping mobility
1305 switch (m_dopingMobilityModel) {
1306 case DopingMobility::Minimos:
1307 UpdateDopingMobilityMinimos();
1308 break;
1309 case DopingMobility::Masetti:
1310 UpdateDopingMobilityMasetti();
1311 break;
1312 default:
1313 std::cerr << m_className << "::UpdateTransportParameters:\n "
1314 << "Unknown doping mobility model. Program bug!\n";
1315 break;
1316 }
1317 }
1318
1319 // Calculate saturation velocity
1320 if (!m_hasUserSaturationVelocity) {
1321 switch (m_saturationVelocityModel) {
1322 case SaturationVelocity::Minimos:
1323 UpdateSaturationVelocityMinimos();
1324 break;
1325 case SaturationVelocity::Canali:
1326 UpdateSaturationVelocityCanali();
1327 break;
1328 case SaturationVelocity::Reggiani:
1329 UpdateSaturationVelocityReggiani();
1330 break;
1331 }
1332 }
1333
1334 // Calculate high field saturation parameters
1335 if (m_highFieldMobilityModel == HighFieldMobility::Canali) {
1336 UpdateHighFieldMobilityCanali();
1337 }
1338
1339 if (m_debug) {
1340 std::cout << m_className << "::UpdateTransportParameters:\n"
1341 << " Low-field mobility [cm2 V-1 ns-1]\n"
1342 << " Electrons: " << m_eMobility << "\n"
1343 << " Holes: " << m_hMobility << "\n";
1344 if (m_highFieldMobilityModel == HighFieldMobility::Constant) {
1345 std::cout << " Mobility is not field-dependent.\n";
1346 } else {
1347 std::cout << " Saturation velocity [cm / ns]\n"
1348 << " Electrons: " << m_eSatVel << "\n"
1349 << " Holes: " << m_hSatVel << "\n";
1350 }
1351 }
1352
1353 if (!ElectronScatteringRates()) return false;
1354 if (!HoleScatteringRates()) return false;
1355
1356 // Reset the collision counters.
1358
1359 return true;
1360}
1361
1362void MediumSilicon::UpdateLatticeMobilityMinimos() {
1363 // References:
1364 // - S. Selberherr, W. Haensch, M. Seavey, J. Slotboom,
1365 // Solid State Electronics 33 (1990), 1425-1436
1366 // - Minimos 6.1 User's Guide (1999)
1367
1368 // Lattice mobilities at 300 K [cm2 / (V ns)]
1369 constexpr double eMu0 = 1.43e-6;
1370 constexpr double hMu0 = 0.46e-6;
1371 // Temperature normalized to 300 K
1372 const double t = m_temperature / 300.;
1373 // Temperature dependence of lattice mobility
1374 m_eLatticeMobility = eMu0 * pow(t, -2.);
1375 m_hLatticeMobility = hMu0 * pow(t, -2.18);
1376}
1377
1378void MediumSilicon::UpdateLatticeMobilitySentaurus() {
1379 // References:
1380 // - C. Lombardi et al.,
1381 // IEEE Trans. CAD 7 (1988), 1164-1171
1382 // - Sentaurus Device User Guide (2007)
1383
1384 // Lattice mobilities at 300 K [cm2 / (V ns)]
1385 constexpr double eMu0 = 1.417e-6;
1386 constexpr double hMu0 = 0.4705e-6;
1387 // Temperature normalized to 300 K
1388 const double t = m_temperature / 300.;
1389 // Temperature dependence of lattice mobility
1390 m_eLatticeMobility = eMu0 * pow(t, -2.5);
1391 m_hLatticeMobility = hMu0 * pow(t, -2.2);
1392}
1393
1394void MediumSilicon::UpdateLatticeMobilityReggiani() {
1395 // Reference:
1396 // - M. A. Omar, L. Reggiani
1397 // Solid State Electronics 30 (1987), 693-697
1398
1399 // Lattice mobilities at 300 K [cm2 / (V ns)]
1400 constexpr double eMu0 = 1.320e-6;
1401 constexpr double hMu0 = 0.460e-6;
1402 // Temperature normalized to 300 K
1403 const double t = m_temperature / 300.;
1404 // Temperature dependence of lattice mobility
1405 m_eLatticeMobility = eMu0 * pow(t, -2.);
1406 m_hLatticeMobility = hMu0 * pow(t, -2.2);
1407}
1408
1409void MediumSilicon::UpdateDopingMobilityMinimos() {
1410 // References:
1411 // - S. Selberherr, W. Haensch, M. Seavey, J. Slotboom,
1412 // Solid State Electronics 33 (1990), 1425-1436
1413 // - Minimos 6.1 User's Guide (1999)
1414
1415 // Mobility reduction due to ionised impurity scattering
1416 // Surface term not taken into account
1417 double eMuMin = 0.080e-6;
1418 double hMuMin = 0.045e-6;
1419 if (m_temperature > 200.) {
1420 const double c0 = pow(m_temperature / 300., -0.45);
1421 eMuMin *= c0;
1422 hMuMin *= c0;
1423 } else {
1424 const double c0 = pow(2. / 3., -0.45) * pow(m_temperature / 200., -0.15);
1425 eMuMin *= c0;
1426 hMuMin *= c0;
1427 }
1428 const double c1 = pow(m_temperature / 300., 3.2);
1429 const double eRefC = 1.12e17 * c1;
1430 const double hRefC = 2.23e17 * c1;
1431 const double alpha = 0.72 * pow(m_temperature / 300., 0.065);
1432 // Assume impurity concentration equal to doping concentration
1433 m_eMobility = eMuMin +
1434 (m_eLatticeMobility - eMuMin) /
1435 (1. + pow(m_dopingConcentration / eRefC, alpha));
1436 m_hMobility = hMuMin +
1437 (m_hLatticeMobility - hMuMin) /
1438 (1. + pow(m_dopingConcentration / hRefC, alpha));
1439}
1440
1441void MediumSilicon::UpdateDopingMobilityMasetti() {
1442 // Reference:
1443 // - G. Masetti, M. Severi, S. Solmi,
1444 // IEEE Trans. Electron Devices 30 (1983), 764-769
1445 // - Sentaurus Device User Guide (2007)
1446 // - Minimos NT User Guide (2004)
1447
1448 if (m_dopingConcentration < 1.e13) {
1449 m_eMobility = m_eLatticeMobility;
1450 m_hMobility = m_hLatticeMobility;
1451 return;
1452 }
1453
1454 // Parameters adopted from Minimos NT User Guide
1455 constexpr double eMuMin1 = 0.0522e-6;
1456 constexpr double eMuMin2 = 0.0522e-6;
1457 constexpr double eMu1 = 0.0434e-6;
1458 constexpr double hMuMin1 = 0.0449e-6;
1459 constexpr double hMuMin2 = 0.;
1460 constexpr double hMu1 = 0.029e-6;
1461 constexpr double eCr = 9.68e16;
1462 constexpr double eCs = 3.42e20;
1463 constexpr double hCr = 2.23e17;
1464 constexpr double hCs = 6.10e20;
1465 constexpr double hPc = 9.23e16;
1466 constexpr double eAlpha = 0.68;
1467 constexpr double eBeta = 2.;
1468 constexpr double hAlpha = 0.719;
1469 constexpr double hBeta = 2.;
1470
1471 m_eMobility = eMuMin1 +
1472 (m_eLatticeMobility - eMuMin2) /
1473 (1. + pow(m_dopingConcentration / eCr, eAlpha)) -
1474 eMu1 / (1. + pow(eCs / m_dopingConcentration, eBeta));
1475
1476 m_hMobility = hMuMin1 * exp(-hPc / m_dopingConcentration) +
1477 (m_hLatticeMobility - hMuMin2) /
1478 (1. + pow(m_dopingConcentration / hCr, hAlpha)) -
1479 hMu1 / (1. + pow(hCs / m_dopingConcentration, hBeta));
1480}
1481
1482void MediumSilicon::UpdateSaturationVelocityMinimos() {
1483 // References:
1484 // - R. Quay, C. Moglestue, V. Palankovski, S. Selberherr,
1485 // Materials Science in Semiconductor Processing 3 (2000), 149-155
1486 // - Minimos NT User Guide (2004)
1487
1488 // Temperature-dependence of saturation velocities [cm / ns]
1489 m_eSatVel = 1.e-2 / (1. + 0.74 * (m_temperature / 300. - 1.));
1490 m_hSatVel = 0.704e-2 / (1. + 0.37 * (m_temperature / 300. - 1.));
1491}
1492
1493void MediumSilicon::UpdateSaturationVelocityCanali() {
1494 // References:
1495 // - C. Canali, G. Majni, R. Minder, G. Ottaviani,
1496 // IEEE Transactions on Electron Devices 22 (1975), 1045-1047
1497 // - Sentaurus Device User Guide (2007)
1498
1499 m_eSatVel = 1.07e-2 * pow(300. / m_temperature, 0.87);
1500 m_hSatVel = 8.37e-3 * pow(300. / m_temperature, 0.52);
1501}
1502
1503void MediumSilicon::UpdateSaturationVelocityReggiani() {
1504 // Reference:
1505 // - M. A. Omar, L. Reggiani
1506 // Solid State Electronics 30 (1987), 693-697
1507
1508 m_eSatVel = 1.470e-2 * sqrt(tanh(150. / m_temperature));
1509 m_hSatVel = 0.916e-2 * sqrt(tanh(300. / m_temperature));
1510}
1511
1512void MediumSilicon::UpdateHighFieldMobilityCanali() {
1513 // References:
1514 // - C. Canali, G. Majni, R. Minder, G. Ottaviani,
1515 // IEEE Transactions on Electron Devices 22 (1975), 1045-1047
1516 // - Sentaurus Device User Guide (2007)
1517
1518 // Temperature dependent exponent in high-field mobility formula
1519 m_eBetaCanali = 1.109 * pow(m_temperature / 300., 0.66);
1520 m_hBetaCanali = 1.213 * pow(m_temperature / 300., 0.17);
1521 m_eBetaCanaliInv = 1. / m_eBetaCanali;
1522 m_hBetaCanaliInv = 1. / m_hBetaCanali;
1523}
1524
1525void MediumSilicon::UpdateImpactIonisationVanOverstraetenDeMan() {
1526
1527 // References:
1528 // - R. van Overstraeten and H. de Man,
1529 // Solid State Electronics 13 (1970), 583-608
1530 // - W. Maes, K. de Meyer and R. van Overstraeten,
1531 // Solid State Electronics 33 (1990), 705-718
1532 // - Sentaurus Device User Guide (2016)
1533
1534 // Temperature dependence as in Sentaurus Device
1535 // Optical phonon energy
1536 constexpr double hbarOmega = 0.063;
1537 // Temperature scaling coefficient
1538 const double gamma =
1539 tanh(hbarOmega / (2. * BoltzmannConstant * 300.)) /
1540 tanh(hbarOmega / (2. * BoltzmannConstant * m_temperature));
1541
1542 // Low field coefficients taken from Maes, de Meyer, van Overstraeten
1543 // eImpactA0 = gamma * 3.318e5;
1544 // eImpactB0 = gamma * 1.135e6;
1545 m_eImpactA0 = gamma * 7.03e5;
1546 m_eImpactB0 = gamma * 1.231e6;
1547 m_eImpactA1 = gamma * 7.03e5;
1548 m_eImpactB1 = gamma * 1.231e6;
1549
1550 m_hImpactA0 = gamma * 1.582e6;
1551 m_hImpactB0 = gamma * 2.036e6;
1552 m_hImpactA1 = gamma * 6.71e5;
1553 m_hImpactB1 = gamma * 1.693e6;
1554}
1555
1556void MediumSilicon::UpdateImpactIonisationGrant() {
1557 // References:
1558 // - W. N. Grant,
1559 // Solid State Electronics 16 (1973), 1189 - 1203
1560 // - Sentaurus Device User Guide (2007)
1561
1562 // Temperature dependence as in Sentaurus Device
1563 // Optical phonon energy
1564 constexpr double hbarOmega = 0.063;
1565 // Temperature scaling coefficient
1566 const double gamma =
1567 tanh(hbarOmega / (2. * BoltzmannConstant * 300.)) /
1568 tanh(hbarOmega / (2. * BoltzmannConstant * m_temperature));
1569
1570 m_eImpactA0 = 2.60e6 * gamma;
1571 m_eImpactB0 = 1.43e6 * gamma;
1572 m_eImpactA1 = 6.20e5 * gamma;
1573 m_eImpactB1 = 1.08e6 * gamma;
1574 m_eImpactA2 = 5.05e5 * gamma;
1575 m_eImpactB2 = 9.90e5 * gamma;
1576
1577 m_hImpactA0 = 2.00e6 * gamma;
1578 m_hImpactB0 = 1.97e6 * gamma;
1579 m_hImpactA1 = 5.60e5 * gamma;
1580 m_hImpactB1 = 1.32e6 * gamma;
1581}
1582
1583bool MediumSilicon::ElectronMobilityMinimos(const double e, double& mu) const {
1584 // Reference:
1585 // - Minimos User's Guide (1999)
1586
1587 if (e < Small) {
1588 mu = 0.;
1589 } else {
1590 const double r = 2 * m_eMobility * e / m_eSatVel;
1591 mu = 2. * m_eMobility / (1. + sqrt(1. + r * r));
1592 }
1593 return true;
1594}
1595
1596bool MediumSilicon::ElectronMobilityCanali(const double e, double& mu) const {
1597 // Reference:
1598 // - Sentaurus Device User Guide (2007)
1599
1600 if (e < Small) {
1601 mu = 0.;
1602 } else {
1603 const double r = m_eMobility * e / m_eSatVel;
1604 mu = m_eMobility / pow(1. + pow(r, m_eBetaCanali), m_eBetaCanaliInv);
1605 }
1606 return true;
1607}
1608
1609bool MediumSilicon::ElectronMobilityReggiani(const double e, double& mu) const {
1610 // Reference:
1611 // - M. A. Omar, L. Reggiani
1612 // Solid State Electronics 30 (1987), 693-697
1613
1614 if (e < Small) {
1615 mu = 0.;
1616 } else {
1617 const double r = m_eMobility * e / m_eSatVel;
1618 constexpr double k = 1. / 1.5;
1619 mu = m_eMobility / pow(1. + pow(r, 1.5), k);
1620 }
1621 return true;
1622}
1623
1624bool MediumSilicon::ElectronImpactIonisationVanOverstraetenDeMan(
1625 const double e, double& alpha) const {
1626 // References:
1627 // - R. van Overstraeten and H. de Man,
1628 // Solid State Electronics 13 (1970), 583-608
1629 // - W. Maes, K. de Meyer and R. van Overstraeten,
1630 // Solid State Electronics 33 (1990), 705-718
1631 // - Sentaurus Device User Guide (2016)
1632
1633 if (e < Small) {
1634 alpha = 0.;
1635 } else if (e < 4e5) {
1636 alpha = m_eImpactA0 * exp(-m_eImpactB0 / e);
1637 } else {
1638 alpha = m_eImpactA1 * exp(-m_eImpactB1 / e);
1639 }
1640 return true;
1641}
1642
1643bool MediumSilicon::ElectronImpactIonisationGrant(const double e,
1644 double& alpha) const {
1645 // Reference:
1646 // - W. N. Grant, Solid State Electronics 16 (1973), 1189 - 1203
1647
1648 if (e < Small) {
1649 alpha = 0.;
1650 } else if (e < 2.4e5) {
1651 alpha = m_eImpactA0 * exp(-m_eImpactB0 / e);
1652 } else if (e < 5.3e5) {
1653 alpha = m_eImpactA1 * exp(-m_eImpactB1 / e);
1654 } else {
1655 alpha = m_eImpactA2 * exp(-m_eImpactB2 / e);
1656 }
1657 return true;
1658}
1659
1660bool MediumSilicon::ElectronImpactIonisationMassey(const double e,
1661 double& alpha) const {
1662
1663 // Reference:
1664 // - D. J. Massey, J. P. R. David, and G. J. Rees,
1665 // IEEE Trans. Electron Devices 53 (2006), 2328 - 2334
1666 constexpr double a = 4.43e5;
1667 constexpr double c = 9.66e5;
1668 constexpr double d = 4.99e2;
1669 const double b = c + d * m_temperature;
1670 if (e < Small) {
1671 alpha = 0.;
1672 } else {
1673 alpha = a * exp(-b / e);
1674 }
1675 return true;
1676}
1677
1678bool MediumSilicon::HoleMobilityMinimos(const double e, double& mu) const {
1679 // Reference:
1680 // - Minimos User's Guide (1999)
1681
1682 if (e < Small) {
1683 mu = 0.;
1684 } else {
1685 mu = m_hMobility / (1. + m_hMobility * e / m_eSatVel);
1686 }
1687 return true;
1688}
1689
1690bool MediumSilicon::HoleMobilityCanali(const double e, double& mu) const {
1691 // Reference:
1692 // - Sentaurus Device User Guide (2007)
1693
1694 if (e < Small) {
1695 mu = 0.;
1696 } else {
1697 const double r = m_hMobility * e / m_hSatVel;
1698 mu = m_hMobility / pow(1. + pow(r, m_hBetaCanali), m_hBetaCanaliInv);
1699 }
1700 return true;
1701}
1702
1703bool MediumSilicon::HoleMobilityReggiani(const double e, double& mu) const {
1704 // Reference:
1705 // - M. A. Omar, L. Reggiani
1706 // Solid State Electronics 30 (1987), 693-697
1707
1708 if (e < Small) {
1709 mu = 0.;
1710 } else {
1711 const double r = m_hMobility * e / m_hSatVel;
1712 mu = m_hMobility / sqrt(1. + r * r);
1713 }
1714 return true;
1715}
1716
1717bool MediumSilicon::HoleImpactIonisationVanOverstraetenDeMan(
1718 const double e, double& alpha) const {
1719 // References:
1720 // - R. van Overstraeten and H. de Man,
1721 // Solid State Electronics 13 (1970), 583-608
1722 // - Sentaurus Device User Guide (2016)
1723
1724 if (e < Small) {
1725 alpha = 0.;
1726 } else if (e < 4e5) {
1727 alpha = m_hImpactA0 * exp(-m_hImpactB0 / e);
1728 } else {
1729 alpha = m_hImpactA1 * exp(-m_hImpactB1 / e);
1730 }
1731 return true;
1732}
1733
1734bool MediumSilicon::HoleImpactIonisationGrant(const double e,
1735 double& alpha) const {
1736 // Reference:
1737 // - W. N. Grant, Solid State Electronics 16 (1973), 1189 - 1203
1738
1739 if (e < Small) {
1740 alpha = 0.;
1741 } else if (e < 5.3e5) {
1742 alpha = m_hImpactA0 * exp(-m_hImpactB0 / e);
1743 } else {
1744 alpha = m_hImpactA1 * exp(-m_hImpactB1 / e);
1745 }
1746 return true;
1747}
1748
1749bool MediumSilicon::HoleImpactIonisationMassey(const double e,
1750 double& alpha) const {
1751
1752 // Reference:
1753 // - D. J. Massey, J. P. R. David, and G. J. Rees,
1754 // IEEE Trans. Electron Devices 53 (2006), 2328 - 2334
1755 constexpr double a = 1.13e6;
1756 constexpr double c = 1.17e6;
1757 constexpr double d = 1.09e3;
1758 const double b = c + d * m_temperature;
1759 if (e < Small) {
1760 alpha = 0.;
1761 } else {
1762 alpha = a * exp(-b / e);
1763 }
1764 return true;
1765}
1766
1767bool MediumSilicon::LoadOpticalData(const std::string& filename) {
1768 // Clear the optical data table.
1769 m_opticalDataEnergies.clear();
1770 m_opticalDataEpsilon.clear();
1771
1772 // Get the path to the data directory.
1773 char* pPath = getenv("GARFIELD_HOME");
1774 if (pPath == 0) {
1775 std::cerr << m_className << "::LoadOpticalData:\n"
1776 << " Environment variable GARFIELD_HOME is not set.\n";
1777 return false;
1778 }
1779 const std::string filepath = std::string(pPath) + "/Data/" + filename;
1780
1781 // Open the file.
1782 std::ifstream infile;
1783 infile.open(filepath.c_str(), std::ios::in);
1784 // Make sure the file could actually be opened.
1785 if (!infile) {
1786 std::cerr << m_className << "::LoadOpticalData:\n"
1787 << " Error opening file " << filename << ".\n";
1788 return false;
1789 }
1790
1791 double lastEnergy = -1.;
1792 double energy, eps1, eps2, loss;
1793 // Read the file line by line.
1794 std::string line;
1795 std::istringstream dataStream;
1796 int i = 0;
1797 while (!infile.eof()) {
1798 ++i;
1799 // Read the next line.
1800 std::getline(infile, line);
1801 // Strip white space from the beginning of the line.
1802 ltrim(line);
1803 if (line.empty()) continue;
1804 // Skip comments.
1805 if (line[0] == '#' || line[0] == '*' || (line[0] == '/' && line[1] == '/'))
1806 continue;
1807 // Extract the values.
1808 dataStream.str(line);
1809 dataStream >> energy >> eps1 >> eps2 >> loss;
1810 if (dataStream.eof()) break;
1811 // Check if the data has been read correctly.
1812 if (infile.fail()) {
1813 std::cerr << m_className << "::LoadOpticalData:\n Error reading file "
1814 << filename << " (line " << i << ").\n";
1815 return false;
1816 }
1817 // Reset the stringstream.
1818 dataStream.str("");
1819 dataStream.clear();
1820 // Make sure the values make sense.
1821 // The table has to be in ascending order
1822 // with respect to the photon energy.
1823 if (energy <= lastEnergy) {
1824 std::cerr << m_className << "::LoadOpticalData:\n Table is not in "
1825 << "monotonically increasing order (line " << i << ").\n"
1826 << " " << lastEnergy << " " << energy << " " << eps1
1827 << " " << eps2 << "\n";
1828 return false;
1829 }
1830 // The imaginary part of the dielectric function has to be positive.
1831 if (eps2 < 0.) {
1832 std::cerr << m_className << "::LoadOpticalData:\n Negative value "
1833 << "of the loss function (line " << i << ").\n";
1834 return false;
1835 }
1836 // Ignore negative photon energies.
1837 if (energy <= 0.) continue;
1838 // Add the values to the list.
1839 m_opticalDataEnergies.emplace_back(energy);
1840 m_opticalDataEpsilon.emplace_back(std::make_pair(eps1, eps2));
1841 lastEnergy = energy;
1842 }
1843
1844 if (m_opticalDataEnergies.empty()) {
1845 std::cerr << m_className << "::LoadOpticalData:\n"
1846 << " Import of data from file " << filepath << "failed.\n"
1847 << " No valid data found.\n";
1848 return false;
1849 }
1850
1851 if (m_debug) {
1852 std::cout << m_className << "::LoadOpticalData:\n Read "
1853 << m_opticalDataEnergies.size() << " values from file "
1854 << filepath << "\n";
1855 }
1856 return true;
1857}
1858
1859bool MediumSilicon::ElectronScatteringRates() {
1860 // Reset the scattering rates
1861 m_cfTotElectronsX.assign(nEnergyStepsXL, 0.);
1862 m_cfTotElectronsL.assign(nEnergyStepsXL, 0.);
1863 m_cfTotElectronsG.assign(nEnergyStepsG, 0.);
1864 m_cfElectronsX.assign(nEnergyStepsXL, std::vector<double>());
1865 m_cfElectronsL.assign(nEnergyStepsXL, std::vector<double>());
1866 m_cfElectronsG.assign(nEnergyStepsG, std::vector<double>());
1867 m_energyLossElectronsX.clear();
1868 m_energyLossElectronsL.clear();
1869 m_energyLossElectronsG.clear();
1870 m_scatTypeElectronsX.clear();
1871 m_scatTypeElectronsL.clear();
1872 m_scatTypeElectronsG.clear();
1873 m_cfNullElectronsX = 0.;
1874 m_cfNullElectronsL = 0.;
1875 m_cfNullElectronsG = 0.;
1876
1877 m_nLevelsX = 0;
1878 m_nLevelsL = 0;
1879 m_nLevelsG = 0;
1880 // Fill the scattering rate tables.
1881 ElectronAcousticScatteringRates();
1882 ElectronOpticalScatteringRates();
1883 ElectronImpurityScatteringRates();
1884 ElectronIntervalleyScatteringRatesXX();
1885 ElectronIntervalleyScatteringRatesXL();
1886 ElectronIntervalleyScatteringRatesLL();
1887 ElectronIntervalleyScatteringRatesXGLG();
1888 ElectronIonisationRatesXL();
1889 ElectronIonisationRatesG();
1890
1891 if (m_debug) {
1892 std::cout << m_className << "::ElectronScatteringRates:\n"
1893 << " " << m_nLevelsX << " X-valley scattering terms\n"
1894 << " " << m_nLevelsL << " L-valley scattering terms\n"
1895 << " " << m_nLevelsG << " higher band scattering terms\n";
1896 }
1897
1898 std::ofstream outfileX;
1899 std::ofstream outfileL;
1900 if (m_cfOutput) {
1901 outfileX.open("ratesX.txt", std::ios::out);
1902 outfileL.open("ratesL.txt", std::ios::out);
1903 }
1904
1905 m_ieMinL = int(m_eMinL / m_eStepXL) + 1;
1906 for (int i = 0; i < nEnergyStepsXL; ++i) {
1907 // Sum up the scattering rates of all processes.
1908 for (int j = m_nLevelsX; j--;) m_cfTotElectronsX[i] += m_cfElectronsX[i][j];
1909 for (int j = m_nLevelsL; j--;) m_cfTotElectronsL[i] += m_cfElectronsL[i][j];
1910
1911 if (m_cfOutput) {
1912 outfileX << i * m_eStepXL << " " << m_cfTotElectronsX[i] << " ";
1913 for (int j = 0; j < m_nLevelsX; ++j) {
1914 outfileX << m_cfElectronsX[i][j] << " ";
1915 }
1916 outfileX << "\n";
1917 outfileL << i * m_eStepXL << " " << m_cfTotElectronsL[i] << " ";
1918 for (int j = 0; j < m_nLevelsL; ++j) {
1919 outfileL << m_cfElectronsL[i][j] << " ";
1920 }
1921 outfileL << "\n";
1922 }
1923
1924 if (m_cfTotElectronsX[i] > m_cfNullElectronsX) {
1925 m_cfNullElectronsX = m_cfTotElectronsX[i];
1926 }
1927 if (m_cfTotElectronsL[i] > m_cfNullElectronsL) {
1928 m_cfNullElectronsL = m_cfTotElectronsL[i];
1929 }
1930
1931 // Make sure the total scattering rate is positive.
1932 if (m_cfTotElectronsX[i] <= 0.) {
1933 std::cerr << m_className << "::ElectronScatteringRates:\n X-valley "
1934 << "scattering rate at " << i * m_eStepXL << " eV <= 0.\n";
1935 return false;
1936 }
1937 // Normalise the rates.
1938 for (int j = 0; j < m_nLevelsX; ++j) {
1939 m_cfElectronsX[i][j] /= m_cfTotElectronsX[i];
1940 if (j > 0) m_cfElectronsX[i][j] += m_cfElectronsX[i][j - 1];
1941 }
1942
1943 // Make sure the total scattering rate is positive.
1944 if (m_cfTotElectronsL[i] <= 0.) {
1945 if (i < m_ieMinL) continue;
1946 std::cerr << m_className << "::ElectronScatteringRates:\n L-valley "
1947 << "scattering rate at " << i * m_eStepXL << " eV <= 0.\n";
1948 return false;
1949 }
1950 // Normalise the rates.
1951 for (int j = 0; j < m_nLevelsL; ++j) {
1952 m_cfElectronsL[i][j] /= m_cfTotElectronsL[i];
1953 if (j > 0) m_cfElectronsL[i][j] += m_cfElectronsL[i][j - 1];
1954 }
1955 }
1956
1957 if (m_cfOutput) {
1958 outfileX.close();
1959 outfileL.close();
1960 }
1961
1962 std::ofstream outfileG;
1963 if (m_cfOutput) {
1964 outfileG.open("ratesG.txt", std::ios::out);
1965 }
1966 m_ieMinG = int(m_eMinG / m_eStepG) + 1;
1967 for (int i = 0; i < nEnergyStepsG; ++i) {
1968 // Sum up the scattering rates of all processes.
1969 for (int j = m_nLevelsG; j--;) m_cfTotElectronsG[i] += m_cfElectronsG[i][j];
1970
1971 if (m_cfOutput) {
1972 outfileG << i * m_eStepG << " " << m_cfTotElectronsG[i] << " ";
1973 for (int j = 0; j < m_nLevelsG; ++j) {
1974 outfileG << m_cfElectronsG[i][j] << " ";
1975 }
1976 outfileG << "\n";
1977 }
1978
1979 if (m_cfTotElectronsG[i] > m_cfNullElectronsG) {
1980 m_cfNullElectronsG = m_cfTotElectronsG[i];
1981 }
1982
1983 // Make sure the total scattering rate is positive.
1984 if (m_cfTotElectronsG[i] <= 0.) {
1985 if (i < m_ieMinG) continue;
1986 std::cerr << m_className << "::ElectronScatteringRates:\n Higher "
1987 << "band scattering rate at " << i * m_eStepG << " eV <= 0.\n";
1988 }
1989 // Normalise the rates.
1990 for (int j = 0; j < m_nLevelsG; ++j) {
1991 m_cfElectronsG[i][j] /= m_cfTotElectronsG[i];
1992 if (j > 0) m_cfElectronsG[i][j] += m_cfElectronsG[i][j - 1];
1993 }
1994 }
1995
1996 if (m_cfOutput) {
1997 outfileG.close();
1998 }
1999
2000 return true;
2001}
2002
2003bool MediumSilicon::ElectronAcousticScatteringRates() {
2004 // Reference:
2005 // - C. Jacoboni and L. Reggiani,
2006 // Rev. Mod. Phys. 55, 645-705
2007
2008 // Mass density [(eV/c2)/cm3]
2009 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2010 // Lattice temperature [eV]
2011 const double kbt = BoltzmannConstant * m_temperature;
2012
2013 // Acoustic phonon intraband scattering
2014 // Acoustic deformation potential [eV]
2015 constexpr double defpot2 = 9. * 9.;
2016 // Longitudinal velocity of sound [cm/ns]
2017 constexpr double u = 9.04e-4;
2018 // Prefactor for acoustic deformation potential scattering
2019 const double cIntra = TwoPi * SpeedOfLight * SpeedOfLight * kbt * defpot2 /
2020 (Hbar * u * u * rho);
2021
2022 // Fill the scattering rate tables.
2023 double en = Small;
2024 for (int i = 0; i < nEnergyStepsXL; ++i) {
2025 const double dosX = GetConductionBandDensityOfStates(en, 0);
2026 const double dosL = GetConductionBandDensityOfStates(en, m_nValleysX);
2027
2028 m_cfElectronsX[i].push_back(cIntra * dosX);
2029 m_cfElectronsL[i].push_back(cIntra * dosL);
2030 en += m_eStepXL;
2031 }
2032
2033 en = Small;
2034 for (int i = 0; i < nEnergyStepsG; ++i) {
2035 const double dosG =
2036 GetConductionBandDensityOfStates(en, m_nValleysX + m_nValleysL);
2037 m_cfElectronsG[i].push_back(cIntra * dosG);
2038 en += m_eStepG;
2039 }
2040
2041 // Assume that energy loss is negligible.
2042 m_energyLossElectronsX.push_back(0.);
2043 m_energyLossElectronsL.push_back(0.);
2044 m_energyLossElectronsG.push_back(0.);
2045 m_scatTypeElectronsX.push_back(ElectronCollisionTypeAcousticPhonon);
2046 m_scatTypeElectronsL.push_back(ElectronCollisionTypeAcousticPhonon);
2047 m_scatTypeElectronsG.push_back(ElectronCollisionTypeAcousticPhonon);
2048 ++m_nLevelsX;
2049 ++m_nLevelsL;
2050 ++m_nLevelsG;
2051
2052 return true;
2053}
2054
2055bool MediumSilicon::ElectronOpticalScatteringRates() {
2056 // Reference:
2057 // - K. Hess (editor),
2058 // Monte Carlo device simulation: full band and beyond
2059 // Chapter 5
2060 // - C. Jacoboni and L. Reggiani,
2061 // Rev. Mod. Phys. 55, 645-705
2062 // - M. Lundstrom,
2063 // Fundamentals of carrier transport
2064
2065 // Mass density [(eV/c2)/cm3]
2066 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2067 // Lattice temperature [eV]
2068 const double kbt = BoltzmannConstant * m_temperature;
2069
2070 // Coupling constant [eV/cm]
2071 constexpr double dtk = 2.2e8;
2072 // Phonon energy [eV]
2073 constexpr double eph = 63.0e-3;
2074 // Phonon cccupation numbers
2075 const double nocc = 1. / (exp(eph / kbt) - 1);
2076 // Prefactors
2077 const double c0 = HbarC * SpeedOfLight * Pi / rho;
2078 double c = c0 * dtk * dtk / eph;
2079
2080 double en = 0.;
2081 // L valleys
2082 /*
2083 for (int i = 0; i < nEnergyStepsXL; ++i) {
2084 // Absorption
2085 if (en > m_eMinL) {
2086 double dos = GetConductionBandDensityOfStates(en + eph, m_nValleysX);
2087 m_cfElectronsL[i].push_back(c * nocc * dos);
2088 } else {
2089 m_cfElectronsL[i].push_back(0.);
2090 }
2091 // Emission
2092 if (en > m_eMinL + eph) {
2093 double dos = GetConductionBandDensityOfStates(en - eph, m_nValleysX);
2094 m_cfElectronsL[i].push_back(c * (nocc + 1) * dos);
2095 } else {
2096 m_cfElectronsL[i].push_back(0.);
2097 }
2098 en += m_eStepXL;
2099 }
2100 //*/
2101
2102 en = 0.;
2103 // Higher band(s)
2104 for (int i = 0; i < nEnergyStepsG; ++i) {
2105 // Absorption
2106 if (en > m_eMinG) {
2107 double dos =
2108 GetConductionBandDensityOfStates(en + eph, m_nValleysX + m_nValleysL);
2109 m_cfElectronsG[i].push_back(c * nocc * dos);
2110 } else {
2111 m_cfElectronsG[i].push_back(0.);
2112 }
2113 // Emission
2114 if (en > m_eMinG + eph) {
2115 double dos =
2116 GetConductionBandDensityOfStates(en - eph, m_nValleysX + m_nValleysL);
2117 m_cfElectronsG[i].push_back(c * (nocc + 1) * dos);
2118 } else {
2119 m_cfElectronsG[i].push_back(0.);
2120 }
2121 en += m_eStepG;
2122 }
2123
2124 // Absorption
2125 // m_energyLossElectronsL.push_back(-eph);
2126 m_energyLossElectronsG.push_back(-eph);
2127 // Emission
2128 // m_energyLossElectronsL.push_back(eph);
2129 m_energyLossElectronsG.push_back(eph);
2130 // m_scatTypeElectronsL.push_back(ElectronCollisionTypeOpticalPhonon);
2131 // m_scatTypeElectronsL.push_back(ElectronCollisionTypeOpticalPhonon);
2132 m_scatTypeElectronsG.push_back(ElectronCollisionTypeOpticalPhonon);
2133 m_scatTypeElectronsG.push_back(ElectronCollisionTypeOpticalPhonon);
2134
2135 // m_nLevelsL += 2;
2136 m_nLevelsG += 2;
2137
2138 return true;
2139}
2140
2141bool MediumSilicon::ElectronIntervalleyScatteringRatesXX() {
2142 // Reference:
2143 // - C. Jacoboni and L. Reggiani,
2144 // Rev. Mod. Phys. 55, 645-705
2145
2146 // Mass density [(eV/c2)/cm3]
2147 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2148 // Lattice temperature [eV]
2149 const double kbt = BoltzmannConstant * m_temperature;
2150
2151 constexpr unsigned int nPhonons = 6;
2152 // f-type scattering: transition between orthogonal axes (multiplicity 4)
2153 // g-type scattering: transition between opposite axes (multiplicity 1)
2154 // Sequence of transitions in the table:
2155 // TA (g) - LA (g) - LO (g) - TA (f) - LA (f) - TO (f)
2156 // Coupling constants [eV/cm]
2157 constexpr double dtk[nPhonons] = {0.5e8, 0.8e8, 1.1e9, 0.3e8, 2.0e8, 2.0e8};
2158 // Phonon energies [eV]
2159 constexpr double eph[nPhonons] = {12.06e-3, 18.53e-3, 62.04e-3,
2160 18.86e-3, 47.39e-3, 59.03e-3};
2161 // Phonon cccupation numbers
2162 double nocc[nPhonons];
2163 // Prefactors
2164 const double c0 = HbarC * SpeedOfLight * Pi / rho;
2165 double c[nPhonons];
2166
2167 for (unsigned int j = 0; j < nPhonons; ++j) {
2168 nocc[j] = 1. / (exp(eph[j] / kbt) - 1);
2169 c[j] = c0 * dtk[j] * dtk[j] / eph[j];
2170 if (j > 2) c[j] *= 4;
2171 }
2172
2173 double en = 0.;
2174 for (int i = 0; i < nEnergyStepsXL; ++i) {
2175 for (unsigned int j = 0; j < nPhonons; ++j) {
2176 // Absorption
2177 double dos = GetConductionBandDensityOfStates(en + eph[j], 0);
2178 m_cfElectronsX[i].push_back(c[j] * nocc[j] * dos);
2179 // Emission
2180 if (en > eph[j]) {
2181 dos = GetConductionBandDensityOfStates(en - eph[j], 0);
2182 m_cfElectronsX[i].push_back(c[j] * (nocc[j] + 1) * dos);
2183 } else {
2184 m_cfElectronsX[i].push_back(0.);
2185 }
2186 }
2187 en += m_eStepXL;
2188 }
2189
2190 for (unsigned int j = 0; j < nPhonons; ++j) {
2191 // Absorption
2192 m_energyLossElectronsX.push_back(-eph[j]);
2193 // Emission
2194 m_energyLossElectronsX.push_back(eph[j]);
2195 if (j <= 2) {
2196 m_scatTypeElectronsX.push_back(ElectronCollisionTypeIntervalleyG);
2197 m_scatTypeElectronsX.push_back(ElectronCollisionTypeIntervalleyG);
2198 } else {
2199 m_scatTypeElectronsX.push_back(ElectronCollisionTypeIntervalleyF);
2200 m_scatTypeElectronsX.push_back(ElectronCollisionTypeIntervalleyF);
2201 }
2202 }
2203
2204 m_nLevelsX += 2 * nPhonons;
2205
2206 return true;
2207}
2208
2209bool MediumSilicon::ElectronIntervalleyScatteringRatesXL() {
2210 // Reference:
2211 // - M. Lundstrom, Fundamentals of carrier transport
2212 // - M. Martin et al.,
2213 // Semicond. Sci. Technol. 8, 1291-1297
2214
2215 // Mass density [(eV/c2)/cm3]
2216 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2217 // Lattice temperature [eV]
2218 const double kbt = BoltzmannConstant * m_temperature;
2219
2220 constexpr unsigned int nPhonons = 4;
2221
2222 // Coupling constants [eV/cm]
2223 constexpr double dtk[nPhonons] = {2.e8, 2.e8, 2.e8, 2.e8};
2224 // Phonon energies [eV]
2225 constexpr double eph[nPhonons] = {58.e-3, 55.e-3, 41.e-3, 17.e-3};
2226 // Number of equivalent valleys
2227 constexpr unsigned int zX = 6;
2228 constexpr unsigned int zL = 8;
2229
2230 // Phonon cccupation numbers
2231 double nocc[nPhonons] = {0.};
2232 // Prefactors
2233 const double c0 = HbarC * SpeedOfLight * Pi / rho;
2234 double c[nPhonons];
2235
2236 for (unsigned int j = 0; j < nPhonons; ++j) {
2237 nocc[j] = 1. / (exp(eph[j] / kbt) - 1);
2238 c[j] = c0 * dtk[j] * dtk[j] / eph[j];
2239 }
2240
2241 double en = 0.;
2242 for (int i = 0; i < nEnergyStepsXL; ++i) {
2243 for (unsigned int j = 0; j < nPhonons; ++j) {
2244 // XL
2245 // Absorption
2246 if (en + eph[j] > m_eMinL) {
2247 double dos = GetConductionBandDensityOfStates(en + eph[j], m_nValleysX);
2248 m_cfElectronsX[i].push_back(zL * c[j] * nocc[j] * dos);
2249 } else {
2250 m_cfElectronsX[i].push_back(0.);
2251 }
2252 // Emission
2253 if (en - eph[j] > m_eMinL) {
2254 double dos = GetConductionBandDensityOfStates(en - eph[j], m_nValleysX);
2255 m_cfElectronsX[i].push_back(zL * c[j] * (nocc[j] + 1) * dos);
2256 } else {
2257 m_cfElectronsX[i].push_back(0.);
2258 }
2259 // LX
2260 if (en > m_eMinL) {
2261 // Absorption
2262 double dos = GetConductionBandDensityOfStates(en + eph[j], 0);
2263 m_cfElectronsL[i].push_back(zX * c[j] * nocc[j] * dos);
2264 // Emission
2265 dos = GetConductionBandDensityOfStates(en - eph[j], 0);
2266 m_cfElectronsL[i].push_back(zX * c[j] * (nocc[j] + 1) * dos);
2267 } else {
2268 m_cfElectronsL[i].push_back(0.);
2269 m_cfElectronsL[i].push_back(0.);
2270 }
2271 }
2272 en += m_eStepXL;
2273 }
2274
2275 for (unsigned int j = 0; j < nPhonons; ++j) {
2276 // Absorption
2277 m_energyLossElectronsX.push_back(-eph[j]);
2278 m_energyLossElectronsL.push_back(-eph[j]);
2279 // Emission
2280 m_energyLossElectronsX.push_back(eph[j]);
2281 m_energyLossElectronsL.push_back(eph[j]);
2282 m_scatTypeElectronsX.push_back(ElectronCollisionTypeInterbandXL);
2283 m_scatTypeElectronsX.push_back(ElectronCollisionTypeInterbandXL);
2284 m_scatTypeElectronsL.push_back(ElectronCollisionTypeInterbandXL);
2285 m_scatTypeElectronsL.push_back(ElectronCollisionTypeInterbandXL);
2286 }
2287
2288 m_nLevelsX += 2 * nPhonons;
2289 m_nLevelsL += 2 * nPhonons;
2290
2291 return true;
2292}
2293
2294bool MediumSilicon::ElectronIntervalleyScatteringRatesLL() {
2295 // Reference:
2296 // - K. Hess (editor),
2297 // Monte Carlo device simulation: full band and beyond
2298 // Chapter 5
2299 // - M. J. Martin et al.,
2300 // Semicond. Sci. Technol. 8, 1291-1297
2301
2302 // Mass density [(eV/c2)/cm3]
2303 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2304 // Lattice temperature [eV]
2305 const double kbt = BoltzmannConstant * m_temperature;
2306
2307 const int nPhonons = 1;
2308 // Coupling constant [eV/cm]
2309 const double dtk[nPhonons] = {2.63e8};
2310 // Phonon energy [eV]
2311 const double eph[nPhonons] = {38.87e-3};
2312 // Phonon cccupation numbers
2313 double nocc[nPhonons];
2314 // Prefactors
2315 const double c0 = HbarC * SpeedOfLight * Pi / rho;
2316 double c[nPhonons];
2317
2318 for (int j = 0; j < nPhonons; ++j) {
2319 nocc[j] = 1. / (exp(eph[j] / kbt) - 1);
2320 c[j] = c0 * dtk[j] * dtk[j] / eph[j];
2321 c[j] *= 7;
2322 }
2323
2324 double en = 0.;
2325 double dos = 0.;
2326 for (int i = 0; i < nEnergyStepsXL; ++i) {
2327 for (int j = 0; j < nPhonons; ++j) {
2328 // Absorption
2329 dos = GetConductionBandDensityOfStates(en + eph[j], m_nValleysX);
2330 m_cfElectronsL[i].push_back(c[j] * nocc[j] * dos);
2331 // Emission
2332 if (en > m_eMinL + eph[j]) {
2333 dos = GetConductionBandDensityOfStates(en - eph[j], m_nValleysX);
2334 m_cfElectronsL[i].push_back(c[j] * (nocc[j] + 1) * dos);
2335 } else {
2336 m_cfElectronsL[i].push_back(0.);
2337 }
2338 }
2339 en += m_eStepXL;
2340 }
2341
2342 for (int j = 0; j < nPhonons; ++j) {
2343 // Absorption
2344 m_energyLossElectronsL.push_back(-eph[j]);
2345 // Emission
2346 m_energyLossElectronsL.push_back(eph[j]);
2347 m_scatTypeElectronsL.push_back(ElectronCollisionTypeIntervalleyF);
2348 m_scatTypeElectronsL.push_back(ElectronCollisionTypeIntervalleyF);
2349 }
2350
2351 m_nLevelsL += 2 * nPhonons;
2352
2353 return true;
2354}
2355
2356bool MediumSilicon::ElectronIntervalleyScatteringRatesXGLG() {
2357 // Reference:
2358 // - K. Hess (editor),
2359 // Monte Carlo device simulation: full band and beyond
2360 // Chapter 5
2361
2362 // Mass density [(eV/c2)/cm3]
2363 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2364 // Lattice temperature [eV]
2365 const double kbt = BoltzmannConstant * m_temperature;
2366
2367 const int nPhonons = 1;
2368
2369 // Coupling constants [eV/cm]
2370 // Average of XG and LG
2371 const double dtk[nPhonons] = {2.43e8};
2372 // Phonon energies [eV]
2373 const double eph[nPhonons] = {37.65e-3};
2374 // Number of equivalent valleys
2375 const int zX = 6;
2376 const int zL = 8;
2377 const int zG = 1;
2378
2379 // Phonon cccupation numbers
2380 double nocc[nPhonons] = {0.};
2381 // Prefactors
2382 const double c0 = HbarC * SpeedOfLight * Pi / rho;
2383 double c[nPhonons];
2384
2385 for (int j = 0; j < nPhonons; ++j) {
2386 nocc[j] = 1. / (exp(eph[j] / kbt) - 1);
2387 c[j] = c0 * dtk[j] * dtk[j] / eph[j];
2388 }
2389
2390 double en = 0.;
2391 double dos = 0.;
2392 // XG, LG
2393 for (int i = 0; i < nEnergyStepsXL; ++i) {
2394 for (int j = 0; j < nPhonons; ++j) {
2395 // Absorption
2396 if (en + eph[j] > m_eMinG) {
2397 dos = GetConductionBandDensityOfStates(en + eph[j],
2398 m_nValleysX + m_nValleysL);
2399 m_cfElectronsX[i].push_back(zG * c[j] * nocc[j] * dos);
2400 m_cfElectronsL[i].push_back(zG * c[j] * nocc[j] * dos);
2401 } else {
2402 m_cfElectronsX[i].push_back(0.);
2403 m_cfElectronsL[i].push_back(0.);
2404 }
2405 // Emission
2406 if (en - eph[j] > m_eMinG) {
2407 dos = GetConductionBandDensityOfStates(en - eph[j],
2408 m_nValleysX + m_nValleysL);
2409 m_cfElectronsX[i].push_back(zG * c[j] * (nocc[j] + 1) * dos);
2410 m_cfElectronsL[i].push_back(zG * c[j] * (nocc[j] + 1) * dos);
2411 } else {
2412 m_cfElectronsX[i].push_back(0.);
2413 m_cfElectronsL[i].push_back(0.);
2414 }
2415 }
2416 en += m_eStepXL;
2417 }
2418
2419 // GX, GL
2420 en = 0.;
2421 double dosX = 0., dosL = 0.;
2422 for (int i = 0; i < nEnergyStepsG; ++i) {
2423 for (int j = 0; j < nPhonons; ++j) {
2424 // Absorption
2425 dosX = GetConductionBandDensityOfStates(en + eph[j], 0);
2426 dosL = GetConductionBandDensityOfStates(en + eph[j], m_nValleysX);
2427 m_cfElectronsG[i].push_back(zX * c[j] * nocc[j] * dosX);
2428 if (en > m_eMinL) {
2429 m_cfElectronsG[i].push_back(zL * c[j] * nocc[j] * dosL);
2430 } else {
2431 m_cfElectronsG[i].push_back(0.);
2432 }
2433 // Emission
2434 dosX = GetConductionBandDensityOfStates(en - eph[j], 0);
2435 dosL = GetConductionBandDensityOfStates(en - eph[j], m_nValleysX);
2436 if (en > eph[j]) {
2437 m_cfElectronsG[i].push_back(zX * c[j] * (nocc[j] + 1) * dosX);
2438 } else {
2439 m_cfElectronsG[i].push_back(0.);
2440 }
2441 if (en - eph[j] > m_eMinL) {
2442 m_cfElectronsG[i].push_back(zL * c[j] * (nocc[j] + 1) * dosL);
2443 } else {
2444 m_cfElectronsG[i].push_back(0.);
2445 }
2446 }
2447 en += m_eStepG;
2448 }
2449
2450 for (int j = 0; j < nPhonons; ++j) {
2451 // Absorption (XL)
2452 m_energyLossElectronsX.push_back(-eph[j]);
2453 m_energyLossElectronsL.push_back(-eph[j]);
2454 // Emission (XL)
2455 m_energyLossElectronsX.push_back(eph[j]);
2456 m_energyLossElectronsL.push_back(eph[j]);
2457 // Absorption (G)
2458 m_energyLossElectronsG.push_back(-eph[j]);
2459 m_energyLossElectronsG.push_back(-eph[j]);
2460 // Emission (G)
2461 m_energyLossElectronsG.push_back(eph[j]);
2462 m_energyLossElectronsG.push_back(eph[j]);
2463
2464 m_scatTypeElectronsX.push_back(ElectronCollisionTypeInterbandXG);
2465 m_scatTypeElectronsX.push_back(ElectronCollisionTypeInterbandXG);
2466 m_scatTypeElectronsL.push_back(ElectronCollisionTypeInterbandLG);
2467 m_scatTypeElectronsL.push_back(ElectronCollisionTypeInterbandLG);
2468
2469 m_scatTypeElectronsG.push_back(ElectronCollisionTypeInterbandXG);
2470 m_scatTypeElectronsG.push_back(ElectronCollisionTypeInterbandLG);
2471 m_scatTypeElectronsG.push_back(ElectronCollisionTypeInterbandXG);
2472 m_scatTypeElectronsG.push_back(ElectronCollisionTypeInterbandLG);
2473 }
2474
2475 m_nLevelsX += 2 * nPhonons;
2476 m_nLevelsL += 2 * nPhonons;
2477 m_nLevelsG += 4 * nPhonons;
2478
2479 return true;
2480}
2481
2482bool MediumSilicon::ElectronIonisationRatesXL() {
2483 // References:
2484 // - E. Cartier, M. V. Fischetti, E. A. Eklund and F. R. McFeely,
2485 // Appl. Phys. Lett 62, 3339-3341
2486 // - DAMOCLES web page: www.research.ibm.com/DAMOCLES
2487
2488 // Coefficients [ns-1]
2489 constexpr double p[3] = {6.25e1, 3.e3, 6.8e5};
2490 // Threshold energies [eV]
2491 constexpr double eth[3] = {1.2, 1.8, 3.45};
2492
2493 double en = 0.;
2494 for (int i = 0; i < nEnergyStepsXL; ++i) {
2495 double fIon = 0.;
2496 if (en > eth[0]) {
2497 fIon += p[0] * (en - eth[0]) * (en - eth[0]);
2498 }
2499 if (en > eth[1]) {
2500 fIon += p[1] * (en - eth[1]) * (en - eth[1]);
2501 }
2502 if (en > eth[2]) {
2503 fIon += p[2] * (en - eth[2]) * (en - eth[2]);
2504 }
2505 m_cfElectronsX[i].push_back(fIon);
2506 m_cfElectronsL[i].push_back(fIon);
2507 en += m_eStepXL;
2508 }
2509
2510 m_energyLossElectronsX.push_back(eth[0]);
2511 m_energyLossElectronsL.push_back(eth[0]);
2512 m_scatTypeElectronsX.push_back(ElectronCollisionTypeIonisation);
2513 m_scatTypeElectronsL.push_back(ElectronCollisionTypeIonisation);
2514 ++m_nLevelsX;
2515 ++m_nLevelsL;
2516
2517 return true;
2518}
2519
2520bool MediumSilicon::ElectronIonisationRatesG() {
2521 // References:
2522 // - E. Cartier, M. V. Fischetti, E. A. Eklund and F. R. McFeely,
2523 // Appl. Phys. Lett 62, 3339-3341
2524 // - S. Tanuma, C. J. Powell and D. R. Penn
2525 // Surf. Interface Anal. (2010)
2526
2527 // Coefficients [ns-1]
2528 constexpr double p[3] = {6.25e1, 3.e3, 6.8e5};
2529 // Threshold energies [eV]
2530 constexpr double eth[3] = {1.2, 1.8, 3.45};
2531
2532 double en = 0.;
2533 for (int i = 0; i < nEnergyStepsG; ++i) {
2534 double fIon = 0.;
2535 if (en > eth[0]) {
2536 fIon += p[0] * (en - eth[0]) * (en - eth[0]);
2537 }
2538 if (en > eth[1]) {
2539 fIon += p[1] * (en - eth[1]) * (en - eth[1]);
2540 }
2541 if (en > eth[2]) {
2542 fIon += p[2] * (en - eth[2]) * (en - eth[2]);
2543 }
2544 if (en >= m_eMinG) {
2545 m_cfElectronsG[i].push_back(fIon);
2546 } else {
2547 m_cfElectronsG[i].push_back(0.);
2548 }
2549 en += m_eStepG;
2550 }
2551
2552 m_energyLossElectronsG.push_back(eth[0]);
2553 m_scatTypeElectronsG.push_back(ElectronCollisionTypeIonisation);
2554 ++m_nLevelsG;
2555
2556 return true;
2557}
2558
2559bool MediumSilicon::ElectronImpurityScatteringRates() {
2560 // Lattice temperature [eV]
2561 const double kbt = BoltzmannConstant * m_temperature;
2562
2563 // Band parameters
2564 // Density of states effective masses
2565 const double mdX =
2566 ElectronMass * pow(m_mLongX * m_mTransX * m_mTransX, 1. / 3.);
2567 const double mdL =
2568 ElectronMass * pow(m_mLongL * m_mTransL * m_mTransL, 1. / 3.);
2569
2570 // Dielectric constant
2571 const double eps = GetDielectricConstant();
2572 // Impurity concentration
2573 const double impurityConcentration = m_dopingConcentration;
2574 if (impurityConcentration < Small) return true;
2575
2576 // Screening length
2577 const double ls = sqrt(eps * kbt / (4 * Pi * FineStructureConstant * HbarC *
2578 impurityConcentration));
2579 const double ebX = 0.5 * HbarC * HbarC / (mdX * ls * ls);
2580 const double ebL = 0.5 * HbarC * HbarC / (mdL * ls * ls);
2581
2582 // Prefactor
2583 // const double c = pow(2., 2.5) * Pi * impurityConcentration *
2584 // pow(FineStructureConstant * HbarC, 2) *
2585 // SpeedOfLight / (eps * eps * sqrt(md) * eb * eb);
2586 // Use momentum-transfer cross-section
2587 const double cX = impurityConcentration * Pi *
2588 pow(FineStructureConstant * HbarC, 2) * SpeedOfLight /
2589 (sqrt(2 * mdX) * eps * eps);
2590 const double cL = impurityConcentration * Pi *
2591 pow(FineStructureConstant * HbarC, 2) * SpeedOfLight /
2592 (sqrt(2 * mdL) * eps * eps);
2593
2594 double en = 0.;
2595 for (int i = 0; i < nEnergyStepsXL; ++i) {
2596 const double gammaX = en * (1. + m_alphaX * en);
2597 const double gammaL = (en - m_eMinL) * (1. + m_alphaL * (en - m_eMinL));
2598 // m_cfElectrons[i][iLevel] = c * sqrt(gamma) * (1. + 2 * alpha * en) /
2599 // (1. + 4. * gamma / eb);
2600 if (gammaX <= 0.) {
2601 m_cfElectronsX[i].push_back(0.);
2602 } else {
2603 const double b = 4 * gammaX / ebX;
2604 m_cfElectronsX[i].push_back((cX / pow(gammaX, 1.5)) *
2605 (log(1. + b) - b / (1. + b)));
2606 }
2607 if (en <= m_eMinL || gammaL <= 0.) {
2608 m_cfElectronsL[i].push_back(0.);
2609 } else {
2610 const double b = 4 * gammaL / ebL;
2611 m_cfElectronsL[i].push_back((cL / pow(gammaL, 1.5)) *
2612 (log(1. + b) - b / (1. + b)));
2613 }
2614 en += m_eStepXL;
2615 }
2616
2617 m_energyLossElectronsX.push_back(0.);
2618 m_energyLossElectronsL.push_back(0.);
2619 m_scatTypeElectronsX.push_back(ElectronCollisionTypeImpurity);
2620 m_scatTypeElectronsL.push_back(ElectronCollisionTypeImpurity);
2621 ++m_nLevelsX;
2622 ++m_nLevelsL;
2623
2624 return true;
2625}
2626
2627bool MediumSilicon::HoleScatteringRates() {
2628 // Reset the scattering rates
2629 m_cfTotHoles.assign(nEnergyStepsV, 0.);
2630 m_cfHoles.assign(nEnergyStepsV, std::vector<double>());
2631 m_energyLossHoles.clear();
2632 m_scatTypeHoles.clear();
2633 m_cfNullHoles = 0.;
2634
2635 m_nLevelsV = 0;
2636 // Fill the scattering rates table
2637 HoleAcousticScatteringRates();
2638 HoleOpticalScatteringRates();
2639 // HoleImpurityScatteringRates();
2640 HoleIonisationRates();
2641
2642 std::ofstream outfile;
2643 if (m_cfOutput) {
2644 outfile.open("ratesV.txt", std::ios::out);
2645 }
2646
2647 for (int i = 0; i < nEnergyStepsV; ++i) {
2648 // Sum up the scattering rates of all processes.
2649 for (int j = m_nLevelsV; j--;) m_cfTotHoles[i] += m_cfHoles[i][j];
2650
2651 if (m_cfOutput) {
2652 outfile << i * m_eStepV << " " << m_cfTotHoles[i] << " ";
2653 for (int j = 0; j < m_nLevelsV; ++j) {
2654 outfile << m_cfHoles[i][j] << " ";
2655 }
2656 outfile << "\n";
2657 }
2658
2659 if (m_cfTotHoles[i] > m_cfNullHoles) {
2660 m_cfNullHoles = m_cfTotHoles[i];
2661 }
2662
2663 // Make sure the total scattering rate is positive.
2664 if (m_cfTotHoles[i] <= 0.) {
2665 std::cerr << m_className << "::HoleScatteringRates:\n"
2666 << " Scattering rate at " << i * m_eStepV << " eV <= 0.\n";
2667 return false;
2668 }
2669 // Normalise the rates.
2670 for (int j = 0; j < m_nLevelsV; ++j) {
2671 m_cfHoles[i][j] /= m_cfTotHoles[i];
2672 if (j > 0) m_cfHoles[i][j] += m_cfHoles[i][j - 1];
2673 }
2674 }
2675
2676 if (m_cfOutput) {
2677 outfile.close();
2678 }
2679
2680 return true;
2681}
2682
2683bool MediumSilicon::HoleAcousticScatteringRates() {
2684 // Reference:
2685 // - C. Jacoboni and L. Reggiani,
2686 // Rev. Mod. Phys. 55, 645-705
2687 // - DAMOCLES web page: www.research.ibm.com/DAMOCLES
2688 // - M. Lundstrom, Fundamentals of carrier transport
2689
2690 // Mass density [(eV/c2)/cm3]
2691 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2692 // Lattice temperature [eV]
2693 const double kbt = BoltzmannConstant * m_temperature;
2694
2695 // Acoustic phonon intraband scattering
2696 // Acoustic deformation potential [eV]
2697 // DAMOCLES: 4.6 eV; Lundstrom: 5 eV
2698 constexpr double defpot2 = 4.6 * 4.6;
2699 // Longitudinal velocity of sound [cm/ns]
2700 constexpr double u = 9.04e-4;
2701 // Prefactor for acoustic deformation potential scattering
2702 const double cIntra = TwoPi * SpeedOfLight * SpeedOfLight * kbt * defpot2 /
2703 (Hbar * u * u * rho);
2704
2705 // Fill the scattering rate tables.
2706 double en = Small;
2707 for (int i = 0; i < nEnergyStepsV; ++i) {
2708 const double dos = GetValenceBandDensityOfStates(en, 0);
2709 m_cfHoles[i].push_back(cIntra * dos);
2710 en += m_eStepV;
2711 }
2712
2713 // Assume that energy loss is negligible.
2714 m_energyLossHoles.push_back(0.);
2715 m_scatTypeHoles.push_back(ElectronCollisionTypeAcousticPhonon);
2716 ++m_nLevelsV;
2717
2718 return true;
2719}
2720
2721bool MediumSilicon::HoleOpticalScatteringRates() {
2722 // Reference:
2723 // - C. Jacoboni and L. Reggiani,
2724 // Rev. Mod. Phys. 55, 645-705
2725 // - DAMOCLES web page: www.research.ibm.com/DAMOCLES
2726 // - M. Lundstrom, Fundamentals of carrier transport
2727
2728 // Mass density [(eV/c2)/cm3]
2729 const double rho = m_density * m_a * AtomicMassUnitElectronVolt;
2730 // Lattice temperature [eV]
2731 const double kbt = BoltzmannConstant * m_temperature;
2732
2733 // Coupling constant [eV/cm]
2734 // DAMOCLES: 6.6, Lundstrom: 6.0
2735 constexpr double dtk = 6.6e8;
2736 // Phonon energy [eV]
2737 constexpr double eph = 63.0e-3;
2738 // Phonon cccupation numbers
2739 const double nocc = 1. / (exp(eph / kbt) - 1);
2740 // Prefactors
2741 const double c0 = HbarC * SpeedOfLight * Pi / rho;
2742 double c = c0 * dtk * dtk / eph;
2743
2744 double en = 0.;
2745 for (int i = 0; i < nEnergyStepsV; ++i) {
2746 // Absorption
2747 double dos = GetValenceBandDensityOfStates(en + eph, 0);
2748 m_cfHoles[i].push_back(c * nocc * dos);
2749 // Emission
2750 if (en > eph) {
2751 dos = GetValenceBandDensityOfStates(en - eph, 0);
2752 m_cfHoles[i].push_back(c * (nocc + 1) * dos);
2753 } else {
2754 m_cfHoles[i].push_back(0.);
2755 }
2756 en += m_eStepV;
2757 }
2758
2759 // Absorption
2760 m_energyLossHoles.push_back(-eph);
2761 // Emission
2762 m_energyLossHoles.push_back(eph);
2763 m_scatTypeHoles.push_back(ElectronCollisionTypeOpticalPhonon);
2764 m_scatTypeHoles.push_back(ElectronCollisionTypeOpticalPhonon);
2765
2766 m_nLevelsV += 2;
2767
2768 return true;
2769}
2770
2771bool MediumSilicon::HoleIonisationRates() {
2772 // References:
2773 // - DAMOCLES web page: www.research.ibm.com/DAMOCLES
2774
2775 // Coefficients [ns-1]
2776 constexpr double p[2] = {2., 1.e3};
2777 // Threshold energies [eV]
2778 constexpr double eth[2] = {1.1, 1.45};
2779 // Exponents
2780 constexpr double b[2] = {6., 4.};
2781
2782 double en = 0.;
2783 for (int i = 0; i < nEnergyStepsV; ++i) {
2784 double fIon = 0.;
2785 if (en > eth[0]) {
2786 fIon += p[0] * pow(en - eth[0], b[0]);
2787 }
2788 if (en > eth[1]) {
2789 fIon += p[1] * pow(en - eth[1], b[1]);
2790 }
2791 m_cfHoles[i].push_back(fIon);
2792 en += m_eStepV;
2793 }
2794
2795 m_energyLossHoles.push_back(eth[0]);
2796 m_scatTypeHoles.push_back(ElectronCollisionTypeIonisation);
2797 ++m_nLevelsV;
2798
2799 return true;
2800}
2801
2803 const int band) {
2804 if (band < 0) {
2805 int iE = int(e / m_eStepDos);
2806 const int nPoints = m_fbDosConduction.size();
2807 if (iE >= nPoints || iE < 0) {
2808 return 0.;
2809 } else if (iE == nPoints - 1) {
2810 return m_fbDosConduction[nPoints - 1];
2811 }
2812
2813 const double dos = m_fbDosConduction[iE] +
2814 (m_fbDosConduction[iE + 1] - m_fbDosConduction[iE]) *
2815 (e / m_eStepDos - iE);
2816 return dos * 1.e21;
2817
2818 } else if (band < m_nValleysX) {
2819 // X valleys
2820 if (e <= 0.) return 0.;
2821 // Density-of-states effective mass (cube)
2822 const double md3 = pow(ElectronMass, 3) * m_mLongX * m_mTransX * m_mTransX;
2823
2824 if (m_fullBandDos) {
2825 if (e < m_eMinL) {
2826 return GetConductionBandDensityOfStates(e, -1) / m_nValleysX;
2827 } else if (e < m_eMinG) {
2828 // Subtract the fraction of the full-band density of states
2829 // attributed to the L valleys.
2830 const double dosX =
2832 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL;
2833 return dosX / m_nValleysX;
2834 } else {
2835 // Subtract the fraction of the full-band density of states
2836 // attributed to the L valleys and the higher bands.
2837 const double dosX =
2839 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL -
2840 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
2841 if (dosX <= 0.) return 0.;
2842 return dosX / m_nValleysX;
2843 }
2844 } else if (m_nonParabolic) {
2845 const double alpha = m_alphaX;
2846 return sqrt(md3 * e * (1. + alpha * e) / 2.) * (1. + 2 * alpha * e) /
2847 (Pi2 * pow(HbarC, 3.));
2848 } else {
2849 return sqrt(md3 * e / 2.) / (Pi2 * pow(HbarC, 3.));
2850 }
2851 } else if (band < m_nValleysX + m_nValleysL) {
2852 // L valleys
2853 if (e <= m_eMinL) return 0.;
2854
2855 // Density-of-states effective mass (cube)
2856 const double md3 = pow(ElectronMass, 3) * m_mLongL * m_mTransL * m_mTransL;
2857 // Non-parabolicity parameter
2858 const double alpha = m_alphaL;
2859
2860 if (m_fullBandDos) {
2861 // Energy up to which the non-parabolic approximation is used.
2862 const double ej = m_eMinL + 0.5;
2863 if (e <= ej) {
2864 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
2865 (1. + 2 * alpha * (e - m_eMinL)) /
2866 (Sqrt2 * Pi2 * pow(HbarC, 3.));
2867 } else {
2868 // Fraction of full-band density of states attributed to L valleys
2869 double fL = sqrt(md3 * (ej - m_eMinL) * (1. + alpha * (ej - m_eMinL))) *
2870 (1. + 2 * alpha * (ej - m_eMinL)) /
2871 (Sqrt2 * Pi2 * pow(HbarC, 3.));
2872 fL = m_nValleysL * fL / GetConductionBandDensityOfStates(ej, -1);
2873
2874 double dosXL = GetConductionBandDensityOfStates(e, -1);
2875 if (e > m_eMinG) {
2876 dosXL -=
2877 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
2878 }
2879 if (dosXL <= 0.) return 0.;
2880 return fL * dosXL / 8.;
2881 }
2882 } else if (m_nonParabolic) {
2883 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
2884 (1. + 2 * alpha * (e - m_eMinL)) / (Sqrt2 * Pi2 * pow(HbarC, 3.));
2885 } else {
2886 return sqrt(md3 * (e - m_eMinL) / 2.) / (Pi2 * pow(HbarC, 3.));
2887 }
2888 } else if (band == m_nValleysX + m_nValleysL) {
2889 // Higher bands
2890 const double ej = 2.7;
2891 if (m_eMinG >= ej) {
2892 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2893 << " Cannot determine higher band density-of-states.\n"
2894 << " Program bug. Check offset energy!\n";
2895 return 0.;
2896 }
2897 if (e < m_eMinG) {
2898 return 0.;
2899 } else if (e < ej) {
2900 // Coexistence of XL and higher bands.
2901 const double dj = GetConductionBandDensityOfStates(ej, -1);
2902 // Assume linear increase of density-of-states.
2903 return dj * (e - m_eMinG) / (ej - m_eMinG);
2904 } else {
2906 }
2907 }
2908
2909 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2910 << " Band index (" << band << ") out of range.\n";
2911 return ElectronMass * sqrt(ElectronMass * e / 2.) / (Pi2 * pow(HbarC, 3.));
2912}
2913
2915 const int band) {
2916 if (band <= 0) {
2917 // Total (full-band) density of states.
2918 const int nPoints = m_fbDosValence.size();
2919 int iE = int(e / m_eStepDos);
2920 if (iE >= nPoints || iE < 0) {
2921 return 0.;
2922 } else if (iE == nPoints - 1) {
2923 return m_fbDosValence[nPoints - 1];
2924 }
2925
2926 const double dos =
2927 m_fbDosValence[iE] +
2928 (m_fbDosValence[iE + 1] - m_fbDosValence[iE]) * (e / m_eStepDos - iE);
2929 return dos * 1.e21;
2930 }
2931
2932 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2933 << " Band index (" << band << ") out of range.\n";
2934 return 0.;
2935}
2936
2937void MediumSilicon::ComputeSecondaries(const double e0, double& ee,
2938 double& eh) {
2939 const int nPointsValence = m_fbDosValence.size();
2940 const int nPointsConduction = m_fbDosConduction.size();
2941 const double widthValenceBand = m_eStepDos * nPointsValence;
2942 const double widthConductionBand = m_eStepDos * nPointsConduction;
2943
2944 bool ok = false;
2945 while (!ok) {
2946 // Sample a hole energy according to the valence band DOS.
2947 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2948 int ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2949 while (RndmUniform() > m_fbDosValence[ih] / m_fbDosMaxV) {
2950 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2951 ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2952 }
2953 // Sample an electron energy according to the conduction band DOS.
2954 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2955 int ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2956 while (RndmUniform() > m_fbDosConduction[ie] / m_fbDosMaxC) {
2957 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2958 ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2959 }
2960 // Calculate the energy of the primary electron.
2961 const double ep = e0 - m_bandGap - eh - ee;
2962 if (ep < Small) continue;
2963 if (ep > 5.) return;
2964 // Check if the primary electron energy is consistent with the DOS.
2965 int ip = std::min(int(ep / m_eStepDos), nPointsConduction - 1);
2966 if (RndmUniform() > m_fbDosConduction[ip] / m_fbDosMaxC) continue;
2967 ok = true;
2968 }
2969}
2970
2971void MediumSilicon::InitialiseDensityOfStates() {
2972 m_eStepDos = 0.1;
2973
2974 m_fbDosValence = {
2975 {0., 1.28083, 2.08928, 2.70763, 3.28095, 3.89162, 4.50547, 5.15043,
2976 5.89314, 6.72667, 7.67768, 8.82725, 10.6468, 12.7003, 13.7457, 14.0263,
2977 14.2731, 14.5527, 14.8808, 15.1487, 15.4486, 15.7675, 16.0519, 16.4259,
2978 16.7538, 17.0589, 17.3639, 17.6664, 18.0376, 18.4174, 18.2334, 16.7552,
2979 15.1757, 14.2853, 13.6516, 13.2525, 12.9036, 12.7203, 12.6104, 12.6881,
2980 13.2862, 14.0222, 14.9366, 13.5084, 9.77808, 6.15266, 3.47839, 2.60183,
2981 2.76747, 3.13985, 3.22524, 3.29119, 3.40868, 3.6118, 3.8464, 4.05776,
2982 4.3046, 4.56219, 4.81553, 5.09909, 5.37616, 5.67297, 6.04611, 6.47252,
2983 6.9256, 7.51254, 8.17923, 8.92351, 10.0309, 11.726, 16.2853, 18.2457,
2984 12.8879, 7.86019, 6.02275, 5.21777, 4.79054, 3.976, 3.11855, 2.46854,
2985 1.65381, 0.830278, 0.217735}};
2986
2987 m_fbDosConduction = {
2988 {0., 1.5114, 2.71026, 3.67114, 4.40173, 5.05025, 5.6849, 6.28358,
2989 6.84628, 7.43859, 8.00204, 8.80658, 9.84885, 10.9579, 12.0302, 13.2051,
2990 14.6948, 16.9879, 18.4492, 18.1933, 17.6747, 16.8135, 15.736, 14.4965,
2991 13.1193, 12.1817, 12.6109, 15.3148, 19.4936, 23.0093, 24.4106, 22.2834,
2992 19.521, 18.9894, 18.8015, 17.9363, 17.0252, 15.9871, 14.8486, 14.3797,
2993 14.2426, 14.3571, 14.7271, 14.681, 14.3827, 14.2789, 14.144, 14.1684,
2994 14.1418, 13.9237, 13.7558, 13.5691, 13.4567, 13.2693, 12.844, 12.4006,
2995 12.045, 11.7729, 11.3607, 11.14, 11.0586, 10.5475, 9.73786, 9.34423,
2996 9.4694, 9.58071, 9.6967, 9.84854, 10.0204, 9.82705, 9.09102, 8.30665,
2997 7.67306, 7.18925, 6.79675, 6.40713, 6.21687, 6.33267, 6.5223, 6.17877,
2998 5.48659, 4.92208, 4.44239, 4.02941, 3.5692, 3.05953, 2.6428, 2.36979,
2999 2.16273, 2.00627, 1.85206, 1.71265, 1.59497, 1.46681, 1.34913, 1.23951,
3000 1.13439, 1.03789, 0.924155, 0.834962, 0.751017}};
3001
3002 auto it = std::max_element(m_fbDosValence.begin(), m_fbDosValence.end());
3003 m_fbDosMaxV = *it;
3004 it = std::max_element(m_fbDosConduction.begin(), m_fbDosConduction.end());
3005 m_fbDosMaxC = *it;
3006}
3007}
bool GetDielectricFunction(const double e, double &eps1, double &eps2, const unsigned int i=0) override
Get the complex dielectric function at a given energy.
void SetLatticeMobilityModelReggiani()
Calculate the lattice mobility using the Reggiani model.
void SetDoping(const char type, const double c)
Set doping concentration [cm-3] and type ('i', 'n', 'p').
void GetDoping(char &type, double &c) const
Retrieve doping concentration.
bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha) override
Ionisation coefficient [cm-1].
void SetHighFieldMobilityModelReggiani()
Parameterize the high-field mobility using the Reggiani model.
double GetElectronCollisionRate(const double e, const int band) override
Collision rate [ns-1] for given electron energy.
bool HoleVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz) override
Drift velocity [cm / ns].
unsigned int GetNumberOfLevels() const
void SetSaturationVelocityModelReggiani()
Calculate the saturation velocities using the Reggiani model.
bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta) override
Attachment coefficient [cm-1].
double GetValenceBandDensityOfStates(const double e, const int band=-1)
bool GetElectronCollision(const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< std::pair< int, double > > &secondaries, int &ndxc, int &band) override
Sample the collision type. Update energy and direction vector.
void SetDopingMobilityModelMasetti()
Use the Masetti model for the doping-dependence of the mobility (default).
void SetTrapDensity(const double n)
Trap density [cm-3], by default set to zero.
void SetHighFieldMobilityModelMinimos()
Parameterize the high-field mobility using the Minimos model.
MediumSilicon()
Constructor.
void SetSaturationVelocity(const double vsate, const double vsath)
Specify the saturation velocities of electrons and holes.
bool ElectronVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz) override
Drift velocity [cm / ns].
void ComputeSecondaries(const double e0, double &ee, double &eh)
double GetConductionBandDensityOfStates(const double e, const int band=0)
bool SetMaxElectronEnergy(const double e)
bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha) override
Ionisation coefficient [cm-1].
void SetHighFieldMobilityModelConstant()
Make the velocity proportional to the electric field (no saturation).
void SetLatticeMobilityModelSentaurus()
Calculate the lattice mobility using the Sentaurus model (default).
void SetSaturationVelocityModelCanali()
Calculate the saturation velocities using the Canali model (default).
void SetHighFieldMobilityModelCanali()
Parameterize the high-field mobility using the Canali model (default).
void GetElectronMomentum(const double e, double &px, double &py, double &pz, int &band) override
void SetImpactIonisationModelVanOverstraetenDeMan()
Calculate α using the van Overstraeten-de Man model (default).
void SetTrappingTime(const double etau, const double htau)
Set time constant for trapping of electrons and holes [ns].
void SetSaturationVelocityModelMinimos()
Calculate the saturation velocities using the Minimos model.
void SetLowFieldMobility(const double mue, const double muh)
Specify the low field values of the electron and hole mobilities.
bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta) override
Attachment coefficient [cm-1].
double GetElectronEnergy(const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0) override
Dispersion relation (energy vs. wave vector)
bool GetOpticalDataRange(double &emin, double &emax, const unsigned int i=0) override
Get the energy range [eV] of the available optical data.
double GetElectronNullCollisionRate(const int band) override
Null-collision rate [ns-1].
unsigned int GetNumberOfElectronCollisions() const
void SetLatticeMobilityModelMinimos()
Calculate the lattice mobility using the Minimos model.
void SetTrapCrossSection(const double ecs, const double hcs)
Trapping cross-sections for electrons and holes.
int GetElectronBandPopulation(const int band)
void SetImpactIonisationModelGrant()
Calculate α using the Grant model.
void SetDopingMobilityModelMinimos()
Use the Minimos model for the doping-dependence of the mobility.
void SetImpactIonisationModelMassey()
Calculate α using the Massey model.
unsigned int GetNumberOfElectronBands() const
Abstract base class for media.
Definition: Medium.hh:13
double m_density
Definition: Medium.hh:527
void SetTemperature(const double t)
Set the temperature [K].
Definition: Medium.cc:71
bool m_microscopic
Definition: Medium.hh:531
virtual bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:534
double m_fano
Definition: Medium.hh:537
virtual bool HoleVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:513
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition: Medium.hh:558
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition: Medium.hh:570
virtual void EnableDrift(const bool on=true)
Switch electron/ion/hole on/off.
Definition: Medium.hh:67
double GetDielectricConstant() const
Get the relative static dielectric constant.
Definition: Medium.hh:42
std::string m_name
Definition: Medium.hh:513
virtual bool ElectronVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Drift velocity [cm / ns].
Definition: Medium.cc:379
virtual void EnablePrimaryIonisation(const bool on=true)
Make the medium ionisable or non-ionisable.
Definition: Medium.hh:69
virtual void SetAtomicNumber(const double z)
Set the effective atomic number.
Definition: Medium.cc:114
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition: Medium.hh:553
void SetDielectricConstant(const double eps)
Set the relative static dielectric constant.
Definition: Medium.cc:91
virtual bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:403
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition: Medium.hh:559
virtual void SetMassDensity(const double rho)
Set the mass density [g/cm3].
Definition: Medium.cc:144
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition: Medium.hh:565
virtual bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:416
std::vector< std::vector< std::vector< double > > > m_hAtt
Definition: Medium.hh:571
virtual void SetAtomicWeight(const double a)
Set the effective atomic weight.
Definition: Medium.cc:124
std::string m_className
Definition: Medium.hh:506
bool m_isChanged
Definition: Medium.hh:540
virtual bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:547
double m_temperature
Definition: Medium.hh:515
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
void ltrim(std::string &line)
Definition: Utilities.hh:10
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:107
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
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314