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