Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
MediumSilicon.hh
Go to the documentation of this file.
1// Solid crystalline silicon
2
3#ifndef G_MEDIUM_SILICON_H
4#define G_MEDIUM_SILICON_H
5
6#include "Medium.hh"
7
8namespace Garfield {
9
10class MediumSilicon : public Medium {
11
12 public:
13 // Constructor
15 // Destructor
17
18 bool IsSemiconductor() const { return true; }
19
20 // Doping concentration [cm-3] and type ('i', 'n', 'p')
21 void SetDoping(const char& type, const double& c);
22 void GetDoping(char& type, double& c) const;
23
24 // Trapping cross-section
25 void SetTrapCrossSection(const double& ecs, const double& hcs);
26 void SetTrapDensity(const double& n);
27 void SetTrappingTime(const double& etau, const double& htau);
28
29 // Electron transport parameters
30 bool ElectronVelocity(const double ex, const double ey, const double ez,
31 const double bx, const double by, const double bz,
32 double& vx, double& vy, double& vz);
33 bool ElectronTownsend(const double ex, const double ey, const double ez,
34 const double bx, const double by, const double bz,
35 double& alpha);
36 bool ElectronAttachment(const double ex, const double ey, const double ez,
37 const double bx, const double by, const double bz,
38 double& eta);
39 // Hole transport parameters
40 bool HoleVelocity(const double ex, const double ey, const double ez,
41 const double bx, const double by, const double bz,
42 double& vx, double& vy, double& vz);
43 bool HoleTownsend(const double ex, const double ey, const double ez,
44 const double bx, const double by, const double bz,
45 double& alpha);
46 bool HoleAttachment(const double ex, const double ey, const double ez,
47 const double bx, const double by, const double bz,
48 double& eta);
49
50 void SetLowFieldMobility(const double mue, const double muh);
54
57
58 void SetSaturationVelocity(const double vsate, const double vsath);
62
67
70
71 // Microscopic transport properties
72 bool SetMaxElectronEnergy(const double e);
73 double GetMaxElectronEnergy() const { return eFinalG; }
74
75 bool Initialise();
76
77 // When enabled, the scattering rates table is written to file
78 // when loaded into memory.
79 void EnableScatteringRateOutput() { useCfOutput = true; }
80 void DisableScatteringRateOutput() { useCfOutput = false; }
81
82 void EnableNonParabolicity() { useNonParabolicity = true; }
83 void DisableNonParabolicity() { useNonParabolicity = false; }
84 void EnableFullBandDensityOfStates() { useFullBandDos = true; }
85 void DisableFullBandDensityOfStates() { useFullBandDos = false; }
86 void EnableAnisotropy() { useAnisotropy = true; }
87 void DisableAnisotropy() { useAnisotropy = false; }
88
89 // Get the electron energy (and its gradient)
90 // for a given (crystal) momentum
91 double GetElectronEnergy(const double px, const double py, const double pz,
92 double& vx, double& vy, double& vz,
93 const int band = 0);
94 // Get the electron (crystal) momentum for a given kinetic energy
95 void GetElectronMomentum(const double e, double& px, double& py, double& pz,
96 int& band);
97
98 // Get the null-collision rate [ns-1]
99 double GetElectronNullCollisionRate(const int band);
100 // Get the (real) collision rate [ns-1] at a given electron energy
101 double GetElectronCollisionRate(const double e, const int band);
102 // Sample the collision type
103 bool GetElectronCollision(const double e, int& type, int& level, double& e1,
104 double& dx, double& dy, double& dz, int& nion,
105 int& ndxc, int& band);
106 int GetNumberOfIonisationProducts() { return nIonisationProducts; }
107 bool GetIonisationProduct(const int i, int& type, double& energy);
108
109 // Density of states
110 double GetConductionBandDensityOfStates(const double e, const int band = 0);
111 double GetValenceBandDensityOfStates(const double e, const int band = -1);
112
113 // Reset the collision counters
115 // Get the total number of electron collisions
117 // Get number of scattering rate terms
118 int GetNumberOfLevels();
119 // Get number of collisions for a specific level
120 int GetNumberOfElectronCollisions(const int level) const;
121
123 int GetElectronBandPopulation(const int band);
124
125 bool GetOpticalDataRange(double& emin, double& emax,
126 const unsigned int& i = 0);
127 bool GetDielectricFunction(const double& e, double& eps1, double& eps2,
128 const unsigned int& i = 0);
129
130 void ComputeSecondaries(const double e0, double& ee, double& eh);
131
132 private:
133 static const int LatticeMobilityModelSentaurus = 0;
134 static const int LatticeMobilityModelMinimos = 1;
135 static const int LatticeMobilityModelReggiani = 2;
136 static const int DopingMobilityModelMinimos = 0;
137 static const int DopingMobilityModelMasetti = 1;
138 static const int SaturationVelocityModelMinimos = 0;
139 static const int SaturationVelocityModelCanali = 1;
140 static const int SaturationVelocityModelReggiani = 2;
141 static const int HighFieldMobilityModelMinimos = 0;
142 static const int HighFieldMobilityModelCanali = 1;
143 static const int HighFieldMobilityModelReggiani = 2;
144 static const int HighFieldMobilityModelConstant = 3;
145 static const int ImpactIonisationModelVanOverstraeten = 0;
146 static const int ImpactIonisationModelGrant = 1;
147
148 double m_bandGap;
149 // Doping
150 char m_dopingType;
151 double m_dopingConcentration;
152
153 // Effective masses
154 // X valleys
155 double mLongX, mTransX;
156 // L valleys
157 double mLongL, mTransL;
158 // Non-parabolicity parameters [1/eV]
159 double alphaX, alphaL;
160 // Lattice mobility
161 double eLatticeMobility, hLatticeMobility;
162 // Low-field mobility
163 double eMobility, hMobility;
164 // High-field mobility parameters
165 double eBetaCanali, hBetaCanali;
166 double eBetaCanaliInv, hBetaCanaliInv;
167 // Saturation velocity
168 double eSatVel, hSatVel;
169 // Hall factor
170 double eHallFactor, hHallFactor;
171
172 // Trapping parameters
173 double eTrapCs, hTrapCs;
174 double eTrapDensity, hTrapDensity;
175 double eTrapTime, hTrapTime;
176 int trappingModel;
177
178 // Impact ionisation parameters
179 double eImpactA0, eImpactA1, eImpactA2;
180 double eImpactB0, eImpactB1, eImpactB2;
181 double hImpactA0, hImpactA1;
182 double hImpactB0, hImpactB1;
183
184 // Models
185 bool m_hasUserMobility;
186 bool m_hasUserSaturationVelocity;
187 int latticeMobilityModel;
188 int dopingMobilityModel;
189 int saturationVelocityModel;
190 int highFieldMobilityModel;
191 int impactIonisationModel;
192
193 // Options
194 bool useCfOutput;
195 bool useNonParabolicity;
196 bool useFullBandDos;
197 bool useAnisotropy;
198
199 // Energy range of scattering rates
200 double eFinalXL, eStepXL;
201 double eFinalG, eStepG;
202 double eFinalV, eStepV;
203 static const int nEnergyStepsXL = 2000;
204 static const int nEnergyStepsG = 2000;
205 static const int nEnergyStepsV = 2000;
206
207 // Number of scattering terms
208 int nLevelsX, nLevelsL, nLevelsG;
209 int nLevelsV;
210 // Number of valleys
211 int nValleysX, nValleysL;
212 // Energy offset
213 double eMinL, eMinG;
214 int ieMinL, ieMinG;
215
216 // Electron scattering rates
217 double cfNullElectronsX, cfNullElectronsL, cfNullElectronsG;
218 std::vector<double> cfTotElectronsX;
219 std::vector<double> cfTotElectronsL;
220 std::vector<double> cfTotElectronsG;
221 std::vector<std::vector<double> > cfElectronsX;
222 std::vector<std::vector<double> > cfElectronsL;
223 std::vector<std::vector<double> > cfElectronsG;
224 std::vector<double> energyLossElectronsX;
225 std::vector<double> energyLossElectronsL;
226 std::vector<double> energyLossElectronsG;
227 // Cross-section type
228 std::vector<int> scatTypeElectronsX;
229 std::vector<int> scatTypeElectronsL;
230 std::vector<int> scatTypeElectronsG;
231
232 // Hole scattering rates
233 double cfNullHoles;
234 std::vector<double> cfTotHoles;
235 std::vector<std::vector<double> > cfHoles;
236 std::vector<double> energyLossHoles;
237 // Cross-section type
238 std::vector<int> scatTypeHoles;
239
240 // Collision counters
241 int nCollElectronAcoustic, nCollElectronOptical;
242 int nCollElectronIntervalley;
243 int nCollElectronImpurity;
244 int nCollElectronIonisation;
245 std::vector<int> nCollElectronDetailed;
246 std::vector<int> nCollElectronBand;
247
248 int nIonisationProducts;
249 struct ionProd {
250 int type;
251 double energy;
252 };
253 std::vector<ionProd> ionProducts;
254
255 // Density of states tables
256 double eStepDos;
257 int nFbDosEntriesValence;
258 int nFbDosEntriesConduction;
259 std::vector<double> fbDosValence;
260 std::vector<double> fbDosConduction;
261 double fbDosMaxV, fbDosMaxC;
262
263 // Optical data
264 bool m_hasOpticalData;
265 std::string opticalDataFile;
266 struct opticalData {
267 // Energy [eV]
268 double energy;
269 // Dielectric function
270 double eps1, eps2;
271 };
272 std::vector<opticalData> opticalDataTable;
273
274 bool UpdateTransportParameters();
275 void UpdateLatticeMobilityMinimos();
276 void UpdateLatticeMobilitySentaurus();
277 void UpdateLatticeMobilityReggiani();
278
279 void UpdateDopingMobilityMinimos();
280 void UpdateDopingMobilityMasetti();
281
282 void UpdateSaturationVelocityMinimos();
283 void UpdateSaturationVelocityCanali();
284 void UpdateSaturationVelocityReggiani();
285
286 void UpdateHighFieldMobilityCanali();
287
288 void UpdateImpactIonisationVanOverstraetenDeMan();
289 void UpdateImpactIonisationGrant();
290
291 bool ElectronMobilityMinimos(const double e, double& mu) const;
292 bool ElectronMobilityCanali(const double e, double& mu) const;
293 bool ElectronMobilityReggiani(const double e, double& mu) const;
294 bool ElectronImpactIonisationVanOverstraetenDeMan(const double e,
295 double& alpha) const;
296 bool ElectronImpactIonisationGrant(const double e, double& alpha) const;
297 bool HoleMobilityMinimos(const double e, double& mu) const;
298 bool HoleMobilityCanali(const double e, double& mu) const;
299 bool HoleMobilityReggiani(const double e, double& mu) const;
300 bool HoleImpactIonisationVanOverstraetenDeMan(const double e,
301 double& alpha) const;
302 bool HoleImpactIonisationGrant(const double e, double& alpha) const;
303
304 bool LoadOpticalData(const std::string filename);
305
306 bool ElectronScatteringRates();
307 bool ElectronAcousticScatteringRates();
308 bool ElectronOpticalScatteringRates();
309 bool ElectronIntervalleyScatteringRatesXX();
310 bool ElectronIntervalleyScatteringRatesXL();
311 bool ElectronIntervalleyScatteringRatesLL();
312 bool ElectronIntervalleyScatteringRatesXGLG();
313 bool ElectronIonisationRatesXL();
314 bool ElectronIonisationRatesG();
315 bool ElectronImpurityScatteringRates();
316
317 bool HoleScatteringRates();
318 bool HoleAcousticScatteringRates();
319 bool HoleOpticalScatteringRates();
320 bool HoleIonisationRates();
321
322 // void ComputeSecondaries(const double e0,
323 // double& ee, double& eh);
324 void InitialiseDensityOfStates();
325};
326}
327
328#endif
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)
void SetTrappingTime(const double &etau, const double &htau)
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)
void GetDoping(char &type, double &c) const
void SetHighFieldMobilityModelReggiani()
bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
void SetSaturationVelocityModelReggiani()
double GetValenceBandDensityOfStates(const double e, const int band=-1)
double GetElectronNullCollisionRate(const int band)
bool GetDielectricFunction(const double &e, double &eps1, double &eps2, const unsigned int &i=0)
bool IsSemiconductor() const
void SetSaturationVelocity(const double vsate, const double vsath)
bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
void ComputeSecondaries(const double e0, double &ee, double &eh)
bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
double GetElectronEnergy(const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
void SetDoping(const char &type, const double &c)
void SetTrapDensity(const double &n)
double GetConductionBandDensityOfStates(const double e, const int band=0)
bool SetMaxElectronEnergy(const double e)
void SetHighFieldMobilityModelConstant()
bool GetIonisationProduct(const int i, int &type, double &energy)
void SetImpactIonisationModelVanOverstraetenDeMan()
bool GetOpticalDataRange(double &emin, double &emax, const unsigned int &i=0)
void SetTrapCrossSection(const double &ecs, const double &hcs)
void SetSaturationVelocityModelMinimos()
void SetLowFieldMobility(const double mue, const double muh)
double GetMaxElectronEnergy() const
double GetElectronCollisionRate(const double e, const int band)
bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
void GetElectronMomentum(const double e, double &px, double &py, double &pz, int &band)
int GetNumberOfElectronCollisions() const
int GetElectronBandPopulation(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)