Garfield++ v2r0
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#ifndef G_MEDIUM_SILICON_H
2#define G_MEDIUM_SILICON_H
3
4#include "Medium.hh"
5
6namespace Garfield {
7/// %Solid crystalline silicon
8
9class MediumSilicon : public Medium {
10
11 public:
12 /// Constructor
14 /// Destructor
15 virtual ~MediumSilicon() {}
16
17 bool IsSemiconductor() const { return true; }
18
19 /// Doping concentration [cm-3] and type ('i', 'n', 'p')
20 void SetDoping(const char type, const double c);
21 void GetDoping(char& type, double& c) const;
22
23 /// Trapping cross-sections for electrons and holes.
24 void SetTrapCrossSection(const double ecs, const double hcs);
25 void SetTrapDensity(const double n);
26 void SetTrappingTime(const double etau, const double htau);
27
28 // Electron transport parameters
29 bool ElectronVelocity(const double ex, const double ey, const double ez,
30 const double bx, const double by, const double bz,
31 double& vx, double& vy, double& vz);
32 bool ElectronTownsend(const double ex, const double ey, const double ez,
33 const double bx, const double by, const double bz,
34 double& alpha);
35 bool ElectronAttachment(const double ex, const double ey, const double ez,
36 const double bx, const double by, const double bz,
37 double& eta);
38 // Hole transport parameters
39 bool HoleVelocity(const double ex, const double ey, const double ez,
40 const double bx, const double by, const double bz,
41 double& vx, double& vy, double& vz);
42 bool HoleTownsend(const double ex, const double ey, const double ez,
43 const double bx, const double by, const double bz,
44 double& alpha);
45 bool HoleAttachment(const double ex, const double ey, const double ez,
46 const double bx, const double by, const double bz,
47 double& eta);
48
49 void SetLowFieldMobility(const double mue, const double muh);
53
56
57 void SetSaturationVelocity(const double vsate, const double vsath);
61
66
69
70
71 // Scaling
72 void SetDiffusionScaling(const double d){
73 diffScale = d;
74 }
75
76 // Microscopic transport properties
77 bool SetMaxElectronEnergy(const double e);
78 double GetMaxElectronEnergy() const { return m_eFinalG; }
79
80 bool Initialise();
81
82 // When enabled, the scattering rates table is written to file
83 // when loaded into memory.
84 void EnableScatteringRateOutput() { m_useCfOutput = true; }
85 void DisableScatteringRateOutput() { m_useCfOutput = false; }
86
87 void EnableNonParabolicity() { m_useNonParabolicity = true; }
88 void DisableNonParabolicity() { m_useNonParabolicity = false; }
89 void EnableFullBandDensityOfStates() { m_useFullBandDos = true; }
90 void DisableFullBandDensityOfStates() { m_useFullBandDos = false; }
91 void EnableAnisotropy() { m_useAnisotropy = true; }
92 void DisableAnisotropy() { m_useAnisotropy = false; }
93
94 // Get the electron energy (and its gradient)
95 // for a given (crystal) momentum
96 double GetElectronEnergy(const double px, const double py, const double pz,
97 double& vx, double& vy, double& vz,
98 const int band = 0);
99 // Get the electron (crystal) momentum for a given kinetic energy
100 void GetElectronMomentum(const double e, double& px, double& py, double& pz,
101 int& band);
102
103 // Get the null-collision rate [ns-1]
104 double GetElectronNullCollisionRate(const int band);
105 // Get the (real) collision rate [ns-1] at a given electron energy
106 double GetElectronCollisionRate(const double e, const int band);
107 // Sample the collision type
108 bool GetElectronCollision(const double e, int& type, int& level, double& e1,
109 double& dx, double& dy, double& dz, int& nion,
110 int& ndxc, int& band);
111 unsigned int GetNumberOfIonisationProducts() const {
112 return m_ionProducts.size();
113 }
114 bool GetIonisationProduct(const unsigned int i, int& type,
115 double& energy) const;
116
117 // Density of states
118 double GetConductionBandDensityOfStates(const double e, const int band = 0);
119 double GetValenceBandDensityOfStates(const double e, const int band = -1);
120
121 // Reset the collision counters
123 // Get the total number of electron collisions
124 unsigned int GetNumberOfElectronCollisions() const;
125 // Get number of scattering rate terms
126 unsigned int GetNumberOfLevels() const;
127 // Get number of collisions for a specific level
128 unsigned int GetNumberOfElectronCollisions(const unsigned int level) const;
129
130 unsigned int GetNumberOfElectronBands() const;
131 int GetElectronBandPopulation(const int band);
132
133 bool GetOpticalDataRange(double& emin, double& emax,
134 const unsigned int i = 0);
135 bool GetDielectricFunction(const double e, double& eps1, double& eps2,
136 const unsigned int i = 0);
137
138 void ComputeSecondaries(const double e0, double& ee, double& eh);
139
140 private:
141 static const int LatticeMobilityModelSentaurus = 0;
142 static const int LatticeMobilityModelMinimos = 1;
143 static const int LatticeMobilityModelReggiani = 2;
144 static const int DopingMobilityModelMinimos = 0;
145 static const int DopingMobilityModelMasetti = 1;
146 static const int SaturationVelocityModelMinimos = 0;
147 static const int SaturationVelocityModelCanali = 1;
148 static const int SaturationVelocityModelReggiani = 2;
149 static const int HighFieldMobilityModelMinimos = 0;
150 static const int HighFieldMobilityModelCanali = 1;
151 static const int HighFieldMobilityModelReggiani = 2;
152 static const int HighFieldMobilityModelConstant = 3;
153 static const int ImpactIonisationModelVanOverstraeten = 0;
154 static const int ImpactIonisationModelGrant = 1;
155
156 // DiffusionScale
157 double diffScale;
158
159 double m_bandGap;
160 // Doping
161 char m_dopingType;
162 double m_dopingConcentration;
163
164 // Effective masses
165 // X valleys
166 double m_mLongX, m_mTransX;
167 // L valleys
168 double m_mLongL, m_mTransL;
169 // Non-parabolicity parameters [1/eV]
170 double m_alphaX, m_alphaL;
171 // Lattice mobility
172 double m_eLatticeMobility, m_hLatticeMobility;
173 // Low-field mobility
174 double m_eMobility, m_hMobility;
175 // High-field mobility parameters
176 double m_eBetaCanali, m_hBetaCanali;
177 double m_eBetaCanaliInv, m_hBetaCanaliInv;
178 // Saturation velocity
179 double m_eSatVel, m_hSatVel;
180 // Hall factor
181 double m_eHallFactor, m_hHallFactor;
182
183 // Trapping parameters
184 double m_eTrapCs, m_hTrapCs;
185 double m_eTrapDensity, m_hTrapDensity;
186 double m_eTrapTime, m_hTrapTime;
187 int m_trappingModel;
188
189 // Impact ionisation parameters
190 double m_eImpactA0, m_eImpactA1, m_eImpactA2;
191 double m_eImpactB0, m_eImpactB1, m_eImpactB2;
192 double m_hImpactA0, m_hImpactA1;
193 double m_hImpactB0, m_hImpactB1;
194
195 // Models
196 bool m_hasUserMobility;
197 bool m_hasUserSaturationVelocity;
198 int m_latticeMobilityModel;
199 int m_dopingMobilityModel;
200 int m_saturationVelocityModel;
201 int m_highFieldMobilityModel;
202 int m_impactIonisationModel;
203
204 // Options
205 bool m_useCfOutput;
206 bool m_useNonParabolicity;
207 bool m_useFullBandDos;
208 bool m_useAnisotropy;
209
210 // Energy range of scattering rates
211 double m_eFinalXL, m_eStepXL;
212 double m_eFinalG, m_eStepG;
213 double m_eFinalV, m_eStepV;
214 static const int nEnergyStepsXL = 2000;
215 static const int nEnergyStepsG = 2000;
216 static const int nEnergyStepsV = 2000;
217
218 // Number of scattering terms
219 int m_nLevelsX, m_nLevelsL, m_nLevelsG;
220 int m_nLevelsV;
221 // Number of valleys
222 int m_nValleysX, m_nValleysL;
223 // Energy offset
224 double m_eMinL, m_eMinG;
225 int m_ieMinL, m_ieMinG;
226
227 // Electron scattering rates
228 double m_cfNullElectronsX;
229 double m_cfNullElectronsL;
230 double m_cfNullElectronsG;
231 std::vector<double> m_cfTotElectronsX;
232 std::vector<double> m_cfTotElectronsL;
233 std::vector<double> m_cfTotElectronsG;
234 std::vector<std::vector<double> > m_cfElectronsX;
235 std::vector<std::vector<double> > m_cfElectronsL;
236 std::vector<std::vector<double> > m_cfElectronsG;
237 std::vector<double> m_energyLossElectronsX;
238 std::vector<double> m_energyLossElectronsL;
239 std::vector<double> m_energyLossElectronsG;
240 // Cross-section type
241 std::vector<int> m_scatTypeElectronsX;
242 std::vector<int> m_scatTypeElectronsL;
243 std::vector<int> m_scatTypeElectronsG;
244
245 // Hole scattering rates
246 double m_cfNullHoles;
247 std::vector<double> m_cfTotHoles;
248 std::vector<std::vector<double> > m_cfHoles;
249 std::vector<double> m_energyLossHoles;
250 // Cross-section type
251 std::vector<int> m_scatTypeHoles;
252
253 // Collision counters
254 unsigned int m_nCollElectronAcoustic;
255 unsigned int m_nCollElectronOptical;
256 unsigned int m_nCollElectronIntervalley;
257 unsigned int m_nCollElectronImpurity;
258 unsigned int m_nCollElectronIonisation;
259 std::vector<unsigned int> m_nCollElectronDetailed;
260 std::vector<unsigned int> m_nCollElectronBand;
261
262 struct ionProd {
263 int type;
264 double energy;
265 };
266 std::vector<ionProd> m_ionProducts;
267
268 // Density of states tables
269 double m_eStepDos;
270 std::vector<double> m_fbDosValence;
271 std::vector<double> m_fbDosConduction;
272 double m_fbDosMaxV, m_fbDosMaxC;
273
274 // Optical data
275 std::string m_opticalDataFile;
276 struct opticalData {
277 // Energy [eV]
278 double energy;
279 // Dielectric function
280 double eps1, eps2;
281 };
282 std::vector<opticalData> m_opticalDataTable;
283
284 bool UpdateTransportParameters();
285 void UpdateLatticeMobilityMinimos();
286 void UpdateLatticeMobilitySentaurus();
287 void UpdateLatticeMobilityReggiani();
288
289 void UpdateDopingMobilityMinimos();
290 void UpdateDopingMobilityMasetti();
291
292 void UpdateSaturationVelocityMinimos();
293 void UpdateSaturationVelocityCanali();
294 void UpdateSaturationVelocityReggiani();
295
296 void UpdateHighFieldMobilityCanali();
297
298 void UpdateImpactIonisationVanOverstraetenDeMan();
299 void UpdateImpactIonisationGrant();
300
301 bool ElectronMobilityMinimos(const double e, double& mu) const;
302 bool ElectronMobilityCanali(const double e, double& mu) const;
303 bool ElectronMobilityReggiani(const double e, double& mu) const;
304 bool ElectronImpactIonisationVanOverstraetenDeMan(const double e,
305 double& alpha) const;
306 bool ElectronImpactIonisationGrant(const double e, double& alpha) const;
307 bool HoleMobilityMinimos(const double e, double& mu) const;
308 bool HoleMobilityCanali(const double e, double& mu) const;
309 bool HoleMobilityReggiani(const double e, double& mu) const;
310 bool HoleImpactIonisationVanOverstraetenDeMan(const double e,
311 double& alpha) const;
312 bool HoleImpactIonisationGrant(const double e, double& alpha) const;
313
314 bool LoadOpticalData(const std::string& filename);
315
316 bool ElectronScatteringRates();
317 bool ElectronAcousticScatteringRates();
318 bool ElectronOpticalScatteringRates();
319 bool ElectronIntervalleyScatteringRatesXX();
320 bool ElectronIntervalleyScatteringRatesXL();
321 bool ElectronIntervalleyScatteringRatesLL();
322 bool ElectronIntervalleyScatteringRatesXGLG();
323 bool ElectronIonisationRatesXL();
324 bool ElectronIonisationRatesG();
325 bool ElectronImpurityScatteringRates();
326
327 bool HoleScatteringRates();
328 bool HoleAcousticScatteringRates();
329 bool HoleOpticalScatteringRates();
330 bool HoleIonisationRates();
331
332 // void ComputeSecondaries(const double e0, double& ee, double& eh);
333 void InitialiseDensityOfStates();
334};
335}
336
337#endif
Solid crystalline silicon
Definition: MediumSilicon.hh:9
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 ~MediumSilicon()
Destructor.
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 SetDoping(const char type, const double c)
Doping concentration [cm-3] and type ('i', 'n', 'p')
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)
unsigned int GetNumberOfLevels() const
void SetSaturationVelocityModelReggiani()
double GetValenceBandDensityOfStates(const double e, const int band=-1)
void SetTrapDensity(const double n)
double GetElectronNullCollisionRate(const int band)
bool GetDielectricFunction(const double e, double &eps1, double &eps2, const unsigned int i=0)
bool IsSemiconductor() const
MediumSilicon()
Constructor.
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)
double GetConductionBandDensityOfStates(const double e, const int band=0)
bool SetMaxElectronEnergy(const double e)
void SetHighFieldMobilityModelConstant()
bool GetOpticalDataRange(double &emin, double &emax, const unsigned int i=0)
void SetImpactIonisationModelVanOverstraetenDeMan()
void SetTrappingTime(const double etau, const double htau)
bool GetIonisationProduct(const unsigned int i, int &type, double &energy) const
void SetSaturationVelocityModelMinimos()
void SetLowFieldMobility(const double mue, const double muh)
double GetMaxElectronEnergy() const
unsigned int GetNumberOfIonisationProducts() const
double GetElectronCollisionRate(const double e, const int band)
unsigned int GetNumberOfElectronCollisions() const
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)
void SetTrapCrossSection(const double ecs, const double hcs)
Trapping cross-sections for electrons and holes.
void SetDiffusionScaling(const double d)
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)
unsigned int GetNumberOfElectronBands() const
Abstract base class for media.
Definition: Medium.hh:11