Garfield++ v1r0
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

#include <MediumSilicon.hh>

+ Inheritance diagram for Garfield::MediumSilicon:

Public Member Functions

 MediumSilicon ()
 
 ~MediumSilicon ()
 
bool IsSemiconductor () const
 
void SetDoping (const char &type, const double &c)
 
void GetDoping (char &type, double &c) const
 
void SetTrapCrossSection (const double &ecs, const double &hcs)
 
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 ()
 
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)
 
int GetNumberOfIonisationProducts ()
 
bool GetIonisationProduct (const int i, int &type, double &energy)
 
double GetConductionBandDensityOfStates (const double e, const int band=0)
 
double GetValenceBandDensityOfStates (const double e, const int band=-1)
 
void ResetCollisionCounters ()
 
int GetNumberOfElectronCollisions () const
 
int GetNumberOfLevels ()
 
int GetNumberOfElectronCollisions (const int level) const
 
int GetNumberOfElectronBands ()
 
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
 
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 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 int GetNumberOfIonisationProducts ()
 
virtual bool GetIonisationProduct (const int i, int &type, double &energy)
 
virtual int GetNumberOfDeexcitationProducts ()
 
virtual bool GetDeexcitationProduct (const int i, double &t, double &s, int &type, double &energy)
 
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 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 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 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 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
 
unsigned int m_nEfields
 
unsigned int m_nBfields
 
unsigned int m_nAngles
 
std::vector< double > eFields
 
std::vector< double > bFields
 
std::vector< double > bAngles
 
bool m_map2d
 
bool m_hasElectronVelocityE
 
bool m_hasElectronVelocityB
 
bool m_hasElectronVelocityExB
 
bool m_hasElectronDiffLong
 
bool m_hasElectronDiffTrans
 
bool m_hasElectronDiffTens
 
bool m_hasElectronTownsend
 
bool m_hasElectronAttachment
 
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< 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_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_intpMobility
 
unsigned int m_intpDissociation
 
- Static Protected Attributes inherited from Garfield::Medium
static int m_idCounter = -1
 

Detailed Description

Definition at line 10 of file MediumSilicon.hh.

Constructor & Destructor Documentation

◆ MediumSilicon()

Garfield::MediumSilicon::MediumSilicon ( )

Definition at line 15 of file MediumSilicon.cc.

16 : Medium(),
17 m_bandGap(1.12),
18 m_dopingType('i'),
19 m_dopingConcentration(0.),
20 mLongX(0.916),
21 mTransX(0.191),
22 mLongL(1.59),
23 mTransL(0.12),
24 alphaX(0.5),
25 alphaL(0.5),
26 eLatticeMobility(1.35e-6),
27 hLatticeMobility(0.45e-6),
28 eMobility(1.35e-6),
29 hMobility(0.45e-6),
30 eBetaCanali(1.109),
31 hBetaCanali(1.213),
32 eBetaCanaliInv(1. / 1.109),
33 hBetaCanaliInv(1. / 1.213),
34 eSatVel(1.02e-2),
35 hSatVel(0.72e-2),
36 eHallFactor(1.15),
37 hHallFactor(0.7),
38 eTrapCs(1.e-15),
39 hTrapCs(1.e-15),
40 eTrapDensity(1.e13),
41 hTrapDensity(1.e13),
42 eTrapTime(0.),
43 hTrapTime(0.),
44 trappingModel(0),
45 eImpactA0(3.318e5),
46 eImpactA1(0.703e6),
47 eImpactA2(0.),
48 eImpactB0(1.135e6),
49 eImpactB1(1.231e6),
50 eImpactB2(0.),
51 hImpactA0(1.582e6),
52 hImpactA1(0.671e6),
53 hImpactB0(2.036e6),
54 hImpactB1(1.693e6),
55 m_hasUserMobility(false),
56 m_hasUserSaturationVelocity(false),
57 latticeMobilityModel(LatticeMobilityModelSentaurus),
58 dopingMobilityModel(DopingMobilityModelMasetti),
59 saturationVelocityModel(SaturationVelocityModelCanali),
60 highFieldMobilityModel(HighFieldMobilityModelCanali),
61 impactIonisationModel(ImpactIonisationModelVanOverstraeten),
62 useCfOutput(false),
63 useNonParabolicity(true),
64 useFullBandDos(true),
65 useAnisotropy(true),
66 eFinalXL(4.),
67 eStepXL(eFinalXL / nEnergyStepsXL),
68 eFinalG(10.),
69 eStepG(eFinalG / nEnergyStepsG),
70 eFinalV(8.5),
71 eStepV(eFinalV / nEnergyStepsV),
72 nLevelsX(0),
73 nLevelsL(0),
74 nLevelsG(0),
75 nLevelsV(0),
76 nValleysX(6),
77 nValleysL(8),
78 eMinL(1.05),
79 eMinG(2.24),
80 ieMinL(0),
81 ieMinG(0),
82 m_hasOpticalData(false),
83 opticalDataFile("OpticalData_Si.txt") {
84
85 m_className = "MediumSilicon";
86 m_name = "Si";
87
88 SetTemperature(300.);
90 SetAtomicNumber(14.);
91 SetAtomicWeight(28.0855);
92 SetMassDensity(2.329);
93
96 m_microscopic = true;
97
98 m_w = 3.6;
99 m_fano = 0.11;
100
101 cfTotElectronsX.clear();
102 cfElectronsX.clear();
103 energyLossElectronsX.clear();
104 scatTypeElectronsX.clear();
105
106 cfTotElectronsL.clear();
107 cfElectronsL.clear();
108 energyLossElectronsL.clear();
109 scatTypeElectronsL.clear();
110
111 cfTotElectronsG.clear();
112 cfElectronsG.clear();
113 energyLossElectronsG.clear();
114 scatTypeElectronsG.clear();
115
116 ieMinL = int(eMinL / eStepXL) + 1;
117 ieMinG = int(eMinG / eStepG) + 1;
118
119 cfTotHoles.clear();
120 cfHoles.clear();
121 energyLossHoles.clear();
122 scatTypeHoles.clear();
123
124 // Load the density of states table.
125 InitialiseDensityOfStates();
126
127 // Initialize the collision counters.
128 nCollElectronAcoustic = nCollElectronOptical = 0;
129 nCollElectronIntervalley = 0;
130 nCollElectronImpurity = 0;
131 nCollElectronIonisation = 0;
132 nCollElectronDetailed.clear();
133 nCollElectronBand.clear();
134
135 nIonisationProducts = 0;
136 ionProducts.clear();
137}
bool m_microscopic
Definition: Medium.hh:309
virtual void SetAtomicWeight(const double &a)
Definition: Medium.cc:178
virtual void SetMassDensity(const double &rho)
Definition: Medium.cc:200
double m_fano
Definition: Medium.hh:313
std::string m_name
Definition: Medium.hh:291
virtual void SetAtomicNumber(const double &z)
Definition: Medium.cc:167
virtual void EnableDrift()
Definition: Medium.hh:52
void SetTemperature(const double &t)
Definition: Medium.cc:117
void SetDielectricConstant(const double &eps)
Definition: Medium.cc:139
virtual void EnablePrimaryIonisation()
Definition: Medium.hh:54
std::string m_className
Definition: Medium.hh:284

◆ ~MediumSilicon()

Garfield::MediumSilicon::~MediumSilicon ( )
inline

Definition at line 16 of file MediumSilicon.hh.

16{}

Member Function Documentation

◆ ComputeSecondaries()

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

Definition at line 3164 of file MediumSilicon.cc.

3165 {
3166
3167 const double widthValenceBand = eStepDos * nFbDosEntriesValence;
3168 const double widthConductionBand = eStepDos * nFbDosEntriesConduction;
3169
3170 bool ok = false;
3171 while (!ok) {
3172 // Sample a hole energy according to the valence band DOS.
3173 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
3174 int ih = std::min(int(eh / eStepDos), nFbDosEntriesValence - 1);
3175 while (RndmUniform() > fbDosValence[ih] / fbDosMaxV) {
3176 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
3177 ih = std::min(int(eh / eStepDos), nFbDosEntriesValence - 1);
3178 }
3179 // Sample an electron energy according to the conduction band DOS.
3180 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
3181 int ie = std::min(int(ee / eStepDos), nFbDosEntriesConduction - 1);
3182 while (RndmUniform() > fbDosConduction[ie] / fbDosMaxC) {
3183 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
3184 ie = std::min(int(ee / eStepDos), nFbDosEntriesConduction - 1);
3185 }
3186 // Calculate the energy of the primary electron.
3187 const double ep = e0 - m_bandGap - eh - ee;
3188 if (ep < Small) continue;
3189 if (ep > 5.) return;
3190 // Check if the primary electron energy is consistent with the DOS.
3191 int ip = std::min(int(ep / eStepDos), nFbDosEntriesConduction - 1);
3192 if (RndmUniform() > fbDosConduction[ip] / fbDosMaxC) continue;
3193 ok = true;
3194 }
3195}
double RndmUniform()
Definition: Random.hh:16
double RndmUniformPos()
Definition: Random.hh:19

Referenced by GetElectronCollision().

◆ DisableAnisotropy()

void Garfield::MediumSilicon::DisableAnisotropy ( )
inline

Definition at line 87 of file MediumSilicon.hh.

87{ useAnisotropy = false; }

◆ DisableFullBandDensityOfStates()

void Garfield::MediumSilicon::DisableFullBandDensityOfStates ( )
inline

Definition at line 85 of file MediumSilicon.hh.

85{ useFullBandDos = false; }

◆ DisableNonParabolicity()

void Garfield::MediumSilicon::DisableNonParabolicity ( )
inline

Definition at line 83 of file MediumSilicon.hh.

83{ useNonParabolicity = false; }

◆ DisableScatteringRateOutput()

void Garfield::MediumSilicon::DisableScatteringRateOutput ( )
inline

Definition at line 80 of file MediumSilicon.hh.

80{ 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 331 of file MediumSilicon.cc.

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

◆ 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 294 of file MediumSilicon.cc.

297 {
298
299 alpha = 0.;
300 if (m_isChanged) {
301 if (!UpdateTransportParameters()) {
302 std::cerr << m_className << "::ElectronTownsend:\n";
303 std::cerr << " Error calculating the transport parameters.\n";
304 return false;
305 }
306 m_isChanged = false;
307 }
308
310 // Interpolation in user table.
311 return Medium::ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
312 }
313
314 const double e = sqrt(ex * ex + ey * ey + ez * ez);
315
316 switch (impactIonisationModel) {
317 case ImpactIonisationModelVanOverstraeten:
318 return ElectronImpactIonisationVanOverstraetenDeMan(e, alpha);
319 break;
320 case ImpactIonisationModelGrant:
321 return ElectronImpactIonisationGrant(e, alpha);
322 break;
323 default:
324 std::cerr << m_className << "::ElectronTownsend:\n";
325 std::cerr << " Unknown model activated. Program bug!\n";
326 break;
327 }
328 return false;
329}
bool m_hasElectronTownsend
Definition: Medium.hh:335
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:542

◆ 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 236 of file MediumSilicon.cc.

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

Referenced by ElectronAttachment().

◆ EnableAnisotropy()

void Garfield::MediumSilicon::EnableAnisotropy ( )
inline

Definition at line 86 of file MediumSilicon.hh.

86{ useAnisotropy = true; }

◆ EnableFullBandDensityOfStates()

void Garfield::MediumSilicon::EnableFullBandDensityOfStates ( )
inline

Definition at line 84 of file MediumSilicon.hh.

84{ useFullBandDos = true; }

◆ EnableNonParabolicity()

void Garfield::MediumSilicon::EnableNonParabolicity ( )
inline

Definition at line 82 of file MediumSilicon.hh.

82{ useNonParabolicity = true; }

◆ EnableScatteringRateOutput()

void Garfield::MediumSilicon::EnableScatteringRateOutput ( )
inline

Definition at line 79 of file MediumSilicon.hh.

79{ useCfOutput = true; }

◆ GetConductionBandDensityOfStates()

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

Definition at line 3031 of file MediumSilicon.cc.

3032 {
3033 if (band < 0) {
3034 int iE = int(e / eStepDos);
3035 if (iE >= nFbDosEntriesConduction || iE < 0) {
3036 return 0.;
3037 } else if (iE == nFbDosEntriesConduction - 1) {
3038 return fbDosConduction[nFbDosEntriesConduction - 1];
3039 }
3040
3041 const double dos =
3042 fbDosConduction[iE] +
3043 (fbDosConduction[iE + 1] - fbDosConduction[iE]) * (e / eStepDos - iE);
3044 return dos * 1.e21;
3045
3046 } else if (band < nValleysX) {
3047 // X valleys
3048 if (e <= 0.) return 0.;
3049 // Density-of-states effective mass (cube)
3050 const double md3 = pow(ElectronMass, 3) * mLongX * mTransX * mTransX;
3051
3052 if (useFullBandDos) {
3053 if (e < eMinL) {
3054 return GetConductionBandDensityOfStates(e, -1) / nValleysX;
3055 } else if (e < eMinG) {
3056 // Subtract the fraction of the full-band density of states
3057 // attributed to the L valleys.
3058 const double dosX =
3060 GetConductionBandDensityOfStates(e, nValleysX) * nValleysL;
3061 return dosX / nValleysX;
3062 } else {
3063 // Subtract the fraction of the full-band density of states
3064 // attributed to the L valleys and the higher bands.
3065 const double dosX =
3067 GetConductionBandDensityOfStates(e, nValleysX) * nValleysL -
3068 GetConductionBandDensityOfStates(e, nValleysX + nValleysL);
3069 if (dosX <= 0.) return 0.;
3070 return dosX / nValleysX;
3071 }
3072 } else if (useNonParabolicity) {
3073 const double alpha = alphaX;
3074 return sqrt(md3 * e * (1. + alpha * e) / 2.) * (1. + 2 * alpha * e) /
3075 (Pi2 * pow(HbarC, 3.));
3076 } else {
3077 return sqrt(md3 * e / 2.) / (Pi2 * pow(HbarC, 3.));
3078 }
3079 } else if (band < nValleysX + nValleysL) {
3080 // L valleys
3081 if (e <= eMinL) return 0.;
3082
3083 // Density-of-states effective mass (cube)
3084 const double md3 = pow(ElectronMass, 3) * mLongL * mTransL * mTransL;
3085 // Non-parabolicity parameter
3086 const double alpha = alphaL;
3087
3088 if (useFullBandDos) {
3089 // Energy up to which the non-parabolic approximation is used.
3090 const double ej = eMinL + 0.5;
3091 if (e <= ej) {
3092 return sqrt(md3 * (e - eMinL) * (1. + alpha * (e - eMinL))) *
3093 (1. + 2 * alpha * (e - eMinL)) / (Sqrt2 * Pi2 * pow(HbarC, 3.));
3094 } else {
3095 // Fraction of full-band density of states attributed to L valleys
3096 double fL = sqrt(md3 * (ej - eMinL) * (1. + alpha * (ej - eMinL))) *
3097 (1. + 2 * alpha * (ej - eMinL)) /
3098 (Sqrt2 * Pi2 * pow(HbarC, 3.));
3099 fL = nValleysL * fL / GetConductionBandDensityOfStates(ej, -1);
3100
3101 double dosXL = GetConductionBandDensityOfStates(e, -1);
3102 if (e > eMinG) {
3103 dosXL -= GetConductionBandDensityOfStates(e, nValleysX + nValleysL);
3104 }
3105 if (dosXL <= 0.) return 0.;
3106 return fL * dosXL / 8.;
3107 }
3108 } else if (useNonParabolicity) {
3109 return sqrt(md3 * (e - eMinL) * (1. + alpha * (e - eMinL))) *
3110 (1. + 2 * alpha * (e - eMinL)) / (Sqrt2 * Pi2 * pow(HbarC, 3.));
3111 } else {
3112 return sqrt(md3 * (e - eMinL) / 2.) / (Pi2 * pow(HbarC, 3.));
3113 }
3114 } else if (band == nValleysX + nValleysL) {
3115 // Higher bands
3116 const double ej = 2.7;
3117 if (eMinG >= ej) {
3118 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n";
3119 std::cerr << " Cannot determine higher band density-of-states.\n";
3120 std::cerr << " Program bug. Check offset energy!\n";
3121 return 0.;
3122 }
3123 if (e < eMinG) {
3124 return 0.;
3125 } else if (e < ej) {
3126 // Coexistence of XL and higher bands.
3127 const double dj = GetConductionBandDensityOfStates(ej, -1);
3128 // Assume linear increase of density-of-states.
3129 return dj * (e - eMinG) / (ej - eMinG);
3130 } else {
3132 }
3133 }
3134
3135 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n";
3136 std::cerr << " Band index (" << band << ") out of range.\n";
3137 return ElectronMass * sqrt(ElectronMass * e / 2.) / (Pi2 * pow(HbarC, 3.));
3138}
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 1399 of file MediumSilicon.cc.

1400 {
1401
1402 if (i != 0) {
1403 std::cerr << m_className << "::GetDielectricFunction:\n";
1404 std::cerr << " Medium has only one component.\n";
1405 return false;
1406 }
1407
1408 // Make sure the optical data table has been loaded.
1409 if (!m_hasOpticalData) {
1410 if (!LoadOpticalData(opticalDataFile)) {
1411 std::cerr << m_className << "::GetDielectricFunction:\n";
1412 std::cerr << " Optical data table could not be loaded.\n";
1413 return false;
1414 }
1415 m_hasOpticalData = true;
1416 }
1417
1418 // Make sure the requested energy is within the range of the table.
1419 const double emin = opticalDataTable[0].energy;
1420 const double emax = opticalDataTable.back().energy;
1421 if (e < emin || e > emax) {
1422 std::cerr << m_className << "::GetDielectricFunction:\n";
1423 std::cerr << " Requested energy (" << e << " eV) "
1424 << " is outside the range of the optical data table.\n";
1425 std::cerr << " " << emin << " < E [eV] < " << emax << "\n";
1426 eps1 = eps2 = 0.;
1427 return false;
1428 }
1429
1430 // Locate the requested energy in the table.
1431 int iLow = 0;
1432 int iUp = opticalDataTable.size() - 1;
1433 int iM;
1434 while (iUp - iLow > 1) {
1435 iM = (iUp + iLow) >> 1;
1436 if (e >= opticalDataTable[iM].energy) {
1437 iLow = iM;
1438 } else {
1439 iUp = iM;
1440 }
1441 }
1442
1443 // Interpolate the real part of dielectric function.
1444 // Use linear interpolation if one of the values is negative,
1445 // Otherwise use log-log interpolation.
1446 const double logX0 = log(opticalDataTable[iLow].energy);
1447 const double logX1 = log(opticalDataTable[iUp].energy);
1448 const double logX = log(e);
1449 if (opticalDataTable[iLow].eps1 <= 0. || opticalDataTable[iUp].eps1 <= 0.) {
1450 eps1 = opticalDataTable[iLow].eps1 +
1451 (e - opticalDataTable[iLow].energy) *
1452 (opticalDataTable[iUp].eps1 - opticalDataTable[iLow].eps1) /
1453 (opticalDataTable[iUp].energy - opticalDataTable[iLow].energy);
1454 } else {
1455 const double logY0 = log(opticalDataTable[iLow].eps1);
1456 const double logY1 = log(opticalDataTable[iUp].eps1);
1457 eps1 = logY0 + (logX - logX0) * (logY1 - logY0) / (logX1 - logX0);
1458 eps1 = exp(eps1);
1459 }
1460
1461 // Interpolate the imaginary part of dielectric function,
1462 // using log-log interpolation.
1463 const double logY0 = log(opticalDataTable[iLow].eps2);
1464 const double logY1 = log(opticalDataTable[iUp].eps2);
1465 eps2 = logY0 + (log(e) - logX0) * (logY1 - logY0) / (logX1 - logX0);
1466 eps2 = exp(eps2);
1467 return true;
1468}
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:376

◆ GetDoping()

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

Definition at line 176 of file MediumSilicon.cc.

176 {
177
178 type = m_dopingType;
179 c = m_dopingConcentration;
180}

◆ GetElectronBandPopulation()

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

Definition at line 1361 of file MediumSilicon.cc.

1361 {
1362
1363 const int nBands = nValleysX + nValleysL + 1;
1364 if (band < 0 || band >= nBands) {
1365 std::cerr << m_className << "::GetElectronBandPopulation:\n";
1366 std::cerr << " Band index (" << band << ") out of range.\n";
1367 return 0;
1368 }
1369 return nCollElectronBand[band];
1370}

◆ 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 928 of file MediumSilicon.cc.

931 {
932
933 if (e > eFinalG) {
934 std::cerr << m_className << "::GetElectronCollision:\n";
935 std::cerr << " Requested electron energy (" << e;
936 std::cerr << " eV) exceeds current energy range (" << eFinalG;
937 std::cerr << " eV).\n";
938 std::cerr << " Increasing energy range to " << 1.05 * e << " eV.\n";
939 SetMaxElectronEnergy(1.05 * e);
940 } else if (e <= 0.) {
941 std::cerr << m_className << "::GetElectronCollision:\n";
942 std::cerr << " Electron energy must be greater than zero.\n";
943 return false;
944 }
945
946 if (m_isChanged) {
947 if (!UpdateTransportParameters()) {
948 std::cerr << m_className << "::GetElectronCollision:\n";
949 std::cerr << " Error calculating the collision rates table.\n";
950 return false;
951 }
952 m_isChanged = false;
953 }
954
955 // Energy loss
956 double loss = 0.;
957 // Sample the scattering process.
958 if (band >= 0 && band < nValleysX) {
959 // X valley
960 // Get the energy interval.
961 int iE = int(e / eStepXL);
962 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
963 if (iE < 0) iE = 0;
964 // Select the scattering process.
965 const double r = RndmUniform();
966 int iLow = 0;
967 int iUp = nLevelsX - 1;
968 if (r <= cfElectronsX[iE][iLow]) {
969 level = iLow;
970 } else if (r >= cfElectronsX[iE][nLevelsX - 1]) {
971 level = iUp;
972 } else {
973 int iMid;
974 while (iUp - iLow > 1) {
975 iMid = (iLow + iUp) >> 1;
976 if (r < cfElectronsX[iE][iMid]) {
977 iUp = iMid;
978 } else {
979 iLow = iMid;
980 }
981 }
982 level = iUp;
983 }
984
985 // Get the collision type.
986 type = scatTypeElectronsX[level];
987 // Fill the collision counters.
988 ++nCollElectronDetailed[level];
989 ++nCollElectronBand[band];
990 if (type == ElectronCollisionTypeAcousticPhonon) {
991 ++nCollElectronAcoustic;
992 } else if (type == ElectronCollisionTypeIntervalleyG) {
993 // Intervalley scattering (g type)
994 ++nCollElectronIntervalley;
995 // Final valley is in opposite direction.
996 switch (band) {
997 case 0:
998 band = 1;
999 break;
1000 case 1:
1001 band = 0;
1002 break;
1003 case 2:
1004 band = 3;
1005 break;
1006 case 3:
1007 band = 2;
1008 break;
1009 case 4:
1010 band = 5;
1011 break;
1012 case 5:
1013 band = 4;
1014 break;
1015 default:
1016 break;
1017 }
1018 } else if (type == ElectronCollisionTypeIntervalleyF) {
1019 // Intervalley scattering (f type)
1020 ++nCollElectronIntervalley;
1021 // Final valley is perpendicular.
1022 switch (band) {
1023 case 0:
1024 case 1:
1025 band = int(RndmUniform() * 4) + 2;
1026 break;
1027 case 2:
1028 case 3:
1029 band = int(RndmUniform() * 4);
1030 if (band > 1) band += 2;
1031 break;
1032 case 4:
1033 case 5:
1034 band = int(RndmUniform() * 4);
1035 break;
1036 }
1037 } else if (type == ElectronCollisionTypeInterbandXL) {
1038 // XL scattering
1039 ++nCollElectronIntervalley;
1040 // Final valley is in L band.
1041 band = nValleysX + int(RndmUniform() * nValleysL);
1042 if (band >= nValleysX + nValleysL) band = nValleysX + nValleysL - 1;
1043 } else if (type == ElectronCollisionTypeInterbandXG) {
1044 ++nCollElectronIntervalley;
1045 band = nValleysX + nValleysL;
1046 } else if (type == ElectronCollisionTypeImpurity) {
1047 ++nCollElectronImpurity;
1048 } else if (type == ElectronCollisionTypeIonisation) {
1049 ++nCollElectronIonisation;
1050 } else {
1051 std::cerr << m_className << "::GetElectronCollision:\n";
1052 std::cerr << " Unexpected collision type (" << type << ").\n";
1053 }
1054
1055 // Get the energy loss.
1056 loss = energyLossElectronsX[level];
1057
1058 } else if (band >= nValleysX && band < nValleysX + nValleysL) {
1059 // L valley
1060 // Get the energy interval.
1061 int iE = int(e / eStepXL);
1062 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
1063 if (iE < ieMinL) iE = ieMinL;
1064 // Select the scattering process.
1065 const double r = RndmUniform();
1066 int iLow = 0;
1067 int iUp = nLevelsL - 1;
1068 if (r <= cfElectronsL[iE][iLow]) {
1069 level = iLow;
1070 } else if (r >= cfElectronsL[iE][nLevelsL - 1]) {
1071 level = iUp;
1072 } else {
1073 int iMid;
1074 while (iUp - iLow > 1) {
1075 iMid = (iLow + iUp) >> 1;
1076 if (r < cfElectronsL[iE][iMid]) {
1077 iUp = iMid;
1078 } else {
1079 iLow = iMid;
1080 }
1081 }
1082 level = iUp;
1083 }
1084
1085 // Get the collision type.
1086 type = scatTypeElectronsL[level];
1087 // Fill the collision counters.
1088 ++nCollElectronDetailed[nLevelsX + level];
1089 ++nCollElectronBand[band];
1090 if (type == ElectronCollisionTypeAcousticPhonon) {
1091 ++nCollElectronAcoustic;
1092 } else if (type == ElectronCollisionTypeOpticalPhonon) {
1093 ++nCollElectronOptical;
1094 } else if (type == ElectronCollisionTypeIntervalleyG ||
1095 type == ElectronCollisionTypeIntervalleyF) {
1096 // Equivalent intervalley scattering
1097 ++nCollElectronIntervalley;
1098 // Randomise the final valley.
1099 band = nValleysX + int(RndmUniform() * nValleysL);
1100 if (band >= nValleysX + nValleysL) band = nValleysX + nValleysL;
1101 } else if (type == ElectronCollisionTypeInterbandXL) {
1102 // LX scattering
1103 ++nCollElectronIntervalley;
1104 // Randomise the final valley.
1105 band = int(RndmUniform() * nValleysX);
1106 if (band >= nValleysX) band = nValleysX - 1;
1107 } else if (type == ElectronCollisionTypeInterbandLG) {
1108 // LG scattering
1109 ++nCollElectronIntervalley;
1110 band = nValleysX + nValleysL;
1111 } else if (type == ElectronCollisionTypeImpurity) {
1112 ++nCollElectronImpurity;
1113 } else if (type == ElectronCollisionTypeIonisation) {
1114 ++nCollElectronIonisation;
1115 } else {
1116 std::cerr << m_className << "::GetElectronCollision:\n";
1117 std::cerr << " Unexpected collision type (" << type << ").\n";
1118 }
1119
1120 // Get the energy loss.
1121 loss = energyLossElectronsL[level];
1122 } else if (band == nValleysX + nValleysL) {
1123 // Higher bands
1124 // Get the energy interval.
1125 int iE = int(e / eStepG);
1126 if (iE >= nEnergyStepsG) iE = nEnergyStepsG - 1;
1127 if (iE < ieMinG) iE = ieMinG;
1128 // Select the scattering process.
1129 const double r = RndmUniform();
1130 int iLow = 0;
1131 int iUp = nLevelsG - 1;
1132 if (r <= cfElectronsG[iE][iLow]) {
1133 level = iLow;
1134 } else if (r >= cfElectronsG[iE][nLevelsG - 1]) {
1135 level = iUp;
1136 } else {
1137 int iMid;
1138 while (iUp - iLow > 1) {
1139 iMid = (iLow + iUp) >> 1;
1140 if (r < cfElectronsG[iE][iMid]) {
1141 iUp = iMid;
1142 } else {
1143 iLow = iMid;
1144 }
1145 }
1146 level = iUp;
1147 }
1148
1149 // Get the collision type.
1150 type = scatTypeElectronsG[level];
1151 // Fill the collision counters.
1152 ++nCollElectronDetailed[nLevelsX + nLevelsL + level];
1153 ++nCollElectronBand[band];
1154 if (type == ElectronCollisionTypeAcousticPhonon) {
1155 ++nCollElectronAcoustic;
1156 } else if (type == ElectronCollisionTypeOpticalPhonon) {
1157 ++nCollElectronOptical;
1158 } else if (type == ElectronCollisionTypeIntervalleyG ||
1159 type == ElectronCollisionTypeIntervalleyF) {
1160 // Equivalent intervalley scattering
1161 ++nCollElectronIntervalley;
1162 } else if (type == ElectronCollisionTypeInterbandXG) {
1163 // GX scattering
1164 ++nCollElectronIntervalley;
1165 // Randomise the final valley.
1166 band = int(RndmUniform() * nValleysX);
1167 if (band >= nValleysX) band = nValleysX - 1;
1168 } else if (type == ElectronCollisionTypeInterbandLG) {
1169 // GL scattering
1170 ++nCollElectronIntervalley;
1171 // Randomise the final valley.
1172 band = nValleysX + int(RndmUniform() * nValleysL);
1173 if (band >= nValleysX + nValleysL) band = nValleysX + nValleysL - 1;
1174 } else if (type == ElectronCollisionTypeIonisation) {
1175 ++nCollElectronIonisation;
1176 } else {
1177 std::cerr << m_className << "::GetElectronCollision:\n";
1178 std::cerr << " Unexpected collision type (" << type << ").\n";
1179 }
1180
1181 // Get the energy loss.
1182 loss = energyLossElectronsG[level];
1183 } else {
1184 std::cerr << m_className << "::GetElectronCollision:\n";
1185 std::cerr << " Band index (" << band << ") out of range.\n";
1186 return false;
1187 }
1188
1189 // Secondaries
1190 nion = ndxc = 0;
1191 // Ionising collision
1192 if (type == ElectronCollisionTypeIonisation) {
1193 double ee = 0., eh = 0.;
1194 ComputeSecondaries(e, ee, eh);
1195 loss = ee + eh + m_bandGap;
1196 ionProducts.clear();
1197 // Add the secondary electron.
1198 ionProd newIonProd;
1199 newIonProd.type = IonProdTypeElectron;
1200 newIonProd.energy = ee;
1201 ionProducts.push_back(newIonProd);
1202 // Add the hole.
1203 newIonProd.type = IonProdTypeHole;
1204 newIonProd.energy = eh;
1205 ionProducts.push_back(newIonProd);
1206 nion = nIonisationProducts = 2;
1207 }
1208
1209 if (e < loss) loss = e - 0.00001;
1210 // Update the energy.
1211 e1 = e - loss;
1212 if (e1 < Small) e1 = Small;
1213
1214 // Update the momentum.
1215 if (band >= 0 && band < nValleysX) {
1216 // X valleys
1217 double pstar = sqrt(2. * ElectronMass * e1);
1218 if (useNonParabolicity) {
1219 const double alpha = alphaX;
1220 pstar *= sqrt(1. + alpha * e1);
1221 }
1222
1223 const double ctheta = 1. - 2. * RndmUniform();
1224 const double stheta = sqrt(1. - ctheta * ctheta);
1225 const double phi = TwoPi * RndmUniform();
1226
1227 if (useAnisotropy) {
1228 const double pl = pstar * sqrt(mLongX);
1229 const double pt = pstar * sqrt(mTransX);
1230 switch (band) {
1231 case 0:
1232 case 1:
1233 // 100
1234 px = pl * ctheta;
1235 py = pt * cos(phi) * stheta;
1236 pz = pt * sin(phi) * stheta;
1237 break;
1238 case 2:
1239 case 3:
1240 // 010
1241 px = pt * sin(phi) * stheta;
1242 py = pl * ctheta;
1243 pz = pt * cos(phi) * stheta;
1244 break;
1245 case 4:
1246 case 5:
1247 // 001
1248 px = pt * cos(phi) * stheta;
1249 py = pt * sin(phi) * stheta;
1250 pz = pl * ctheta;
1251 break;
1252 default:
1253 return false;
1254 break;
1255 }
1256 } else {
1257 pstar *= sqrt(3. / (1. / mLongX + 2. / mTransX));
1258 px = pstar * cos(phi) * stheta;
1259 py = pstar * sin(phi) * stheta;
1260 pz = pstar * ctheta;
1261 }
1262 return true;
1263
1264 } else if (band >= nValleysX && band < nValleysX + nValleysL) {
1265 // L valleys
1266 double pstar = sqrt(2. * ElectronMass * (e1 - eMinL));
1267 if (useNonParabolicity) {
1268 const double alpha = alphaL;
1269 pstar *= sqrt(1. + alpha * (e1 - eMinL));
1270 }
1271
1272 const double ctheta = 1. - 2. * RndmUniform();
1273 const double stheta = sqrt(1. - ctheta * ctheta);
1274 const double phi = TwoPi * RndmUniform();
1275
1276 pstar *= sqrt(3. / (1. / mLongL + 2. / mTransL));
1277 px = pstar * cos(phi) * stheta;
1278 py = pstar * sin(phi) * stheta;
1279 pz = pstar * ctheta;
1280 return true;
1281
1282 } else {
1283 double pstar = sqrt(2. * ElectronMass * e1);
1284
1285 const double ctheta = 1. - 2. * RndmUniform();
1286 const double stheta = sqrt(1. - ctheta * ctheta);
1287 const double phi = TwoPi * RndmUniform();
1288
1289 px = pstar * cos(phi) * stheta;
1290 py = pstar * sin(phi) * stheta;
1291 pz = pstar * ctheta;
1292 return true;
1293 }
1294
1295 std::cerr << m_className << "::GetElectronCollision:\n";
1296 std::cerr << " Band index (" << band << ") out of range.\n";
1297 e1 = e;
1298 type = 0;
1299 return false;
1300}
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:431
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:383
void ComputeSecondaries(const double e0, double &ee, double &eh)
bool SetMaxElectronEnergy(const double e)

◆ GetElectronCollisionRate()

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

Reimplemented from Garfield::Medium.

Definition at line 874 of file MediumSilicon.cc.

874 {
875
876 if (e <= 0.) {
877 std::cerr << m_className << "::GetElectronCollisionRate:\n";
878 std::cerr << " Electron energy must be positive.\n";
879 return 0.;
880 }
881
882 if (e > eFinalG) {
883 std::cerr << m_className << "::GetElectronCollisionRate:\n";
884 std::cerr << " Collision rate at " << e << " eV (band " << band
885 << ") is not included"
886 << " in the current table.\n";
887 std::cerr << " Increasing energy range to " << 1.05 * e << " eV.\n";
888 SetMaxElectronEnergy(1.05 * e);
889 }
890
891 if (m_isChanged) {
892 if (!UpdateTransportParameters()) {
893 std::cerr << m_className << "::GetElectronCollisionRate:\n";
894 std::cerr << " Error calculating the collision rates table.\n";
895 return 0.;
896 }
897 m_isChanged = false;
898 }
899
900 if (band >= 0 && band < nValleysX) {
901 int iE = int(e / eStepXL);
902 if (iE >= nEnergyStepsXL)
903 iE = nEnergyStepsXL - 1;
904 else if (iE < 0)
905 iE = 0;
906 return cfTotElectronsX[iE];
907 } else if (band >= nValleysX && band < nValleysX + nValleysL) {
908 int iE = int(e / eStepXL);
909 if (iE >= nEnergyStepsXL)
910 iE = nEnergyStepsXL - 1;
911 else if (iE < ieMinL)
912 iE = ieMinL;
913 return cfTotElectronsL[iE];
914 } else if (band == nValleysX + nValleysL) {
915 int iE = int(e / eStepG);
916 if (iE >= nEnergyStepsG)
917 iE = nEnergyStepsG - 1;
918 else if (iE < ieMinG)
919 iE = ieMinG;
920 return cfTotElectronsG[iE];
921 }
922
923 std::cerr << m_className << "::GetElectronCollisionRate:\n";
924 std::cerr << " Band index (" << band << ") out of range.\n";
925 return 0.;
926}

◆ 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 643 of file MediumSilicon.cc.

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

◆ GetElectronMomentum()

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

Reimplemented from Garfield::Medium.

Definition at line 729 of file MediumSilicon.cc.

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

◆ GetElectronNullCollisionRate()

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

Reimplemented from Garfield::Medium.

Definition at line 851 of file MediumSilicon.cc.

851 {
852
853 if (m_isChanged) {
854 if (!UpdateTransportParameters()) {
855 std::cerr << m_className << "::GetElectronNullCollisionRate:\n";
856 std::cerr << " Error calculating the collision rates table.\n";
857 return 0.;
858 }
859 m_isChanged = false;
860 }
861
862 if (band >= 0 && band < nValleysX) {
863 return cfNullElectronsX;
864 } else if (band >= nValleysX && band < nValleysX + nValleysL) {
865 return cfNullElectronsL;
866 } else if (band == nValleysX + nValleysL) {
867 return cfNullElectronsG;
868 }
869 std::cerr << m_className << "::GetElectronNullCollisionRate:\n";
870 std::cerr << " Band index (" << band << ") out of range.\n";
871 return 0.;
872}

◆ GetIonisationProduct()

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

Reimplemented from Garfield::Medium.

Definition at line 1302 of file MediumSilicon.cc.

1303 {
1304
1305 if (i < 0 || i >= nIonisationProducts) {
1306 std::cerr << m_className << "::GetIonisationProduct:\n";
1307 std::cerr << " Index (" << i << ") out of range.\n";
1308 return false;
1309 }
1310
1311 type = ionProducts[i].type;
1312 energy = ionProducts[i].energy;
1313 return true;
1314}

◆ GetMaxElectronEnergy()

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

Definition at line 73 of file MediumSilicon.hh.

73{ return eFinalG; }

◆ GetNumberOfElectronBands()

int Garfield::MediumSilicon::GetNumberOfElectronBands ( )

Definition at line 1355 of file MediumSilicon.cc.

1355 {
1356
1357 const int nBands = nValleysX + nValleysL + 1;
1358 return nBands;
1359}

◆ GetNumberOfElectronCollisions() [1/2]

int Garfield::MediumSilicon::GetNumberOfElectronCollisions ( ) const

Definition at line 1330 of file MediumSilicon.cc.

1330 {
1331
1332 return nCollElectronAcoustic + nCollElectronOptical +
1333 nCollElectronIntervalley + nCollElectronImpurity +
1334 nCollElectronIonisation;
1335}

◆ GetNumberOfElectronCollisions() [2/2]

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

Definition at line 1343 of file MediumSilicon.cc.

1343 {
1344
1345 const int nLevels = nLevelsX + nLevelsL + nLevelsG;
1346 if (level < 0 || level >= nLevels) {
1347 std::cerr << m_className << "::GetNumberOfElectronCollisions:\n";
1348 std::cerr << " Requested scattering rate term (" << level
1349 << ") does not exist.\n";
1350 return 0;
1351 }
1352 return nCollElectronDetailed[level];
1353}

◆ GetNumberOfIonisationProducts()

int Garfield::MediumSilicon::GetNumberOfIonisationProducts ( )
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 106 of file MediumSilicon.hh.

106{ return nIonisationProducts; }

◆ GetNumberOfLevels()

int Garfield::MediumSilicon::GetNumberOfLevels ( )

Definition at line 1337 of file MediumSilicon.cc.

1337 {
1338
1339 const int nLevels = nLevelsX + nLevelsL + nLevelsG;
1340 return nLevels;
1341}

◆ GetOpticalDataRange()

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

Reimplemented from Garfield::Medium.

Definition at line 1372 of file MediumSilicon.cc.

1373 {
1374
1375 if (i != 0) {
1376 std::cerr << m_className << "::GetOpticalDataRange:\n";
1377 std::cerr << " Medium has only one component.\n";
1378 }
1379
1380 // Make sure the optical data table has been loaded.
1381 if (!m_hasOpticalData) {
1382 if (!LoadOpticalData(opticalDataFile)) {
1383 std::cerr << m_className << "::GetOpticalDataRange:\n";
1384 std::cerr << " Optical data table could not be loaded.\n";
1385 return false;
1386 }
1387 m_hasOpticalData = true;
1388 }
1389
1390 emin = opticalDataTable[0].energy;
1391 emax = opticalDataTable.back().energy;
1392 if (m_debug) {
1393 std::cout << m_className << "::GetOpticalDataRange:\n";
1394 std::cout << " " << emin << " < E [eV] < " << emax << "\n";
1395 }
1396 return true;
1397}

◆ GetValenceBandDensityOfStates()

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

Definition at line 3140 of file MediumSilicon.cc.

3141 {
3142
3143 if (band <= 0) {
3144 // Total (full-band) density of states.
3145
3146 int iE = int(e / eStepDos);
3147 if (iE >= nFbDosEntriesValence || iE < 0) {
3148 return 0.;
3149 } else if (iE == nFbDosEntriesValence - 1) {
3150 return fbDosValence[nFbDosEntriesValence - 1];
3151 }
3152
3153 const double dos =
3154 fbDosValence[iE] +
3155 (fbDosValence[iE + 1] - fbDosValence[iE]) * (e / eStepDos - iE);
3156 return dos * 1.e21;
3157 }
3158
3159 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n";
3160 std::cerr << " Band index (" << band << ") out of range.\n";
3161 return 0.;
3162}

◆ 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 464 of file MediumSilicon.cc.

467 {
468
469 eta = 0.;
470 if (m_isChanged) {
471 if (!UpdateTransportParameters()) {
472 std::cerr << m_className << "::HoleAttachment:\n";
473 std::cerr << " Error calculating the transport parameters.\n";
474 return false;
475 }
476 m_isChanged = false;
477 }
478
480 // Interpolation in user table.
481 return Medium::HoleAttachment(ex, ey, ez, bx, by, bz, eta);
482 }
483
484 switch (trappingModel) {
485 case 0:
486 eta = hTrapCs * hTrapDensity;
487 break;
488 case 1:
489 double vx, vy, vz;
490 HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
491 eta = hTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
492 if (eta > 0.) eta = 1. / eta;
493 break;
494 default:
495 std::cerr << m_className << "::HoleAttachment:\n";
496 std::cerr << " Unknown model activated. Program bug!\n";
497 return false;
498 break;
499 }
500
501 return true;
502}
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:350
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:1144

◆ 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 427 of file MediumSilicon.cc.

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

◆ 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 371 of file MediumSilicon.cc.

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

Referenced by HoleAttachment().

◆ Initialise()

bool Garfield::MediumSilicon::Initialise ( )

Definition at line 1470 of file MediumSilicon.cc.

1470 {
1471
1472 if (!m_isChanged) {
1473 if (m_debug) {
1474 std::cerr << m_className << "::Initialise:\n";
1475 std::cerr << " Nothing changed.\n";
1476 }
1477 return true;
1478 }
1479 if (!UpdateTransportParameters()) {
1480 std::cerr << m_className << "::Initialise:\n";
1481 std::cerr << " Error preparing transport parameters/"
1482 << "calculating collision rates.\n";
1483 return false;
1484 }
1485 return true;
1486}

◆ IsSemiconductor()

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

Reimplemented from Garfield::Medium.

Definition at line 18 of file MediumSilicon.hh.

18{ return true; }

◆ ResetCollisionCounters()

void Garfield::MediumSilicon::ResetCollisionCounters ( )

Definition at line 1316 of file MediumSilicon.cc.

1316 {
1317
1318 nCollElectronAcoustic = nCollElectronOptical = 0;
1319 nCollElectronIntervalley = 0;
1320 nCollElectronImpurity = 0;
1321 nCollElectronIonisation = 0;
1322 const int nLevels = nLevelsX + nLevelsL + nLevelsG;
1323 nCollElectronDetailed.resize(nLevels);
1324 for (int j = nLevels; j--;) nCollElectronDetailed[j] = 0;
1325 const int nBands = nValleysX + nValleysL + 1;
1326 nCollElectronBand.resize(nBands);
1327 for (int j = nBands; j--;) nCollElectronBand[j] = 0;
1328}

◆ SetDoping()

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

Definition at line 139 of file MediumSilicon.cc.

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

◆ SetDopingMobilityModelMasetti()

void Garfield::MediumSilicon::SetDopingMobilityModelMasetti ( )

Definition at line 546 of file MediumSilicon.cc.

546 {
547
548 dopingMobilityModel = DopingMobilityModelMasetti;
549 m_hasUserMobility = false;
550 m_isChanged = true;
551}

◆ SetDopingMobilityModelMinimos()

void Garfield::MediumSilicon::SetDopingMobilityModelMinimos ( )

Definition at line 539 of file MediumSilicon.cc.

539 {
540
541 dopingMobilityModel = DopingMobilityModelMinimos;
542 m_hasUserMobility = false;
543 m_isChanged = true;
544}

◆ SetHighFieldMobilityModelCanali()

void Garfield::MediumSilicon::SetHighFieldMobilityModelCanali ( )

Definition at line 596 of file MediumSilicon.cc.

596 {
597
598 highFieldMobilityModel = HighFieldMobilityModelCanali;
599 m_isChanged = true;
600}

◆ SetHighFieldMobilityModelConstant()

void Garfield::MediumSilicon::SetHighFieldMobilityModelConstant ( )

Definition at line 608 of file MediumSilicon.cc.

608 {
609
610 highFieldMobilityModel = HighFieldMobilityModelConstant;
611}

◆ SetHighFieldMobilityModelMinimos()

void Garfield::MediumSilicon::SetHighFieldMobilityModelMinimos ( )

Definition at line 590 of file MediumSilicon.cc.

590 {
591
592 highFieldMobilityModel = HighFieldMobilityModelMinimos;
593 m_isChanged = true;
594}

◆ SetHighFieldMobilityModelReggiani()

void Garfield::MediumSilicon::SetHighFieldMobilityModelReggiani ( )

Definition at line 602 of file MediumSilicon.cc.

602 {
603
604 highFieldMobilityModel = HighFieldMobilityModelReggiani;
605 m_isChanged = true;
606}

◆ SetImpactIonisationModelGrant()

void Garfield::MediumSilicon::SetImpactIonisationModelGrant ( )

Definition at line 619 of file MediumSilicon.cc.

619 {
620
621 impactIonisationModel = ImpactIonisationModelGrant;
622 m_isChanged = true;
623}

◆ SetImpactIonisationModelVanOverstraetenDeMan()

void Garfield::MediumSilicon::SetImpactIonisationModelVanOverstraetenDeMan ( )

Definition at line 613 of file MediumSilicon.cc.

613 {
614
615 impactIonisationModel = ImpactIonisationModelVanOverstraeten;
616 m_isChanged = true;
617}

◆ SetLatticeMobilityModelMinimos()

void Garfield::MediumSilicon::SetLatticeMobilityModelMinimos ( )

Definition at line 518 of file MediumSilicon.cc.

518 {
519
520 latticeMobilityModel = LatticeMobilityModelMinimos;
521 m_hasUserMobility = false;
522 m_isChanged = true;
523}

◆ SetLatticeMobilityModelReggiani()

void Garfield::MediumSilicon::SetLatticeMobilityModelReggiani ( )

Definition at line 532 of file MediumSilicon.cc.

532 {
533
534 latticeMobilityModel = LatticeMobilityModelReggiani;
535 m_hasUserMobility = false;
536 m_isChanged = true;
537}

◆ SetLatticeMobilityModelSentaurus()

void Garfield::MediumSilicon::SetLatticeMobilityModelSentaurus ( )

Definition at line 525 of file MediumSilicon.cc.

525 {
526
527 latticeMobilityModel = LatticeMobilityModelSentaurus;
528 m_hasUserMobility = false;
529 m_isChanged = true;
530}

◆ SetLowFieldMobility()

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

Definition at line 504 of file MediumSilicon.cc.

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

◆ SetMaxElectronEnergy()

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

Definition at line 625 of file MediumSilicon.cc.

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

Referenced by GetElectronCollision(), and GetElectronCollisionRate().

◆ SetSaturationVelocity()

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

Definition at line 553 of file MediumSilicon.cc.

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

◆ SetSaturationVelocityModelCanali()

void Garfield::MediumSilicon::SetSaturationVelocityModelCanali ( )

Definition at line 576 of file MediumSilicon.cc.

576 {
577
578 saturationVelocityModel = SaturationVelocityModelCanali;
579 m_hasUserSaturationVelocity = false;
580 m_isChanged = true;
581}

◆ SetSaturationVelocityModelMinimos()

void Garfield::MediumSilicon::SetSaturationVelocityModelMinimos ( )

Definition at line 569 of file MediumSilicon.cc.

569 {
570
571 saturationVelocityModel = SaturationVelocityModelMinimos;
572 m_hasUserSaturationVelocity = false;
573 m_isChanged = true;
574}

◆ SetSaturationVelocityModelReggiani()

void Garfield::MediumSilicon::SetSaturationVelocityModelReggiani ( )

Definition at line 583 of file MediumSilicon.cc.

583 {
584
585 saturationVelocityModel = SaturationVelocityModelReggiani;
586 m_hasUserSaturationVelocity = false;
587 m_isChanged = true;
588}

◆ SetTrapCrossSection()

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

Definition at line 182 of file MediumSilicon.cc.

182 {
183
184 if (ecs < 0.) {
185 std::cerr << m_className << "::SetTrapCrossSection:\n";
186 std::cerr << " Capture cross-section [cm2] must positive.\n";
187 } else {
188 eTrapCs = ecs;
189 }
190
191 if (hcs < 0.) {
192 std::cerr << m_className << "::SetTrapCrossSection:\n";
193 std::cerr << " Capture cross-section [cm2] must be positive.n";
194 } else {
195 hTrapCs = hcs;
196 }
197
198 trappingModel = 0;
199 m_isChanged = true;
200}

◆ SetTrapDensity()

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

Definition at line 202 of file MediumSilicon.cc.

202 {
203
204 if (n < 0.) {
205 std::cerr << m_className << "::SetTrapDensity:\n";
206 std::cerr << " Trap density [cm-3] must be greater than zero.\n";
207 } else {
208 eTrapDensity = n;
209 hTrapDensity = n;
210 }
211
212 trappingModel = 0;
213 m_isChanged = true;
214}

◆ SetTrappingTime()

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

Definition at line 216 of file MediumSilicon.cc.

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

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