3#ifndef G_MEDIUM_SILICON_H
4#define G_MEDIUM_SILICON_H
21 void SetDoping(
const char& type,
const double& c);
22 void GetDoping(
char& type,
double& c)
const;
31 const double bx,
const double by,
const double bz,
32 double& vx,
double& vy,
double& vz);
34 const double bx,
const double by,
const double bz,
37 const double bx,
const double by,
const double bz,
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,
46 bool HoleAttachment(
const double ex,
const double ey,
const double ez,
47 const double bx,
const double by,
const double bz,
92 double& vx,
double& vy,
double& vz,
104 double& dx,
double& dy,
double& dz,
int& nion,
105 int& ndxc,
int& band);
126 const unsigned int& i = 0);
128 const unsigned int& i = 0);
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;
151 double m_dopingConcentration;
155 double mLongX, mTransX;
157 double mLongL, mTransL;
159 double alphaX, alphaL;
161 double eLatticeMobility, hLatticeMobility;
163 double eMobility, hMobility;
165 double eBetaCanali, hBetaCanali;
166 double eBetaCanaliInv, hBetaCanaliInv;
168 double eSatVel, hSatVel;
170 double eHallFactor, hHallFactor;
173 double eTrapCs, hTrapCs;
174 double eTrapDensity, hTrapDensity;
175 double eTrapTime, hTrapTime;
179 double eImpactA0, eImpactA1, eImpactA2;
180 double eImpactB0, eImpactB1, eImpactB2;
181 double hImpactA0, hImpactA1;
182 double hImpactB0, hImpactB1;
185 bool m_hasUserMobility;
186 bool m_hasUserSaturationVelocity;
187 int latticeMobilityModel;
188 int dopingMobilityModel;
189 int saturationVelocityModel;
190 int highFieldMobilityModel;
191 int impactIonisationModel;
195 bool useNonParabolicity;
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;
208 int nLevelsX, nLevelsL, nLevelsG;
211 int nValleysX, nValleysL;
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;
228 std::vector<int> scatTypeElectronsX;
229 std::vector<int> scatTypeElectronsL;
230 std::vector<int> scatTypeElectronsG;
234 std::vector<double> cfTotHoles;
235 std::vector<std::vector<double> > cfHoles;
236 std::vector<double> energyLossHoles;
238 std::vector<int> scatTypeHoles;
241 int nCollElectronAcoustic, nCollElectronOptical;
242 int nCollElectronIntervalley;
243 int nCollElectronImpurity;
244 int nCollElectronIonisation;
245 std::vector<int> nCollElectronDetailed;
246 std::vector<int> nCollElectronBand;
248 int nIonisationProducts;
253 std::vector<ionProd> ionProducts;
257 int nFbDosEntriesValence;
258 int nFbDosEntriesConduction;
259 std::vector<double> fbDosValence;
260 std::vector<double> fbDosConduction;
261 double fbDosMaxV, fbDosMaxC;
264 bool m_hasOpticalData;
265 std::string opticalDataFile;
272 std::vector<opticalData> opticalDataTable;
274 bool UpdateTransportParameters();
275 void UpdateLatticeMobilityMinimos();
276 void UpdateLatticeMobilitySentaurus();
277 void UpdateLatticeMobilityReggiani();
279 void UpdateDopingMobilityMinimos();
280 void UpdateDopingMobilityMasetti();
282 void UpdateSaturationVelocityMinimos();
283 void UpdateSaturationVelocityCanali();
284 void UpdateSaturationVelocityReggiani();
286 void UpdateHighFieldMobilityCanali();
288 void UpdateImpactIonisationVanOverstraetenDeMan();
289 void UpdateImpactIonisationGrant();
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;
304 bool LoadOpticalData(
const std::string filename);
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();
317 bool HoleScatteringRates();
318 bool HoleAcousticScatteringRates();
319 bool HoleOpticalScatteringRates();
320 bool HoleIonisationRates();
324 void InitialiseDensityOfStates();
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 SetLatticeMobilityModelReggiani()
void EnableNonParabolicity()
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)
void DisableFullBandDensityOfStates()
void SetDopingMobilityModelMasetti()
double GetElectronNullCollisionRate(const int band)
bool GetDielectricFunction(const double &e, double &eps1, double &eps2, const unsigned int &i=0)
void SetHighFieldMobilityModelMinimos()
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)
void ResetCollisionCounters()
int GetNumberOfElectronBands()
void DisableScatteringRateOutput()
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 SetLatticeMobilityModelSentaurus()
void SetSaturationVelocityModelCanali()
void SetHighFieldMobilityModelCanali()
void SetImpactIonisationModelVanOverstraetenDeMan()
bool GetOpticalDataRange(double &emin, double &emax, const unsigned int &i=0)
void SetTrapCrossSection(const double &ecs, const double &hcs)
void SetSaturationVelocityModelMinimos()
void EnableFullBandDensityOfStates()
void SetLowFieldMobility(const double mue, const double muh)
int GetNumberOfIonisationProducts()
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 DisableNonParabolicity()
void SetLatticeMobilityModelMinimos()
void GetElectronMomentum(const double e, double &px, double &py, double &pz, int &band)
int GetNumberOfElectronCollisions() const
int GetElectronBandPopulation(const int band)
void SetImpactIonisationModelGrant()
bool GetElectronCollision(const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, int &nion, int &ndxc, int &band)
void SetDopingMobilityModelMinimos()
void EnableScatteringRateOutput()