Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::MediumSilicon Class Reference

Solid crystalline silicon More...

#include <MediumSilicon.hh>

+ Inheritance diagram for Garfield::MediumSilicon:

Public Member Functions

 MediumSilicon ()
 Constructor.
 
virtual ~MediumSilicon ()
 Destructor.
 
bool IsSemiconductor () const
 
void SetDoping (const char type, const double c)
 Doping concentration [cm-3] and type ('i', 'n', 'p')
 
void GetDoping (char &type, double &c) const
 
void SetTrapCrossSection (const double ecs, const double hcs)
 Trapping cross-sections for electrons and holes.
 
void SetTrapDensity (const double n)
 
void SetTrappingTime (const double etau, const double htau)
 
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)
 
bool ElectronTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
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)
 
bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
void SetLowFieldMobility (const double mue, const double muh)
 
void SetLatticeMobilityModelMinimos ()
 
void SetLatticeMobilityModelSentaurus ()
 
void SetLatticeMobilityModelReggiani ()
 
void SetDopingMobilityModelMinimos ()
 
void SetDopingMobilityModelMasetti ()
 
void SetSaturationVelocity (const double vsate, const double vsath)
 
void SetSaturationVelocityModelMinimos ()
 
void SetSaturationVelocityModelCanali ()
 
void SetSaturationVelocityModelReggiani ()
 
void SetHighFieldMobilityModelMinimos ()
 
void SetHighFieldMobilityModelCanali ()
 
void SetHighFieldMobilityModelReggiani ()
 
void SetHighFieldMobilityModelConstant ()
 
void SetImpactIonisationModelVanOverstraetenDeMan ()
 
void SetImpactIonisationModelGrant ()
 
void SetDiffusionScaling (const double d)
 
bool SetMaxElectronEnergy (const double e)
 
double GetMaxElectronEnergy () const
 
bool Initialise ()
 
void EnableScatteringRateOutput ()
 
void DisableScatteringRateOutput ()
 
void EnableNonParabolicity ()
 
void DisableNonParabolicity ()
 
void EnableFullBandDensityOfStates ()
 
void DisableFullBandDensityOfStates ()
 
void EnableAnisotropy ()
 
void DisableAnisotropy ()
 
double GetElectronEnergy (const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
 
void GetElectronMomentum (const double e, double &px, double &py, double &pz, int &band)
 
double GetElectronNullCollisionRate (const int band)
 
double GetElectronCollisionRate (const double e, const int band)
 
bool GetElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, int &nion, int &ndxc, int &band)
 
unsigned int GetNumberOfIonisationProducts () const
 
bool GetIonisationProduct (const unsigned int i, int &type, double &energy) const
 
double GetConductionBandDensityOfStates (const double e, const int band=0)
 
double GetValenceBandDensityOfStates (const double e, const int band=-1)
 
void ResetCollisionCounters ()
 
unsigned int GetNumberOfElectronCollisions () const
 
unsigned int GetNumberOfLevels () const
 
unsigned int GetNumberOfElectronCollisions (const unsigned int level) const
 
unsigned int GetNumberOfElectronBands () const
 
int GetElectronBandPopulation (const int band)
 
bool GetOpticalDataRange (double &emin, double &emax, const unsigned int i=0)
 
bool GetDielectricFunction (const double e, double &eps1, double &eps2, const unsigned int i=0)
 
void ComputeSecondaries (const double e0, double &ee, double &eh)
 
- Public Member Functions inherited from Garfield::Medium
 Medium ()
 
virtual ~Medium ()
 
int GetId () const
 
const std::string & GetName () const
 
virtual bool IsGas () const
 
virtual bool IsSemiconductor () const
 
void SetTemperature (const double t)
 
double GetTemperature () const
 
void SetPressure (const double p)
 
double GetPressure () const
 
void SetDielectricConstant (const double eps)
 
double GetDielectricConstant () const
 
unsigned int GetNumberOfComponents () const
 
virtual void GetComponent (const unsigned int i, std::string &label, double &f)
 
virtual void SetAtomicNumber (const double z)
 
virtual double GetAtomicNumber () const
 
virtual void SetAtomicWeight (const double a)
 
virtual double GetAtomicWeight () const
 
virtual void SetNumberDensity (const double n)
 
virtual double GetNumberDensity () const
 
virtual void SetMassDensity (const double rho)
 
virtual double GetMassDensity () const
 
virtual void EnableDrift ()
 
void DisableDrift ()
 
virtual void EnablePrimaryIonisation ()
 
void DisablePrimaryIonisation ()
 
bool IsDriftable () const
 
bool IsMicroscopic () const
 
bool IsIonisable () const
 
void SetW (const double w)
 
double GetW ()
 
void SetFanoFactor (const double f)
 
double GetFanoFactor ()
 
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)
 
virtual bool ElectronDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
 
virtual bool ElectronDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double cov[3][3])
 
virtual bool ElectronTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
virtual bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
virtual bool ElectronLorentzAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
 
virtual double GetElectronEnergy (const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
 
virtual void GetElectronMomentum (const double e, double &px, double &py, double &pz, int &band)
 
virtual double GetElectronNullCollisionRate (const int band=0)
 
virtual double GetElectronCollisionRate (const double e, const int band=0)
 
virtual bool GetElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, int &nion, int &ndxc, int &band)
 
virtual unsigned int GetNumberOfIonisationProducts () const
 
virtual bool GetIonisationProduct (const unsigned int i, int &type, double &energy) const
 
virtual unsigned int GetNumberOfDeexcitationProducts () const
 
virtual bool GetDeexcitationProduct (const unsigned int i, double &t, double &s, int &type, double &energy) const
 
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)
 
virtual bool HoleDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
 
virtual bool HoleDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double cov[3][3])
 
virtual bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
virtual bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
virtual bool IonVelocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
 
virtual bool IonDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
 
virtual bool IonDissociation (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &diss)
 
void SetFieldGrid (double emin, double emax, int ne, bool logE, double bmin=0., double bmax=0., int nb=1, double amin=0., double amax=0., int na=1)
 
void SetFieldGrid (const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles)
 
void GetFieldGrid (std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles)
 
bool GetElectronVelocityE (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetElectronVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetElectronVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetElectronLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 
bool GetElectronTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 
bool GetElectronTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &alpha)
 
bool GetElectronAttachment (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &eta)
 
bool GetElectronLorentzAngle (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &lor)
 
bool GetHoleVelocityE (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetHoleVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetHoleVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetHoleLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 
bool GetHoleTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 
bool GetHoleTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &alpha)
 
bool GetHoleAttachment (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &eta)
 
bool GetIonMobility (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &mu)
 
bool GetIonLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 
bool GetIonTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 
bool GetIonDissociation (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &diss)
 
void ResetElectronVelocity ()
 
void ResetElectronDiffusion ()
 
void ResetElectronTownsend ()
 
void ResetElectronAttachment ()
 
void ResetElectronLorentzAngle ()
 
void ResetHoleVelocity ()
 
void ResetHoleDiffusion ()
 
void ResetHoleTownsend ()
 
void ResetHoleAttachment ()
 
void ResetIonMobility ()
 
void ResetIonDiffusion ()
 
void ResetIonDissociation ()
 
bool SetIonMobility (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double mu)
 
bool SetIonMobility (const std::vector< double > &fields, const std::vector< double > &mobilities)
 
void SetExtrapolationMethodVelocity (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodDiffusion (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodTownsend (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodAttachment (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodIonMobility (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodIonDissociation (const std::string &extrLow, const std::string &extrHigh)
 
void SetInterpolationMethodVelocity (const unsigned int intrp)
 
void SetInterpolationMethodDiffusion (const unsigned int intrp)
 
void SetInterpolationMethodTownsend (const unsigned int intrp)
 
void SetInterpolationMethodAttachment (const unsigned int intrp)
 
void SetInterpolationMethodIonMobility (const unsigned int intrp)
 
void SetInterpolationMethodIonDissociation (const unsigned int intrp)
 
virtual double ScaleElectricField (const double e) const
 
virtual double UnScaleElectricField (const double e) const
 
virtual double ScaleVelocity (const double v) const
 
virtual double ScaleDiffusion (const double d) const
 
virtual double ScaleDiffusionTensor (const double d) const
 
virtual double ScaleTownsend (const double alpha) const
 
virtual double ScaleAttachment (const double eta) const
 
virtual double ScaleLorentzAngle (const double lor) const
 
virtual double ScaleDissociation (const double diss) const
 
virtual bool GetOpticalDataRange (double &emin, double &emax, const unsigned int i=0)
 
virtual bool GetDielectricFunction (const double e, double &eps1, double &eps2, const unsigned int i=0)
 
virtual bool GetPhotoAbsorptionCrossSection (const double e, double &sigma, const unsigned int i=0)
 
virtual double GetPhotonCollisionRate (const double e)
 
virtual bool GetPhotonCollision (const double e, int &type, int &level, double &e1, double &ctheta, int &nsec, double &esec)
 
void EnableDebugging ()
 
void DisableDebugging ()
 

Additional Inherited Members

- Protected Member Functions inherited from Garfield::Medium
double GetAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const double e, const double b) const
 
double Interpolate1D (const double e, const std::vector< double > &table, const std::vector< double > &fields, const unsigned int intpMeth, const int jExtr, const int iExtr)
 
bool GetExtrapolationIndex (std::string extrStr, unsigned int &extrNb)
 
void CloneTable (std::vector< std::vector< std::vector< double > > > &tab, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int intp, const unsigned int extrLow, const unsigned int extrHigh, const double init, const std::string &label)
 
void CloneTensor (std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const unsigned int n, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int intp, const unsigned int extrLow, const unsigned int extrHigh, const double init, const std::string &label)
 
void InitParamArrays (const unsigned int eRes, const unsigned int bRes, const unsigned int aRes, std::vector< std::vector< std::vector< double > > > &tab, const double val)
 
void InitParamTensor (const unsigned int eRes, const unsigned int bRes, const unsigned int aRes, const unsigned int tRes, std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const double val)
 
- Protected Attributes inherited from Garfield::Medium
std::string m_className
 
int m_id
 
std::string m_name
 
double m_temperature
 
double m_pressure
 
double m_epsilon
 
unsigned int m_nComponents
 
double m_z
 
double m_a
 
double m_density
 
bool m_driftable
 
bool m_microscopic
 
bool m_ionisable
 
double m_w
 
double m_fano
 
bool m_isChanged
 
bool m_debug
 
std::vector< double > m_eFields
 
std::vector< double > m_bFields
 
std::vector< double > m_bAngles
 
bool m_map2d
 
bool m_hasElectronVelocityE
 
bool m_hasElectronVelocityB
 
bool m_hasElectronVelocityExB
 
bool m_hasElectronDiffLong
 
bool m_hasElectronDiffTrans
 
bool m_hasElectronDiffTens
 
bool m_hasElectronAttachment
 
bool m_hasElectronLorentzAngle
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityE
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityExB
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityB
 
std::vector< std::vector< std::vector< double > > > tabElectronDiffLong
 
std::vector< std::vector< std::vector< double > > > tabElectronDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabElectronTownsend
 
std::vector< std::vector< std::vector< double > > > tabElectronAttachment
 
std::vector< std::vector< std::vector< double > > > tabElectronLorentzAngle
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabElectronDiffTens
 
bool m_hasHoleVelocityE
 
bool m_hasHoleVelocityB
 
bool m_hasHoleVelocityExB
 
bool m_hasHoleDiffLong
 
bool m_hasHoleDiffTrans
 
bool m_hasHoleDiffTens
 
bool m_hasHoleTownsend
 
bool m_hasHoleAttachment
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityE
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityExB
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityB
 
std::vector< std::vector< std::vector< double > > > tabHoleDiffLong
 
std::vector< std::vector< std::vector< double > > > tabHoleDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabHoleTownsend
 
std::vector< std::vector< std::vector< double > > > tabHoleAttachment
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabHoleDiffTens
 
bool m_hasIonMobility
 
bool m_hasIonDiffLong
 
bool m_hasIonDiffTrans
 
bool m_hasIonDissociation
 
std::vector< std::vector< std::vector< double > > > tabIonMobility
 
std::vector< std::vector< std::vector< double > > > tabIonDiffLong
 
std::vector< std::vector< std::vector< double > > > tabIonDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabIonDissociation
 
int thrElectronTownsend
 
int thrElectronAttachment
 
int thrHoleTownsend
 
int thrHoleAttachment
 
int thrIonDissociation
 
unsigned int m_extrLowVelocity
 
unsigned int m_extrHighVelocity
 
unsigned int m_extrLowDiffusion
 
unsigned int m_extrHighDiffusion
 
unsigned int m_extrLowTownsend
 
unsigned int m_extrHighTownsend
 
unsigned int m_extrLowAttachment
 
unsigned int m_extrHighAttachment
 
unsigned int m_extrLowLorentzAngle
 
unsigned int m_extrHighLorentzAngle
 
unsigned int m_extrLowMobility
 
unsigned int m_extrHighMobility
 
unsigned int m_extrLowDissociation
 
unsigned int m_extrHighDissociation
 
unsigned int m_intpVelocity
 
unsigned int m_intpDiffusion
 
unsigned int m_intpTownsend
 
unsigned int m_intpAttachment
 
unsigned int m_intpLorentzAngle
 
unsigned int m_intpMobility
 
unsigned int m_intpDissociation
 
- Static Protected Attributes inherited from Garfield::Medium
static int m_idCounter = -1
 

Detailed Description

Solid crystalline silicon

Definition at line 9 of file MediumSilicon.hh.

Constructor & Destructor Documentation

◆ MediumSilicon()

Garfield::MediumSilicon::MediumSilicon ( )

Constructor.

Definition at line 14 of file MediumSilicon.cc.

15 : Medium(),
16 diffScale(1.0),
17 m_bandGap(1.12),
18 m_dopingType('i'),
19 m_dopingConcentration(0.),
20 m_mLongX(0.916),
21 m_mTransX(0.191),
22 m_mLongL(1.59),
23 m_mTransL(0.12),
24 m_alphaX(0.5),
25 m_alphaL(0.5),
26 m_eLatticeMobility(1.35e-6),
27 m_hLatticeMobility(0.45e-6),
28 m_eMobility(1.35e-6),
29 m_hMobility(0.45e-6),
30 m_eBetaCanali(1.109),
31 m_hBetaCanali(1.213),
32 m_eBetaCanaliInv(1. / 1.109),
33 m_hBetaCanaliInv(1. / 1.213),
34 m_eSatVel(1.02e-2),
35 m_hSatVel(0.72e-2),
36 m_eHallFactor(1.15),
37 m_hHallFactor(0.7),
38 m_eTrapCs(1.e-15),
39 m_hTrapCs(1.e-15),
40 m_eTrapDensity(1.e13),
41 m_hTrapDensity(1.e13),
42 m_eTrapTime(0.),
43 m_hTrapTime(0.),
44 m_trappingModel(0),
45 m_eImpactA0(3.318e5),
46 m_eImpactA1(0.703e6),
47 m_eImpactA2(0.),
48 m_eImpactB0(1.135e6),
49 m_eImpactB1(1.231e6),
50 m_eImpactB2(0.),
51 m_hImpactA0(1.582e6),
52 m_hImpactA1(0.671e6),
53 m_hImpactB0(2.036e6),
54 m_hImpactB1(1.693e6),
55 m_hasUserMobility(false),
56 m_hasUserSaturationVelocity(false),
57 m_latticeMobilityModel(LatticeMobilityModelSentaurus),
58 m_dopingMobilityModel(DopingMobilityModelMasetti),
59 m_saturationVelocityModel(SaturationVelocityModelCanali),
60 m_highFieldMobilityModel(HighFieldMobilityModelCanali),
61 m_impactIonisationModel(ImpactIonisationModelVanOverstraeten),
62 m_useCfOutput(false),
63 m_useNonParabolicity(true),
64 m_useFullBandDos(true),
65 m_useAnisotropy(true),
66 m_eFinalXL(4.),
67 m_eStepXL(m_eFinalXL / nEnergyStepsXL),
68 m_eFinalG(10.),
69 m_eStepG(m_eFinalG / nEnergyStepsG),
70 m_eFinalV(8.5),
71 m_eStepV(m_eFinalV / nEnergyStepsV),
72 m_nLevelsX(0),
73 m_nLevelsL(0),
74 m_nLevelsG(0),
75 m_nLevelsV(0),
76 m_nValleysX(6),
77 m_nValleysL(8),
78 m_eMinL(1.05),
79 m_eMinG(2.24),
80 m_ieMinL(0),
81 m_ieMinG(0),
82 m_opticalDataFile("OpticalData_Si.txt") {
83
84 m_className = "MediumSilicon";
85 m_name = "Si";
86
87 SetTemperature(300.);
89 SetAtomicNumber(14.);
90 SetAtomicWeight(28.0855);
91 SetMassDensity(2.329);
92
95 m_microscopic = true;
96
97 m_w = 3.6;
98 m_fano = 0.11;
99
100 m_cfTotElectronsX.clear();
101 m_cfElectronsX.clear();
102 m_energyLossElectronsX.clear();
103 m_scatTypeElectronsX.clear();
104
105 m_cfTotElectronsL.clear();
106 m_cfElectronsL.clear();
107 m_energyLossElectronsL.clear();
108 m_scatTypeElectronsL.clear();
109
110 m_cfTotElectronsG.clear();
111 m_cfElectronsG.clear();
112 m_energyLossElectronsG.clear();
113 m_scatTypeElectronsG.clear();
114
115 m_ieMinL = int(m_eMinL / m_eStepXL) + 1;
116 m_ieMinG = int(m_eMinG / m_eStepG) + 1;
117
118 m_cfTotHoles.clear();
119 m_cfHoles.clear();
120 m_energyLossHoles.clear();
121 m_scatTypeHoles.clear();
122
123 // Load the density of states table.
124 InitialiseDensityOfStates();
125
126 // Initialize the collision counters.
127 m_nCollElectronAcoustic = m_nCollElectronOptical = 0;
128 m_nCollElectronIntervalley = 0;
129 m_nCollElectronImpurity = 0;
130 m_nCollElectronIonisation = 0;
131 m_nCollElectronDetailed.clear();
132 m_nCollElectronBand.clear();
133
134 m_ionProducts.clear();
135}
void SetTemperature(const double t)
Definition: Medium.cc:104
bool m_microscopic
Definition: Medium.hh:319
double m_fano
Definition: Medium.hh:323
std::string m_name
Definition: Medium.hh:301
virtual void SetAtomicNumber(const double z)
Definition: Medium.cc:154
void SetDielectricConstant(const double eps)
Definition: Medium.cc:126
virtual void EnableDrift()
Definition: Medium.hh:52
virtual void SetMassDensity(const double rho)
Definition: Medium.cc:187
virtual void EnablePrimaryIonisation()
Definition: Medium.hh:54
virtual void SetAtomicWeight(const double a)
Definition: Medium.cc:165
std::string m_className
Definition: Medium.hh:294

◆ ~MediumSilicon()

virtual Garfield::MediumSilicon::~MediumSilicon ( )
inlinevirtual

Destructor.

Definition at line 15 of file MediumSilicon.hh.

15{}

Member Function Documentation

◆ ComputeSecondaries()

void Garfield::MediumSilicon::ComputeSecondaries ( const double  e0,
double &  ee,
double &  eh 
)

Definition at line 3182 of file MediumSilicon.cc.

3183 {
3184
3185 const int nPointsValence = m_fbDosValence.size();
3186 const int nPointsConduction = m_fbDosConduction.size();
3187 const double widthValenceBand = m_eStepDos * nPointsValence;
3188 const double widthConductionBand = m_eStepDos * nPointsConduction;
3189
3190 bool ok = false;
3191 while (!ok) {
3192 // Sample a hole energy according to the valence band DOS.
3193 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
3194 int ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
3195 while (RndmUniform() > m_fbDosValence[ih] / m_fbDosMaxV) {
3196 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
3197 ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
3198 }
3199 // Sample an electron energy according to the conduction band DOS.
3200 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
3201 int ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
3202 while (RndmUniform() > m_fbDosConduction[ie] / m_fbDosMaxC) {
3203 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
3204 ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
3205 }
3206 // Calculate the energy of the primary electron.
3207 const double ep = e0 - m_bandGap - eh - ee;
3208 if (ep < Small) continue;
3209 if (ep > 5.) return;
3210 // Check if the primary electron energy is consistent with the DOS.
3211 int ip = std::min(int(ep / m_eStepDos), nPointsConduction - 1);
3212 if (RndmUniform() > m_fbDosConduction[ip] / m_fbDosMaxC) continue;
3213 ok = true;
3214 }
3215}
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition: Random.hh:17

Referenced by GetElectronCollision().

◆ DisableAnisotropy()

void Garfield::MediumSilicon::DisableAnisotropy ( )
inline

Definition at line 92 of file MediumSilicon.hh.

92{ m_useAnisotropy = false; }

◆ DisableFullBandDensityOfStates()

void Garfield::MediumSilicon::DisableFullBandDensityOfStates ( )
inline

Definition at line 90 of file MediumSilicon.hh.

90{ m_useFullBandDos = false; }

◆ DisableNonParabolicity()

void Garfield::MediumSilicon::DisableNonParabolicity ( )
inline

Definition at line 88 of file MediumSilicon.hh.

88{ m_useNonParabolicity = false; }

◆ DisableScatteringRateOutput()

void Garfield::MediumSilicon::DisableScatteringRateOutput ( )
inline

Definition at line 85 of file MediumSilicon.hh.

85{ m_useCfOutput = false; }

◆ ElectronAttachment()

bool Garfield::MediumSilicon::ElectronAttachment ( const double  ex,
const double  ey,
const double  ez,
const double  bx,
const double  by,
const double  bz,
double &  eta 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 329 of file MediumSilicon.cc.

332 {
333
334 eta = 0.;
335 if (m_isChanged) {
336 if (!UpdateTransportParameters()) {
337 std::cerr << m_className << "::ElectronAttachment:\n"
338 << " Error calculating the transport parameters.\n";
339 return false;
340 }
341 m_isChanged = false;
342 }
343
345 // Interpolation in user table.
346 return Medium::ElectronAttachment(ex, ey, ez, bx, by, bz, eta);
347 }
348
349 switch (m_trappingModel) {
350 case 0:
351 eta = m_eTrapCs * m_eTrapDensity;
352 break;
353 case 1:
354 double vx, vy, vz;
355 ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
356 eta = m_eTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
357 if (eta > 0.) eta = 1. / eta;
358 break;
359 default:
360 std::cerr << m_className << "::ElectronAttachment:\n"
361 << " Unknown model. Program bug!\n";
362 return false;
363 break;
364 }
365
366 return true;
367}
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)
virtual bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Definition: Medium.cc:544
bool m_isChanged
Definition: Medium.hh:326
bool m_hasElectronAttachment
Definition: Medium.hh:341
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

◆ ElectronTownsend()

bool Garfield::MediumSilicon::ElectronTownsend ( const double  ex,
const double  ey,
const double  ez,
const double  bx,
const double  by,
const double  bz,
double &  alpha 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 292 of file MediumSilicon.cc.

295 {
296
297 alpha = 0.;
298 if (m_isChanged) {
299 if (!UpdateTransportParameters()) {
300 std::cerr << m_className << "::ElectronTownsend:\n"
301 << " Error calculating the transport parameters.\n";
302 return false;
303 }
304 m_isChanged = false;
305 }
306
307 if (!tabElectronTownsend.empty()) {
308 // Interpolation in user table.
309 return Medium::ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
310 }
311
312 const double e = sqrt(ex * ex + ey * ey + ez * ez);
313
314 switch (m_impactIonisationModel) {
315 case ImpactIonisationModelVanOverstraeten:
316 return ElectronImpactIonisationVanOverstraetenDeMan(e, alpha);
317 break;
318 case ImpactIonisationModelGrant:
319 return ElectronImpactIonisationGrant(e, alpha);
320 break;
321 default:
322 std::cerr << m_className << "::ElectronTownsend:\n"
323 << " Unknown model. Program bug!\n";
324 break;
325 }
326 return false;
327}
virtual bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Definition: Medium.cc:490
std::vector< std::vector< std::vector< double > > > tabElectronTownsend
Definition: Medium.hh:348

◆ ElectronVelocity()

bool Garfield::MediumSilicon::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 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 234 of file MediumSilicon.cc.

237 {
238
239 vx = vy = vz = 0.;
240 if (m_isChanged) {
241 if (!UpdateTransportParameters()) {
242 std::cerr << m_className << "::ElectronVelocity:\n"
243 << " Error calculating the transport parameters.\n";
244 return false;
245 }
246 m_isChanged = false;
247 }
248
250 // Interpolation in user table.
251 return Medium::ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
252 }
253
254 const double e = sqrt(ex * ex + ey * ey + ez * ez);
255
256 // Calculate the mobility
257 double mu;
258 switch (m_highFieldMobilityModel) {
259 case HighFieldMobilityModelMinimos:
260 ElectronMobilityMinimos(e, mu);
261 break;
262 case HighFieldMobilityModelCanali:
263 ElectronMobilityCanali(e, mu);
264 break;
265 case HighFieldMobilityModelReggiani:
266 ElectronMobilityReggiani(e, mu);
267 break;
268 default:
269 mu = m_eMobility;
270 break;
271 }
272 mu = -mu;
273
274 const double b = sqrt(bx * bx + by * by + bz * bz);
275 if (b < Small) {
276 vx = mu * ex;
277 vy = mu * ey;
278 vz = mu * ez;
279 } else {
280 // Hall mobility
281 const double muH = m_eHallFactor * mu;
282 const double eb = bx * ex + by * ey + bz * ez;
283 const double nom = 1. + pow(muH * b, 2);
284 // Compute the drift velocity using the Langevin equation.
285 vx = mu * (ex + muH * (ey * bz - ez * by) + muH * muH * bx * eb) / nom;
286 vy = mu * (ey + muH * (ez * bx - ex * bz) + muH * muH * by * eb) / nom;
287 vz = mu * (ez + muH * (ex * by - ey * bx) + muH * muH * bz * eb) / nom;
288 }
289 return true;
290}
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)
Definition: Medium.cc:204
bool m_hasElectronVelocityE
Definition: Medium.hh:339
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337

Referenced by ElectronAttachment().

◆ EnableAnisotropy()

void Garfield::MediumSilicon::EnableAnisotropy ( )
inline

Definition at line 91 of file MediumSilicon.hh.

91{ m_useAnisotropy = true; }

◆ EnableFullBandDensityOfStates()

void Garfield::MediumSilicon::EnableFullBandDensityOfStates ( )
inline

Definition at line 89 of file MediumSilicon.hh.

89{ m_useFullBandDos = true; }

◆ EnableNonParabolicity()

void Garfield::MediumSilicon::EnableNonParabolicity ( )
inline

Definition at line 87 of file MediumSilicon.hh.

87{ m_useNonParabolicity = true; }

◆ EnableScatteringRateOutput()

void Garfield::MediumSilicon::EnableScatteringRateOutput ( )
inline

Definition at line 84 of file MediumSilicon.hh.

84{ m_useCfOutput = true; }

◆ GetConductionBandDensityOfStates()

double Garfield::MediumSilicon::GetConductionBandDensityOfStates ( const double  e,
const int  band = 0 
)

Definition at line 3048 of file MediumSilicon.cc.

3049 {
3050 if (band < 0) {
3051 int iE = int(e / m_eStepDos);
3052 const int nPoints = m_fbDosConduction.size();
3053 if (iE >= nPoints || iE < 0) {
3054 return 0.;
3055 } else if (iE == nPoints - 1) {
3056 return m_fbDosConduction[nPoints - 1];
3057 }
3058
3059 const double dos =
3060 m_fbDosConduction[iE] +
3061 (m_fbDosConduction[iE + 1] - m_fbDosConduction[iE]) * (e / m_eStepDos - iE);
3062 return dos * 1.e21;
3063
3064 } else if (band < m_nValleysX) {
3065 // X valleys
3066 if (e <= 0.) return 0.;
3067 // Density-of-states effective mass (cube)
3068 const double md3 = pow(ElectronMass, 3) * m_mLongX * m_mTransX * m_mTransX;
3069
3070 if (m_useFullBandDos) {
3071 if (e < m_eMinL) {
3072 return GetConductionBandDensityOfStates(e, -1) / m_nValleysX;
3073 } else if (e < m_eMinG) {
3074 // Subtract the fraction of the full-band density of states
3075 // attributed to the L valleys.
3076 const double dosX =
3078 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL;
3079 return dosX / m_nValleysX;
3080 } else {
3081 // Subtract the fraction of the full-band density of states
3082 // attributed to the L valleys and the higher bands.
3083 const double dosX =
3085 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL -
3086 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
3087 if (dosX <= 0.) return 0.;
3088 return dosX / m_nValleysX;
3089 }
3090 } else if (m_useNonParabolicity) {
3091 const double alpha = m_alphaX;
3092 return sqrt(md3 * e * (1. + alpha * e) / 2.) * (1. + 2 * alpha * e) /
3093 (Pi2 * pow(HbarC, 3.));
3094 } else {
3095 return sqrt(md3 * e / 2.) / (Pi2 * pow(HbarC, 3.));
3096 }
3097 } else if (band < m_nValleysX + m_nValleysL) {
3098 // L valleys
3099 if (e <= m_eMinL) return 0.;
3100
3101 // Density-of-states effective mass (cube)
3102 const double md3 = pow(ElectronMass, 3) * m_mLongL * m_mTransL * m_mTransL;
3103 // Non-parabolicity parameter
3104 const double alpha = m_alphaL;
3105
3106 if (m_useFullBandDos) {
3107 // Energy up to which the non-parabolic approximation is used.
3108 const double ej = m_eMinL + 0.5;
3109 if (e <= ej) {
3110 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
3111 (1. + 2 * alpha * (e - m_eMinL)) / (Sqrt2 * Pi2 * pow(HbarC, 3.));
3112 } else {
3113 // Fraction of full-band density of states attributed to L valleys
3114 double fL = sqrt(md3 * (ej - m_eMinL) * (1. + alpha * (ej - m_eMinL))) *
3115 (1. + 2 * alpha * (ej - m_eMinL)) /
3116 (Sqrt2 * Pi2 * pow(HbarC, 3.));
3117 fL = m_nValleysL * fL / GetConductionBandDensityOfStates(ej, -1);
3118
3119 double dosXL = GetConductionBandDensityOfStates(e, -1);
3120 if (e > m_eMinG) {
3121 dosXL -= GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
3122 }
3123 if (dosXL <= 0.) return 0.;
3124 return fL * dosXL / 8.;
3125 }
3126 } else if (m_useNonParabolicity) {
3127 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
3128 (1. + 2 * alpha * (e - m_eMinL)) / (Sqrt2 * Pi2 * pow(HbarC, 3.));
3129 } else {
3130 return sqrt(md3 * (e - m_eMinL) / 2.) / (Pi2 * pow(HbarC, 3.));
3131 }
3132 } else if (band == m_nValleysX + m_nValleysL) {
3133 // Higher bands
3134 const double ej = 2.7;
3135 if (m_eMinG >= ej) {
3136 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
3137 << " Cannot determine higher band density-of-states.\n"
3138 << " Program bug. Check offset energy!\n";
3139 return 0.;
3140 }
3141 if (e < m_eMinG) {
3142 return 0.;
3143 } else if (e < ej) {
3144 // Coexistence of XL and higher bands.
3145 const double dj = GetConductionBandDensityOfStates(ej, -1);
3146 // Assume linear increase of density-of-states.
3147 return dj * (e - m_eMinG) / (ej - m_eMinG);
3148 } else {
3150 }
3151 }
3152
3153 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
3154 << " Band index (" << band << ") out of range.\n";
3155 return ElectronMass * sqrt(ElectronMass * e / 2.) / (Pi2 * pow(HbarC, 3.));
3156}
double GetConductionBandDensityOfStates(const double e, const int band=0)

Referenced by GetConductionBandDensityOfStates(), and GetElectronMomentum().

◆ GetDielectricFunction()

bool Garfield::MediumSilicon::GetDielectricFunction ( const double  e,
double &  eps1,
double &  eps2,
const unsigned int  i = 0 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 1363 of file MediumSilicon.cc.

1364 {
1365
1366 if (i != 0) {
1367 std::cerr << m_className << "::GetDielectricFunction:\n";
1368 std::cerr << " Medium has only one component.\n";
1369 return false;
1370 }
1371
1372 // Make sure the optical data table has been loaded.
1373 if (m_opticalDataTable.empty()) {
1374 if (!LoadOpticalData(m_opticalDataFile)) {
1375 std::cerr << m_className << "::GetDielectricFunction:\n";
1376 std::cerr << " Optical data table could not be loaded.\n";
1377 return false;
1378 }
1379 }
1380
1381 // Make sure the requested energy is within the range of the table.
1382 const double emin = m_opticalDataTable.front().energy;
1383 const double emax = m_opticalDataTable.back().energy;
1384 if (e < emin || e > emax) {
1385 std::cerr << m_className << "::GetDielectricFunction:\n"
1386 << " Requested energy (" << e << " eV) "
1387 << " is outside the range of the optical data table.\n"
1388 << " " << emin << " < E [eV] < " << emax << "\n";
1389 eps1 = eps2 = 0.;
1390 return false;
1391 }
1392
1393 // Locate the requested energy in the table.
1394 int iLow = 0;
1395 int iUp = m_opticalDataTable.size() - 1;
1396 int iM;
1397 while (iUp - iLow > 1) {
1398 iM = (iUp + iLow) >> 1;
1399 if (e >= m_opticalDataTable[iM].energy) {
1400 iLow = iM;
1401 } else {
1402 iUp = iM;
1403 }
1404 }
1405
1406 // Interpolate the real part of dielectric function.
1407 // Use linear interpolation if one of the values is negative,
1408 // Otherwise use log-log interpolation.
1409 const double logX0 = log(m_opticalDataTable[iLow].energy);
1410 const double logX1 = log(m_opticalDataTable[iUp].energy);
1411 const double logX = log(e);
1412 if (m_opticalDataTable[iLow].eps1 <= 0. ||
1413 m_opticalDataTable[iUp].eps1 <= 0.) {
1414 eps1 = m_opticalDataTable[iLow].eps1 +
1415 (e - m_opticalDataTable[iLow].energy) *
1416 (m_opticalDataTable[iUp].eps1 - m_opticalDataTable[iLow].eps1) /
1417 (m_opticalDataTable[iUp].energy - m_opticalDataTable[iLow].energy);
1418 } else {
1419 const double logY0 = log(m_opticalDataTable[iLow].eps1);
1420 const double logY1 = log(m_opticalDataTable[iUp].eps1);
1421 eps1 = logY0 + (logX - logX0) * (logY1 - logY0) / (logX1 - logX0);
1422 eps1 = exp(eps1);
1423 }
1424
1425 // Interpolate the imaginary part of dielectric function,
1426 // using log-log interpolation.
1427 const double logY0 = log(m_opticalDataTable[iLow].eps2);
1428 const double logY1 = log(m_opticalDataTable[iUp].eps2);
1429 eps2 = logY0 + (log(e) - logX0) * (logY1 - logY0) / (logX1 - logX0);
1430 eps2 = exp(eps2);
1431 return true;
1432}
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377

◆ GetDoping()

void Garfield::MediumSilicon::GetDoping ( char &  type,
double &  c 
) const

Definition at line 174 of file MediumSilicon.cc.

174 {
175
176 type = m_dopingType;
177 c = m_dopingConcentration;
178}

◆ GetElectronBandPopulation()

int Garfield::MediumSilicon::GetElectronBandPopulation ( const int  band)

Definition at line 1326 of file MediumSilicon.cc.

1326 {
1327
1328 const int nBands = m_nValleysX + m_nValleysL + 1;
1329 if (band < 0 || band >= nBands) {
1330 std::cerr << m_className << "::GetElectronBandPopulation:\n";
1331 std::cerr << " Band index (" << band << ") out of range.\n";
1332 return 0;
1333 }
1334 return m_nCollElectronBand[band];
1335}

◆ GetElectronCollision()

bool Garfield::MediumSilicon::GetElectronCollision ( const double  e,
int &  type,
int &  level,
double &  e1,
double &  dx,
double &  dy,
double &  dz,
int &  nion,
int &  ndxc,
int &  band 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 911 of file MediumSilicon.cc.

914 {
915
916 if (e > m_eFinalG) {
917 std::cerr << m_className << "::GetElectronCollision:\n"
918 << " Requested electron energy (" << e << " eV) exceeds the "
919 << "current energy range (" << m_eFinalG << " eV).\n"
920 << " Increasing energy range to " << 1.05 * e << " eV.\n";
921 SetMaxElectronEnergy(1.05 * e);
922 } else if (e <= 0.) {
923 std::cerr << m_className << "::GetElectronCollision:\n";
924 std::cerr << " Electron energy must be greater than zero.\n";
925 return false;
926 }
927
928 if (m_isChanged) {
929 if (!UpdateTransportParameters()) {
930 std::cerr << m_className << "::GetElectronCollision:\n";
931 std::cerr << " Error calculating the collision rates table.\n";
932 return false;
933 }
934 m_isChanged = false;
935 }
936
937 // Energy loss
938 double loss = 0.;
939 // Sample the scattering process.
940 if (band >= 0 && band < m_nValleysX) {
941 // X valley
942 // Get the energy interval.
943 int iE = int(e / m_eStepXL);
944 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
945 if (iE < 0) iE = 0;
946 // Select the scattering process.
947 const double r = RndmUniform();
948 int iLow = 0;
949 int iUp = m_nLevelsX - 1;
950 if (r <= m_cfElectronsX[iE][iLow]) {
951 level = iLow;
952 } else if (r >= m_cfElectronsX[iE][m_nLevelsX - 1]) {
953 level = iUp;
954 } else {
955 int iMid;
956 while (iUp - iLow > 1) {
957 iMid = (iLow + iUp) >> 1;
958 if (r < m_cfElectronsX[iE][iMid]) {
959 iUp = iMid;
960 } else {
961 iLow = iMid;
962 }
963 }
964 level = iUp;
965 }
966
967 // Get the collision type.
968 type = m_scatTypeElectronsX[level];
969 // Fill the collision counters.
970 ++m_nCollElectronDetailed[level];
971 ++m_nCollElectronBand[band];
972 if (type == ElectronCollisionTypeAcousticPhonon) {
973 ++m_nCollElectronAcoustic;
974 } else if (type == ElectronCollisionTypeIntervalleyG) {
975 // Intervalley scattering (g type)
976 ++m_nCollElectronIntervalley;
977 // Final valley is in opposite direction.
978 switch (band) {
979 case 0:
980 band = 1;
981 break;
982 case 1:
983 band = 0;
984 break;
985 case 2:
986 band = 3;
987 break;
988 case 3:
989 band = 2;
990 break;
991 case 4:
992 band = 5;
993 break;
994 case 5:
995 band = 4;
996 break;
997 default:
998 break;
999 }
1000 } else if (type == ElectronCollisionTypeIntervalleyF) {
1001 // Intervalley scattering (f type)
1002 ++m_nCollElectronIntervalley;
1003 // Final valley is perpendicular.
1004 switch (band) {
1005 case 0:
1006 case 1:
1007 band = int(RndmUniform() * 4) + 2;
1008 break;
1009 case 2:
1010 case 3:
1011 band = int(RndmUniform() * 4);
1012 if (band > 1) band += 2;
1013 break;
1014 case 4:
1015 case 5:
1016 band = int(RndmUniform() * 4);
1017 break;
1018 }
1019 } else if (type == ElectronCollisionTypeInterbandXL) {
1020 // XL scattering
1021 ++m_nCollElectronIntervalley;
1022 // Final valley is in L band.
1023 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
1024 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysX + m_nValleysL - 1;
1025 } else if (type == ElectronCollisionTypeInterbandXG) {
1026 ++m_nCollElectronIntervalley;
1027 band = m_nValleysX + m_nValleysL;
1028 } else if (type == ElectronCollisionTypeImpurity) {
1029 ++m_nCollElectronImpurity;
1030 } else if (type == ElectronCollisionTypeIonisation) {
1031 ++m_nCollElectronIonisation;
1032 } else {
1033 std::cerr << m_className << "::GetElectronCollision:\n";
1034 std::cerr << " Unexpected collision type (" << type << ").\n";
1035 }
1036
1037 // Get the energy loss.
1038 loss = m_energyLossElectronsX[level];
1039
1040 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
1041 // L valley
1042 // Get the energy interval.
1043 int iE = int(e / m_eStepXL);
1044 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
1045 if (iE < m_ieMinL) iE = m_ieMinL;
1046 // Select the scattering process.
1047 const double r = RndmUniform();
1048 int iLow = 0;
1049 int iUp = m_nLevelsL - 1;
1050 if (r <= m_cfElectronsL[iE][iLow]) {
1051 level = iLow;
1052 } else if (r >= m_cfElectronsL[iE][m_nLevelsL - 1]) {
1053 level = iUp;
1054 } else {
1055 int iMid;
1056 while (iUp - iLow > 1) {
1057 iMid = (iLow + iUp) >> 1;
1058 if (r < m_cfElectronsL[iE][iMid]) {
1059 iUp = iMid;
1060 } else {
1061 iLow = iMid;
1062 }
1063 }
1064 level = iUp;
1065 }
1066
1067 // Get the collision type.
1068 type = m_scatTypeElectronsL[level];
1069 // Fill the collision counters.
1070 ++m_nCollElectronDetailed[m_nLevelsX + level];
1071 ++m_nCollElectronBand[band];
1072 if (type == ElectronCollisionTypeAcousticPhonon) {
1073 ++m_nCollElectronAcoustic;
1074 } else if (type == ElectronCollisionTypeOpticalPhonon) {
1075 ++m_nCollElectronOptical;
1076 } else if (type == ElectronCollisionTypeIntervalleyG ||
1077 type == ElectronCollisionTypeIntervalleyF) {
1078 // Equivalent intervalley scattering
1079 ++m_nCollElectronIntervalley;
1080 // Randomise the final valley.
1081 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
1082 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysX + m_nValleysL;
1083 } else if (type == ElectronCollisionTypeInterbandXL) {
1084 // LX scattering
1085 ++m_nCollElectronIntervalley;
1086 // Randomise the final valley.
1087 band = int(RndmUniform() * m_nValleysX);
1088 if (band >= m_nValleysX) band = m_nValleysX - 1;
1089 } else if (type == ElectronCollisionTypeInterbandLG) {
1090 // LG scattering
1091 ++m_nCollElectronIntervalley;
1092 band = m_nValleysX + m_nValleysL;
1093 } else if (type == ElectronCollisionTypeImpurity) {
1094 ++m_nCollElectronImpurity;
1095 } else if (type == ElectronCollisionTypeIonisation) {
1096 ++m_nCollElectronIonisation;
1097 } else {
1098 std::cerr << m_className << "::GetElectronCollision:\n";
1099 std::cerr << " Unexpected collision type (" << type << ").\n";
1100 }
1101
1102 // Get the energy loss.
1103 loss = m_energyLossElectronsL[level];
1104 } else if (band == m_nValleysX + m_nValleysL) {
1105 // Higher bands
1106 // Get the energy interval.
1107 int iE = int(e / m_eStepG);
1108 if (iE >= nEnergyStepsG) iE = nEnergyStepsG - 1;
1109 if (iE < m_ieMinG) iE = m_ieMinG;
1110 // Select the scattering process.
1111 const double r = RndmUniform();
1112 int iLow = 0;
1113 int iUp = m_nLevelsG - 1;
1114 if (r <= m_cfElectronsG[iE][iLow]) {
1115 level = iLow;
1116 } else if (r >= m_cfElectronsG[iE][m_nLevelsG - 1]) {
1117 level = iUp;
1118 } else {
1119 int iMid;
1120 while (iUp - iLow > 1) {
1121 iMid = (iLow + iUp) >> 1;
1122 if (r < m_cfElectronsG[iE][iMid]) {
1123 iUp = iMid;
1124 } else {
1125 iLow = iMid;
1126 }
1127 }
1128 level = iUp;
1129 }
1130
1131 // Get the collision type.
1132 type = m_scatTypeElectronsG[level];
1133 // Fill the collision counters.
1134 ++m_nCollElectronDetailed[m_nLevelsX + m_nLevelsL + level];
1135 ++m_nCollElectronBand[band];
1136 if (type == ElectronCollisionTypeAcousticPhonon) {
1137 ++m_nCollElectronAcoustic;
1138 } else if (type == ElectronCollisionTypeOpticalPhonon) {
1139 ++m_nCollElectronOptical;
1140 } else if (type == ElectronCollisionTypeIntervalleyG ||
1141 type == ElectronCollisionTypeIntervalleyF) {
1142 // Equivalent intervalley scattering
1143 ++m_nCollElectronIntervalley;
1144 } else if (type == ElectronCollisionTypeInterbandXG) {
1145 // GX scattering
1146 ++m_nCollElectronIntervalley;
1147 // Randomise the final valley.
1148 band = int(RndmUniform() * m_nValleysX);
1149 if (band >= m_nValleysX) band = m_nValleysX - 1;
1150 } else if (type == ElectronCollisionTypeInterbandLG) {
1151 // GL scattering
1152 ++m_nCollElectronIntervalley;
1153 // Randomise the final valley.
1154 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
1155 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysX + m_nValleysL - 1;
1156 } else if (type == ElectronCollisionTypeIonisation) {
1157 ++m_nCollElectronIonisation;
1158 } else {
1159 std::cerr << m_className << "::GetElectronCollision:\n";
1160 std::cerr << " Unexpected collision type (" << type << ").\n";
1161 }
1162
1163 // Get the energy loss.
1164 loss = m_energyLossElectronsG[level];
1165 } else {
1166 std::cerr << m_className << "::GetElectronCollision:\n";
1167 std::cerr << " Band index (" << band << ") out of range.\n";
1168 return false;
1169 }
1170
1171 // Secondaries
1172 nion = ndxc = 0;
1173 // Ionising collision
1174 if (type == ElectronCollisionTypeIonisation) {
1175 double ee = 0., eh = 0.;
1176 ComputeSecondaries(e, ee, eh);
1177 loss = ee + eh + m_bandGap;
1178 m_ionProducts.clear();
1179 // Add the secondary electron.
1180 ionProd newIonProd;
1181 newIonProd.type = IonProdTypeElectron;
1182 newIonProd.energy = ee;
1183 m_ionProducts.push_back(newIonProd);
1184 // Add the hole.
1185 newIonProd.type = IonProdTypeHole;
1186 newIonProd.energy = eh;
1187 m_ionProducts.push_back(newIonProd);
1188 nion = 2;
1189 }
1190
1191 if (e < loss) loss = e - 0.00001;
1192 // Update the energy.
1193 e1 = e - loss;
1194 if (e1 < Small) e1 = Small;
1195
1196 // Update the momentum.
1197 if (band >= 0 && band < m_nValleysX) {
1198 // X valleys
1199 double pstar = sqrt(2. * ElectronMass * e1);
1200 if (m_useNonParabolicity) {
1201 const double alpha = m_alphaX;
1202 pstar *= sqrt(1. + alpha * e1);
1203 }
1204
1205 const double ctheta = 1. - 2. * RndmUniform();
1206 const double stheta = sqrt(1. - ctheta * ctheta);
1207 const double phi = TwoPi * RndmUniform();
1208
1209 if (m_useAnisotropy) {
1210 const double pl = pstar * sqrt(m_mLongX);
1211 const double pt = pstar * sqrt(m_mTransX);
1212 switch (band) {
1213 case 0:
1214 case 1:
1215 // 100
1216 px = pl * ctheta;
1217 py = pt * cos(phi) * stheta;
1218 pz = pt * sin(phi) * stheta;
1219 break;
1220 case 2:
1221 case 3:
1222 // 010
1223 px = pt * sin(phi) * stheta;
1224 py = pl * ctheta;
1225 pz = pt * cos(phi) * stheta;
1226 break;
1227 case 4:
1228 case 5:
1229 // 001
1230 px = pt * cos(phi) * stheta;
1231 py = pt * sin(phi) * stheta;
1232 pz = pl * ctheta;
1233 break;
1234 default:
1235 return false;
1236 break;
1237 }
1238 } else {
1239 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
1240 px = pstar * cos(phi) * stheta;
1241 py = pstar * sin(phi) * stheta;
1242 pz = pstar * ctheta;
1243 }
1244 return true;
1245
1246 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
1247 // L valleys
1248 double pstar = sqrt(2. * ElectronMass * (e1 - m_eMinL));
1249 if (m_useNonParabolicity) {
1250 const double alpha = m_alphaL;
1251 pstar *= sqrt(1. + alpha * (e1 - m_eMinL));
1252 }
1253 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
1254 RndmDirection(px, py, pz, pstar);
1255 return true;
1256 } else {
1257 const double pstar = sqrt(2. * ElectronMass * e1);
1258 RndmDirection(px, py, pz, pstar);
1259 return true;
1260 }
1261
1262 std::cerr << m_className << "::GetElectronCollision:\n";
1263 std::cerr << " Band index (" << band << ") out of range.\n";
1264 e1 = e;
1265 type = 0;
1266 return false;
1267}
void ComputeSecondaries(const double e0, double &ee, double &eh)
bool SetMaxElectronEnergy(const double e)
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:106
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384

◆ GetElectronCollisionRate()

double Garfield::MediumSilicon::GetElectronCollisionRate ( const double  e,
const int  band 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 858 of file MediumSilicon.cc.

858 {
859
860 if (e <= 0.) {
861 std::cerr << m_className << "::GetElectronCollisionRate:\n"
862 << " Electron energy must be positive.\n";
863 return 0.;
864 }
865
866 if (e > m_eFinalG) {
867 std::cerr << m_className << "::GetElectronCollisionRate:\n"
868 << " Collision rate at " << e << " eV (band " << band
869 << ") is not included in the current table.\n"
870 << " Increasing energy range to " << 1.05 * e << " eV.\n";
871 SetMaxElectronEnergy(1.05 * e);
872 }
873
874 if (m_isChanged) {
875 if (!UpdateTransportParameters()) {
876 std::cerr << m_className << "::GetElectronCollisionRate:\n"
877 << " Error calculating the collision rates table.\n";
878 return 0.;
879 }
880 m_isChanged = false;
881 }
882
883 if (band >= 0 && band < m_nValleysX) {
884 int iE = int(e / m_eStepXL);
885 if (iE >= nEnergyStepsXL)
886 iE = nEnergyStepsXL - 1;
887 else if (iE < 0)
888 iE = 0;
889 return m_cfTotElectronsX[iE];
890 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
891 int iE = int(e / m_eStepXL);
892 if (iE >= nEnergyStepsXL)
893 iE = nEnergyStepsXL - 1;
894 else if (iE < m_ieMinL)
895 iE = m_ieMinL;
896 return m_cfTotElectronsL[iE];
897 } else if (band == m_nValleysX + m_nValleysL) {
898 int iE = int(e / m_eStepG);
899 if (iE >= nEnergyStepsG)
900 iE = nEnergyStepsG - 1;
901 else if (iE < m_ieMinG)
902 iE = m_ieMinG;
903 return m_cfTotElectronsG[iE];
904 }
905
906 std::cerr << m_className << "::GetElectronCollisionRate:\n"
907 << " Band index (" << band << ") out of range.\n";
908 return 0.;
909}

◆ GetElectronEnergy()

double Garfield::MediumSilicon::GetElectronEnergy ( const double  px,
const double  py,
const double  pz,
double &  vx,
double &  vy,
double &  vz,
const int  band = 0 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 641 of file MediumSilicon.cc.

643 {
644
645 // Effective masses
646 double mx = ElectronMass, my = ElectronMass, mz = ElectronMass;
647 // Energy offset
648 double e0 = 0.;
649 if (band >= 0 && band < m_nValleysX) {
650 // X valley
651 if (m_useAnisotropy) {
652 switch (band) {
653 case 0:
654 case 1:
655 // X 100, -100
656 mx *= m_mLongX;
657 my *= m_mTransX;
658 mz *= m_mTransX;
659 break;
660 case 2:
661 case 3:
662 // X 010, 0-10
663 mx *= m_mTransX;
664 my *= m_mLongX;
665 mz *= m_mTransX;
666 break;
667 case 4:
668 case 5:
669 // X 001, 00-1
670 mx *= m_mTransX;
671 my *= m_mTransX;
672 mz *= m_mLongX;
673 break;
674 default:
675 std::cerr << m_className << "::GetElectronEnergy:\n"
676 << " Unexpected band index " << band << "!\n";
677 break;
678 }
679 } else {
680 // Conduction effective mass
681 const double mc = 3. / (1. / m_mLongX + 2. / m_mTransX);
682 mx *= mc;
683 my *= mc;
684 mz *= mc;
685 }
686 } else if (band < m_nValleysX + m_nValleysL) {
687 // L valley, isotropic approximation
688 e0 = m_eMinL;
689 // Effective mass
690 const double mc = 3. / (1. / m_mLongL + 2. / m_mTransL);
691 mx *= mc;
692 my *= mc;
693 mz *= mc;
694 } else if (band == m_nValleysX + m_nValleysL) {
695 // Higher band(s)
696 }
697
698 if (m_useNonParabolicity) {
699 // Non-parabolicity parameter
700 double alpha = 0.;
701 if (band < m_nValleysX) {
702 // X valley
703 alpha = m_alphaX;
704 } else if (band < m_nValleysX + m_nValleysL) {
705 // L valley
706 alpha = m_alphaL;
707 }
708
709 const double p2 = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
710 if (alpha > 0.) {
711 const double e = 0.5 * (sqrt(1. + 4 * alpha * p2) - 1.) / alpha;
712 const double a = SpeedOfLight / (1. + 2 * alpha * e);
713 vx = a * px / mx;
714 vy = a * py / my;
715 vz = a * pz / mz;
716 return e0 + e;
717 }
718 }
719
720 const double e = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
721 vx = SpeedOfLight * px / mx;
722 vy = SpeedOfLight * py / my;
723 vz = SpeedOfLight * pz / mz;
724 return e0 + e;
725}

◆ GetElectronMomentum()

void Garfield::MediumSilicon::GetElectronMomentum ( const double  e,
double &  px,
double &  py,
double &  pz,
int &  band 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 727 of file MediumSilicon.cc.

728 {
729
730 int nBands = m_nValleysX;
731 if (e > m_eMinL) nBands += m_nValleysL;
732 if (e > m_eMinG) ++nBands;
733
734 // If the band index is out of range, choose one at random.
735 if (band < 0 || band > m_nValleysX + m_nValleysL ||
736 (e < m_eMinL || band >= m_nValleysX) ||
737 (e < m_eMinG || band == m_nValleysX + m_nValleysL)) {
738 if (e < m_eMinL) {
739 band = int(m_nValleysX * RndmUniform());
740 if (band >= m_nValleysX) band = m_nValleysX - 1;
741 } else {
742 const double dosX = GetConductionBandDensityOfStates(e, 0);
743 const double dosL = GetConductionBandDensityOfStates(e, m_nValleysX);
744 const double dosG =
745 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
746 const double dosSum = m_nValleysX * dosX + m_nValleysL * dosL + dosG;
747 if (dosSum < Small) {
748 band = m_nValleysX + m_nValleysL;
749 } else {
750 const double r = RndmUniform() * dosSum;
751 if (r < dosX) {
752 band = int(m_nValleysX * RndmUniform());
753 if (band >= m_nValleysX) band = m_nValleysX - 1;
754 } else if (r < dosX + dosL) {
755 band = m_nValleysX + int(m_nValleysL * RndmUniform());
756 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysL - 1;
757 } else {
758 band = m_nValleysX + m_nValleysL;
759 }
760 }
761 }
762 if (m_debug) {
763 std::cout << m_className << "::GetElectronMomentum:\n"
764 << " Randomised band index: " << band << "\n";
765 }
766 }
767 if (band < m_nValleysX) {
768 // X valleys
769 double pstar = sqrt(2. * ElectronMass * e);
770 if (m_useNonParabolicity) {
771 const double alpha = m_alphaX;
772 pstar *= sqrt(1. + alpha * e);
773 }
774
775 const double ctheta = 1. - 2. * RndmUniform();
776 const double stheta = sqrt(1. - ctheta * ctheta);
777 const double phi = TwoPi * RndmUniform();
778
779 if (m_useAnisotropy) {
780 const double pl = pstar * sqrt(m_mLongX);
781 const double pt = pstar * sqrt(m_mTransX);
782 switch (band) {
783 case 0:
784 case 1:
785 // 100
786 px = pl * ctheta;
787 py = pt * cos(phi) * stheta;
788 pz = pt * sin(phi) * stheta;
789 break;
790 case 2:
791 case 3:
792 // 010
793 px = pt * sin(phi) * stheta;
794 py = pl * ctheta;
795 pz = pt * cos(phi) * stheta;
796 break;
797 case 4:
798 case 5:
799 // 001
800 px = pt * cos(phi) * stheta;
801 py = pt * sin(phi) * stheta;
802 pz = pl * ctheta;
803 break;
804 default:
805 // Other band; should not occur.
806 std::cerr << m_className << "::GetElectronMomentum:\n"
807 << " Unexpected band index (" << band << ").\n";
808 px = pstar * stheta * cos(phi);
809 py = pstar * stheta * sin(phi);
810 pz = pstar * ctheta;
811 break;
812 }
813 } else {
814 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
815 px = pstar * cos(phi) * stheta;
816 py = pstar * sin(phi) * stheta;
817 pz = pstar * ctheta;
818 }
819 } else if (band < m_nValleysX + m_nValleysL) {
820 // L valleys
821 double pstar = sqrt(2. * ElectronMass * (e - m_eMinL));
822 if (m_useNonParabolicity) {
823 const double alpha = m_alphaL;
824 pstar *= sqrt(1. + alpha * (e - m_eMinL));
825 }
826 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
827 RndmDirection(px, py, pz, pstar);
828 } else if (band == m_nValleysX + m_nValleysL) {
829 // Higher band
830 const double pstar = sqrt(2. * ElectronMass * e);
831 RndmDirection(px, py, pz, pstar);
832 }
833}

◆ GetElectronNullCollisionRate()

double Garfield::MediumSilicon::GetElectronNullCollisionRate ( const int  band)
virtual

Reimplemented from Garfield::Medium.

Definition at line 835 of file MediumSilicon.cc.

835 {
836
837 if (m_isChanged) {
838 if (!UpdateTransportParameters()) {
839 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
840 << " Error calculating the collision rates table.\n";
841 return 0.;
842 }
843 m_isChanged = false;
844 }
845
846 if (band >= 0 && band < m_nValleysX) {
847 return m_cfNullElectronsX;
848 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
849 return m_cfNullElectronsL;
850 } else if (band == m_nValleysX + m_nValleysL) {
851 return m_cfNullElectronsG;
852 }
853 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
854 << " Band index (" << band << ") out of range.\n";
855 return 0.;
856}

◆ GetIonisationProduct()

bool Garfield::MediumSilicon::GetIonisationProduct ( const unsigned int  i,
int &  type,
double &  energy 
) const
virtual

Reimplemented from Garfield::Medium.

Definition at line 1269 of file MediumSilicon.cc.

1270 {
1271
1272 if (i >= m_ionProducts.size()) {
1273 std::cerr << m_className << "::GetIonisationProduct:\n"
1274 << " Index (" << i << ") out of range.\n";
1275 return false;
1276 }
1277
1278 type = m_ionProducts[i].type;
1279 energy = m_ionProducts[i].energy;
1280 return true;
1281}

◆ GetMaxElectronEnergy()

double Garfield::MediumSilicon::GetMaxElectronEnergy ( ) const
inline

Definition at line 78 of file MediumSilicon.hh.

78{ return m_eFinalG; }

◆ GetNumberOfElectronBands()

unsigned int Garfield::MediumSilicon::GetNumberOfElectronBands ( ) const

Definition at line 1321 of file MediumSilicon.cc.

1321 {
1322
1323 return m_nValleysX + m_nValleysL + 1;
1324}

◆ GetNumberOfElectronCollisions() [1/2]

unsigned int Garfield::MediumSilicon::GetNumberOfElectronCollisions ( ) const

Definition at line 1297 of file MediumSilicon.cc.

1297 {
1298
1299 return m_nCollElectronAcoustic + m_nCollElectronOptical +
1300 m_nCollElectronIntervalley + m_nCollElectronImpurity +
1301 m_nCollElectronIonisation;
1302}

◆ GetNumberOfElectronCollisions() [2/2]

unsigned int Garfield::MediumSilicon::GetNumberOfElectronCollisions ( const unsigned int  level) const

Definition at line 1310 of file MediumSilicon.cc.

1310 {
1311
1312 const unsigned int nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1313 if (level >= nLevels) {
1314 std::cerr << m_className << "::GetNumberOfElectronCollisions:\n"
1315 << " Scattering rate term (" << level << ") does not exist.\n";
1316 return 0;
1317 }
1318 return m_nCollElectronDetailed[level];
1319}

◆ GetNumberOfIonisationProducts()

unsigned int Garfield::MediumSilicon::GetNumberOfIonisationProducts ( ) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 111 of file MediumSilicon.hh.

111 {
112 return m_ionProducts.size();
113 }

◆ GetNumberOfLevels()

unsigned int Garfield::MediumSilicon::GetNumberOfLevels ( ) const

Definition at line 1304 of file MediumSilicon.cc.

1304 {
1305
1306 return m_nLevelsX + m_nLevelsL + m_nLevelsG;
1307}

◆ GetOpticalDataRange()

bool Garfield::MediumSilicon::GetOpticalDataRange ( double &  emin,
double &  emax,
const unsigned int  i = 0 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 1337 of file MediumSilicon.cc.

1338 {
1339
1340 if (i != 0) {
1341 std::cerr << m_className << "::GetOpticalDataRange:\n";
1342 std::cerr << " Medium has only one component.\n";
1343 }
1344
1345 // Make sure the optical data table has been loaded.
1346 if (m_opticalDataTable.empty()) {
1347 if (!LoadOpticalData(m_opticalDataFile)) {
1348 std::cerr << m_className << "::GetOpticalDataRange:\n";
1349 std::cerr << " Optical data table could not be loaded.\n";
1350 return false;
1351 }
1352 }
1353
1354 emin = m_opticalDataTable[0].energy;
1355 emax = m_opticalDataTable.back().energy;
1356 if (m_debug) {
1357 std::cout << m_className << "::GetOpticalDataRange:\n"
1358 << " " << emin << " < E [eV] < " << emax << "\n";
1359 }
1360 return true;
1361}

◆ GetValenceBandDensityOfStates()

double Garfield::MediumSilicon::GetValenceBandDensityOfStates ( const double  e,
const int  band = -1 
)

Definition at line 3158 of file MediumSilicon.cc.

3159 {
3160
3161 if (band <= 0) {
3162 // Total (full-band) density of states.
3163 const int nPoints = m_fbDosValence.size();
3164 int iE = int(e / m_eStepDos);
3165 if (iE >= nPoints || iE < 0) {
3166 return 0.;
3167 } else if (iE == nPoints - 1) {
3168 return m_fbDosValence[nPoints - 1];
3169 }
3170
3171 const double dos =
3172 m_fbDosValence[iE] +
3173 (m_fbDosValence[iE + 1] - m_fbDosValence[iE]) * (e / m_eStepDos - iE);
3174 return dos * 1.e21;
3175 }
3176
3177 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
3178 << " Band index (" << band << ") out of range.\n";
3179 return 0.;
3180}

◆ HoleAttachment()

bool Garfield::MediumSilicon::HoleAttachment ( const double  ex,
const double  ey,
const double  ez,
const double  bx,
const double  by,
const double  bz,
double &  eta 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 462 of file MediumSilicon.cc.

465 {
466
467 eta = 0.;
468 if (m_isChanged) {
469 if (!UpdateTransportParameters()) {
470 std::cerr << m_className << "::HoleAttachment:\n"
471 << " Error calculating the transport parameters.\n";
472 return false;
473 }
474 m_isChanged = false;
475 }
476
478 // Interpolation in user table.
479 return Medium::HoleAttachment(ex, ey, ez, bx, by, bz, eta);
480 }
481
482 switch (m_trappingModel) {
483 case 0:
484 eta = m_hTrapCs * m_hTrapDensity;
485 break;
486 case 1:
487 double vx, vy, vz;
488 HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
489 eta = m_hTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
490 if (eta > 0.) eta = 1. / eta;
491 break;
492 default:
493 std::cerr << m_className << "::HoleAttachment:\n"
494 << " Unknown model. Program bug!\n";
495 return false;
496 break;
497 }
498
499 return true;
500}
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)
bool m_hasHoleAttachment
Definition: Medium.hh:358
virtual bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Definition: Medium.cc:1023

◆ HoleTownsend()

bool Garfield::MediumSilicon::HoleTownsend ( const double  ex,
const double  ey,
const double  ez,
const double  bx,
const double  by,
const double  bz,
double &  alpha 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 425 of file MediumSilicon.cc.

428 {
429
430 alpha = 0.;
431 if (m_isChanged) {
432 if (!UpdateTransportParameters()) {
433 std::cerr << m_className << "::HoleTownsend:\n"
434 << " Error calculating the transport parameters.\n";
435 return false;
436 }
437 m_isChanged = false;
438 }
439
440 if (m_hasHoleTownsend) {
441 // Interpolation in user table.
442 return Medium::HoleTownsend(ex, ey, ez, bx, by, bz, alpha);
443 }
444
445 const double e = sqrt(ex * ex + ey * ey + ez * ez);
446
447 switch (m_impactIonisationModel) {
448 case ImpactIonisationModelVanOverstraeten:
449 return HoleImpactIonisationVanOverstraetenDeMan(e, alpha);
450 break;
451 case ImpactIonisationModelGrant:
452 return HoleImpactIonisationGrant(e, alpha);
453 break;
454 default:
455 std::cerr << m_className << "::HoleTownsend:\n"
456 << " Unknown model. Program bug!\n";
457 break;
458 }
459 return false;
460}
virtual bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Definition: Medium.cc:971
bool m_hasHoleTownsend
Definition: Medium.hh:358

◆ HoleVelocity()

bool Garfield::MediumSilicon::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 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 369 of file MediumSilicon.cc.

372 {
373
374 vx = vy = vz = 0.;
375 if (m_isChanged) {
376 if (!UpdateTransportParameters()) {
377 std::cerr << m_className << "::HoleVelocity:\n"
378 << " Error calculating the transport parameters.\n";
379 return false;
380 }
381 m_isChanged = false;
382 }
383
384 if (m_hasHoleVelocityE) {
385 // Interpolation in user table.
386 return Medium::HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
387 }
388
389 const double e = sqrt(ex * ex + ey * ey + ez * ez);
390
391 // Calculate the mobility
392 double mu;
393 switch (m_highFieldMobilityModel) {
394 case HighFieldMobilityModelMinimos:
395 HoleMobilityMinimos(e, mu);
396 break;
397 case HighFieldMobilityModelCanali:
398 HoleMobilityCanali(e, mu);
399 break;
400 case HighFieldMobilityModelReggiani:
401 HoleMobilityReggiani(e, mu);
402 break;
403 default:
404 mu = m_hMobility;
405 }
406
407 const double b = sqrt(bx * bx + by * by + bz * bz);
408 if (b < Small) {
409 vx = mu * ex;
410 vy = mu * ey;
411 vz = mu * ez;
412 } else {
413 // Hall mobility
414 const double muH = m_hHallFactor * mu;
415 const double eb = bx * ex + by * ey + bz * ez;
416 const double nom = 1. + pow(muH * b, 2);
417 // Compute the drift velocity using the Langevin equation.
418 vx = mu * (ex + muH * (ey * bz - ez * by) + muH * muH * bx * eb) / nom;
419 vy = mu * (ey + muH * (ez * bx - ex * bz) + muH * muH * by * eb) / nom;
420 vz = mu * (ez + muH * (ex * by - ey * bx) + muH * muH * bz * eb) / nom;
421 }
422 return true;
423}
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)
Definition: Medium.cc:704
bool m_hasHoleVelocityE
Definition: Medium.hh:356

Referenced by HoleAttachment().

◆ Initialise()

bool Garfield::MediumSilicon::Initialise ( )

Definition at line 1434 of file MediumSilicon.cc.

1434 {
1435
1436 if (!m_isChanged) {
1437 if (m_debug) {
1438 std::cerr << m_className << "::Initialise:\n Nothing changed.\n";
1439 }
1440 return true;
1441 }
1442 if (!UpdateTransportParameters()) {
1443 std::cerr << m_className << "::Initialise: Error preparing "
1444 << "transport parameters/calculating collision rates.\n";
1445 return false;
1446 }
1447 return true;
1448}

◆ IsSemiconductor()

bool Garfield::MediumSilicon::IsSemiconductor ( ) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 17 of file MediumSilicon.hh.

17{ return true; }

◆ ResetCollisionCounters()

void Garfield::MediumSilicon::ResetCollisionCounters ( )

Definition at line 1283 of file MediumSilicon.cc.

1283 {
1284
1285 m_nCollElectronAcoustic = m_nCollElectronOptical = 0;
1286 m_nCollElectronIntervalley = 0;
1287 m_nCollElectronImpurity = 0;
1288 m_nCollElectronIonisation = 0;
1289 const int nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1290 m_nCollElectronDetailed.resize(nLevels);
1291 for (int j = nLevels; j--;) m_nCollElectronDetailed[j] = 0;
1292 const int nBands = m_nValleysX + m_nValleysL + 1;
1293 m_nCollElectronBand.resize(nBands);
1294 for (int j = nBands; j--;) m_nCollElectronBand[j] = 0;
1295}

◆ SetDiffusionScaling()

void Garfield::MediumSilicon::SetDiffusionScaling ( const double  d)
inline

Definition at line 72 of file MediumSilicon.hh.

72 {
73 diffScale = d;
74 }

◆ SetDoping()

void Garfield::MediumSilicon::SetDoping ( const char  type,
const double  c 
)

Doping concentration [cm-3] and type ('i', 'n', 'p')

Definition at line 137 of file MediumSilicon.cc.

137 {
138
139 if (toupper(type) == 'N') {
140 m_dopingType = 'n';
141 if (c > Small) {
142 m_dopingConcentration = c;
143 } else {
144 std::cerr << m_className << "::SetDoping:\n"
145 << " Doping concentration must be greater than zero.\n"
146 << " Using default value for n-type silicon "
147 << "(10^12 cm-3) instead.\n";
148 m_dopingConcentration = 1.e12;
149 }
150 } else if (toupper(type) == 'P') {
151 m_dopingType = 'p';
152 if (c > Small) {
153 m_dopingConcentration = c;
154 } else {
155 std::cerr << m_className << "::SetDoping:\n"
156 << " Doping concentration must be greater than zero.\n"
157 << " Using default value for p-type silicon "
158 << "(10^18 cm-3) instead.\n";
159 m_dopingConcentration = 1.e18;
160 }
161 } else if (toupper(type) == 'I') {
162 m_dopingType = 'i';
163 m_dopingConcentration = 0.;
164 } else {
165 std::cerr << m_className << "::SetDoping:\n"
166 << " Unknown dopant type (" << type << ").\n"
167 << " Available types are n, p and i (intrinsic).\n";
168 return;
169 }
170
171 m_isChanged = true;
172}

◆ SetDopingMobilityModelMasetti()

void Garfield::MediumSilicon::SetDopingMobilityModelMasetti ( )

Definition at line 544 of file MediumSilicon.cc.

544 {
545
546 m_dopingMobilityModel = DopingMobilityModelMasetti;
547 m_hasUserMobility = false;
548 m_isChanged = true;
549}

◆ SetDopingMobilityModelMinimos()

void Garfield::MediumSilicon::SetDopingMobilityModelMinimos ( )

Definition at line 537 of file MediumSilicon.cc.

537 {
538
539 m_dopingMobilityModel = DopingMobilityModelMinimos;
540 m_hasUserMobility = false;
541 m_isChanged = true;
542}

◆ SetHighFieldMobilityModelCanali()

void Garfield::MediumSilicon::SetHighFieldMobilityModelCanali ( )

Definition at line 594 of file MediumSilicon.cc.

594 {
595
596 m_highFieldMobilityModel = HighFieldMobilityModelCanali;
597 m_isChanged = true;
598}

◆ SetHighFieldMobilityModelConstant()

void Garfield::MediumSilicon::SetHighFieldMobilityModelConstant ( )

Definition at line 606 of file MediumSilicon.cc.

606 {
607
608 m_highFieldMobilityModel = HighFieldMobilityModelConstant;
609}

◆ SetHighFieldMobilityModelMinimos()

void Garfield::MediumSilicon::SetHighFieldMobilityModelMinimos ( )

Definition at line 588 of file MediumSilicon.cc.

588 {
589
590 m_highFieldMobilityModel = HighFieldMobilityModelMinimos;
591 m_isChanged = true;
592}

◆ SetHighFieldMobilityModelReggiani()

void Garfield::MediumSilicon::SetHighFieldMobilityModelReggiani ( )

Definition at line 600 of file MediumSilicon.cc.

600 {
601
602 m_highFieldMobilityModel = HighFieldMobilityModelReggiani;
603 m_isChanged = true;
604}

◆ SetImpactIonisationModelGrant()

void Garfield::MediumSilicon::SetImpactIonisationModelGrant ( )

Definition at line 617 of file MediumSilicon.cc.

617 {
618
619 m_impactIonisationModel = ImpactIonisationModelGrant;
620 m_isChanged = true;
621}

◆ SetImpactIonisationModelVanOverstraetenDeMan()

void Garfield::MediumSilicon::SetImpactIonisationModelVanOverstraetenDeMan ( )

Definition at line 611 of file MediumSilicon.cc.

611 {
612
613 m_impactIonisationModel = ImpactIonisationModelVanOverstraeten;
614 m_isChanged = true;
615}

◆ SetLatticeMobilityModelMinimos()

void Garfield::MediumSilicon::SetLatticeMobilityModelMinimos ( )

Definition at line 516 of file MediumSilicon.cc.

516 {
517
518 m_latticeMobilityModel = LatticeMobilityModelMinimos;
519 m_hasUserMobility = false;
520 m_isChanged = true;
521}

◆ SetLatticeMobilityModelReggiani()

void Garfield::MediumSilicon::SetLatticeMobilityModelReggiani ( )

Definition at line 530 of file MediumSilicon.cc.

530 {
531
532 m_latticeMobilityModel = LatticeMobilityModelReggiani;
533 m_hasUserMobility = false;
534 m_isChanged = true;
535}

◆ SetLatticeMobilityModelSentaurus()

void Garfield::MediumSilicon::SetLatticeMobilityModelSentaurus ( )

Definition at line 523 of file MediumSilicon.cc.

523 {
524
525 m_latticeMobilityModel = LatticeMobilityModelSentaurus;
526 m_hasUserMobility = false;
527 m_isChanged = true;
528}

◆ SetLowFieldMobility()

void Garfield::MediumSilicon::SetLowFieldMobility ( const double  mue,
const double  muh 
)

Definition at line 502 of file MediumSilicon.cc.

502 {
503
504 if (mue <= 0. || muh <= 0.) {
505 std::cerr << m_className << "::SetLowFieldMobility:\n"
506 << " Mobility must be greater than zero.\n";
507 return;
508 }
509
510 m_eMobility = mue;
511 m_hMobility = muh;
512 m_hasUserMobility = true;
513 m_isChanged = true;
514}

◆ SetMaxElectronEnergy()

bool Garfield::MediumSilicon::SetMaxElectronEnergy ( const double  e)

Definition at line 623 of file MediumSilicon.cc.

623 {
624
625 if (e <= m_eMinG + Small) {
626 std::cerr << m_className << "::SetMaxElectronEnergy:\n"
627 << " Requested upper electron energy limit (" << e
628 << " eV) is too small.\n";
629 return false;
630 }
631
632 m_eFinalG = e;
633 // Determine the energy interval size.
634 m_eStepG = m_eFinalG / nEnergyStepsG;
635
636 m_isChanged = true;
637
638 return true;
639}

Referenced by GetElectronCollision(), and GetElectronCollisionRate().

◆ SetSaturationVelocity()

void Garfield::MediumSilicon::SetSaturationVelocity ( const double  vsate,
const double  vsath 
)

Definition at line 551 of file MediumSilicon.cc.

552 {
553
554 if (vsate <= 0. || vsath <= 0.) {
555 std::cout << m_className << "::SetSaturationVelocity:\n"
556 << " Restoring default values.\n";
557 m_hasUserSaturationVelocity = false;
558 } else {
559 m_eSatVel = vsate;
560 m_hSatVel = vsath;
561 m_hasUserSaturationVelocity = true;
562 }
563
564 m_isChanged = true;
565}

◆ SetSaturationVelocityModelCanali()

void Garfield::MediumSilicon::SetSaturationVelocityModelCanali ( )

Definition at line 574 of file MediumSilicon.cc.

574 {
575
576 m_saturationVelocityModel = SaturationVelocityModelCanali;
577 m_hasUserSaturationVelocity = false;
578 m_isChanged = true;
579}

◆ SetSaturationVelocityModelMinimos()

void Garfield::MediumSilicon::SetSaturationVelocityModelMinimos ( )

Definition at line 567 of file MediumSilicon.cc.

567 {
568
569 m_saturationVelocityModel = SaturationVelocityModelMinimos;
570 m_hasUserSaturationVelocity = false;
571 m_isChanged = true;
572}

◆ SetSaturationVelocityModelReggiani()

void Garfield::MediumSilicon::SetSaturationVelocityModelReggiani ( )

Definition at line 581 of file MediumSilicon.cc.

581 {
582
583 m_saturationVelocityModel = SaturationVelocityModelReggiani;
584 m_hasUserSaturationVelocity = false;
585 m_isChanged = true;
586}

◆ SetTrapCrossSection()

void Garfield::MediumSilicon::SetTrapCrossSection ( const double  ecs,
const double  hcs 
)

Trapping cross-sections for electrons and holes.

Definition at line 180 of file MediumSilicon.cc.

180 {
181
182 if (ecs < 0.) {
183 std::cerr << m_className << "::SetTrapCrossSection:\n"
184 << " Capture cross-section [cm2] must non-negative.\n";
185 } else {
186 m_eTrapCs = ecs;
187 }
188
189 if (hcs < 0.) {
190 std::cerr << m_className << "::SetTrapCrossSection:\n"
191 << " Capture cross-section [cm2] must be non-negative.n";
192 } else {
193 m_hTrapCs = hcs;
194 }
195
196 m_trappingModel = 0;
197 m_isChanged = true;
198}

◆ SetTrapDensity()

void Garfield::MediumSilicon::SetTrapDensity ( const double  n)

Definition at line 200 of file MediumSilicon.cc.

200 {
201
202 if (n < 0.) {
203 std::cerr << m_className << "::SetTrapDensity:\n"
204 << " Trap density [cm-3] must be non-negative.\n";
205 } else {
206 m_eTrapDensity = n;
207 m_hTrapDensity = n;
208 }
209
210 m_trappingModel = 0;
211 m_isChanged = true;
212}

◆ SetTrappingTime()

void Garfield::MediumSilicon::SetTrappingTime ( const double  etau,
const double  htau 
)

Definition at line 214 of file MediumSilicon.cc.

214 {
215
216 if (etau <= 0.) {
217 std::cerr << m_className << "::SetTrappingTime:\n"
218 << " Trapping time [ns-1] must be positive.\n";
219 } else {
220 m_eTrapTime = etau;
221 }
222
223 if (htau <= 0.) {
224 std::cerr << m_className << "::SetTrappingTime:\n"
225 << " Trapping time [ns-1] must be positive.\n";
226 } else {
227 m_hTrapTime = htau;
228 }
229
230 m_trappingModel = 1;
231 m_isChanged = true;
232}

The documentation for this class was generated from the following files: