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

Solid crystalline silicon More...

#include <MediumSilicon.hh>

+ Inheritance diagram for Garfield::MediumSilicon:

Public Member Functions

 MediumSilicon ()
 Constructor.
 
virtual ~MediumSilicon ()
 Destructor.
 
bool IsSemiconductor () const override
 Is this medium a semiconductor?
 
void SetDoping (const char type, const double c)
 Set doping concentration [cm-3] and type ('i', 'n', 'p').
 
void GetDoping (char &type, double &c) const
 Retrieve doping concentration.
 
void SetTrapCrossSection (const double ecs, const double hcs)
 Trapping cross-sections for electrons and holes.
 
void SetTrapDensity (const double n)
 Trap density [cm-3], by default set to zero.
 
void SetTrappingTime (const double etau, const double htau)
 Set time constant for trapping of electrons and holes [ns].
 
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) override
 Drift velocity [cm / ns].
 
bool ElectronTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha) override
 Ionisation coefficient [cm-1].
 
bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta) override
 Attachment coefficient [cm-1].
 
double ElectronMobility () override
 Low-field mobility [cm2 V-1 ns-1].
 
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) override
 Drift velocity [cm / ns].
 
bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha) override
 Ionisation coefficient [cm-1].
 
bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta) override
 Attachment coefficient [cm-1].
 
double HoleMobility () override
 Low-field mobility [cm2 V-1 ns-1].
 
void SetLowFieldMobility (const double mue, const double muh)
 Specify the low field values of the electron and hole mobilities.
 
void SetLatticeMobilityModelMinimos ()
 Calculate the lattice mobility using the Minimos model.
 
void SetLatticeMobilityModelSentaurus ()
 Calculate the lattice mobility using the Sentaurus model (default).
 
void SetLatticeMobilityModelReggiani ()
 Calculate the lattice mobility using the Reggiani model.
 
void SetDopingMobilityModelMinimos ()
 Use the Minimos model for the doping-dependence of the mobility.
 
void SetDopingMobilityModelMasetti ()
 Use the Masetti model for the doping-dependence of the mobility (default).
 
void SetSaturationVelocity (const double vsate, const double vsath)
 Specify the saturation velocities of electrons and holes.
 
void SetSaturationVelocityModelMinimos ()
 Calculate the saturation velocities using the Minimos model.
 
void SetSaturationVelocityModelCanali ()
 Calculate the saturation velocities using the Canali model (default).
 
void SetSaturationVelocityModelReggiani ()
 Calculate the saturation velocities using the Reggiani model.
 
void SetHighFieldMobilityModelMinimos ()
 Parameterize the high-field mobility using the Minimos model.
 
void SetHighFieldMobilityModelCanali ()
 Parameterize the high-field mobility using the Canali model (default).
 
void SetHighFieldMobilityModelReggiani ()
 Parameterize the high-field mobility using the Reggiani model.
 
void SetHighFieldMobilityModelConstant ()
 Make the velocity proportional to the electric field (no saturation).
 
void SetImpactIonisationModelVanOverstraetenDeMan ()
 Calculate α using the van Overstraeten-de Man model (default).
 
void SetImpactIonisationModelGrant ()
 Calculate α using the Grant model.
 
void SetImpactIonisationModelMassey ()
 Calculate α using the Massey model.
 
void SetDiffusionScaling (const double d)
 Apply a scaling factor to the diffusion coefficients.
 
bool SetMaxElectronEnergy (const double e)
 
double GetMaxElectronEnergy () const
 
bool Initialise ()
 
void EnableScatteringRateOutput (const bool on=true)
 
void EnableNonParabolicity (const bool on=true)
 
void EnableFullBandDensityOfStates (const bool on=true)
 
void EnableAnisotropy (const bool on=true)
 
double GetElectronEnergy (const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0) override
 Dispersion relation (energy vs. wave vector)
 
void GetElectronMomentum (const double e, double &px, double &py, double &pz, int &band) override
 
double GetElectronNullCollisionRate (const int band) override
 Null-collision rate [ns-1].
 
double GetElectronCollisionRate (const double e, const int band) override
 Collision rate [ns-1] for given electron energy.
 
bool GetElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< std::pair< int, double > > &secondaries, int &ndxc, int &band) override
 Sample the collision type. Update energy and direction vector.
 
double GetConductionBandDensityOfStates (const double e, const int band=0)
 
double GetValenceBandDensityOfStates (const double e, const int band=-1)
 
void ResetCollisionCounters ()
 
unsigned int GetNumberOfElectronCollisions () const
 
unsigned int GetNumberOfLevels () const
 
unsigned int GetNumberOfElectronCollisions (const unsigned int level) const
 
unsigned int GetNumberOfElectronBands () const
 
int GetElectronBandPopulation (const int band)
 
bool GetOpticalDataRange (double &emin, double &emax, const unsigned int i=0) override
 Get the energy range [eV] of the available optical data.
 
bool GetDielectricFunction (const double e, double &eps1, double &eps2, const unsigned int i=0) override
 Get the complex dielectric function at a given energy.
 
void ComputeSecondaries (const double e0, double &ee, double &eh)
 
- Public Member Functions inherited from Garfield::Medium
 Medium ()
 Constructor.
 
virtual ~Medium ()
 Destructor.
 
int GetId () const
 Return the id number of the class instance.
 
const std::string & GetName () const
 Get the medium name/identifier.
 
virtual bool IsGas () const
 Is this medium a gas?
 
virtual bool IsSemiconductor () const
 Is this medium a semiconductor?
 
virtual bool IsConductor () const
 Is this medium a conductor?
 
void SetTemperature (const double t)
 Set the temperature [K].
 
double GetTemperature () const
 Get the temperature [K].
 
void SetPressure (const double p)
 
double GetPressure () const
 
void SetDielectricConstant (const double eps)
 Set the relative static dielectric constant.
 
double GetDielectricConstant () const
 Get the relative static dielectric constant.
 
unsigned int GetNumberOfComponents () const
 Get number of components of the medium.
 
virtual void GetComponent (const unsigned int i, std::string &label, double &f)
 Get the name and fraction of a given component.
 
virtual void SetAtomicNumber (const double z)
 Set the effective atomic number.
 
virtual double GetAtomicNumber () const
 Get the effective atomic number.
 
virtual void SetAtomicWeight (const double a)
 Set the effective atomic weight.
 
virtual double GetAtomicWeight () const
 Get the effective atomic weight.
 
virtual void SetNumberDensity (const double n)
 Set the number density [cm-3].
 
virtual double GetNumberDensity () const
 Get the number density [cm-3].
 
virtual void SetMassDensity (const double rho)
 Set the mass density [g/cm3].
 
virtual double GetMassDensity () const
 Get the mass density [g/cm3].
 
virtual void EnableDrift (const bool on=true)
 Switch electron/ion/hole on/off.
 
virtual void EnablePrimaryIonisation (const bool on=true)
 Make the medium ionisable or non-ionisable.
 
bool IsDriftable () const
 Is charge carrier transport enabled in this medium?
 
bool IsMicroscopic () const
 Does the medium have electron scattering rates?
 
bool IsIonisable () const
 Is charge deposition by charged particles/photon enabled in this medium?
 
void SetW (const double w)
 Set the W value (average energy to produce an electron/ion or e/h pair).
 
double GetW ()
 Get the W value.
 
void SetFanoFactor (const double f)
 Set the Fano factor.
 
double GetFanoFactor ()
 Get the Fano factor.
 
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)
 Drift velocity [cm / ns].
 
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)
 Longitudinal and transverse diffusion coefficients [cm1/2].
 
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)
 Ionisation coefficient [cm-1].
 
virtual bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 Attachment coefficient [cm-1].
 
virtual bool ElectronLorentzAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
 Lorentz angle.
 
virtual double ElectronMobility ()
 Low-field mobility [cm2 V-1 ns-1].
 
virtual double GetElectronEnergy (const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
 Dispersion relation (energy vs. wave vector)
 
virtual void GetElectronMomentum (const double e, double &px, double &py, double &pz, int &band)
 
virtual double GetElectronNullCollisionRate (const int band=0)
 Null-collision rate [ns-1].
 
virtual double GetElectronCollisionRate (const double e, const int band=0)
 Collision rate [ns-1] for given electron energy.
 
virtual bool GetElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< std::pair< int, double > > &secondaries, int &ndxc, int &band)
 Sample the collision type. Update energy and direction vector.
 
virtual unsigned int GetNumberOfDeexcitationProducts () const
 
virtual bool GetDeexcitationProduct (const unsigned int i, double &t, double &s, int &type, double &energy) const
 
virtual bool HoleVelocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
 Drift velocity [cm / ns].
 
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)
 Longitudinal and transverse diffusion coefficients [cm1/2].
 
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])
 Diffusion tensor.
 
virtual bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 Ionisation coefficient [cm-1].
 
virtual bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 Attachment coefficient [cm-1].
 
virtual double HoleMobility ()
 Low-field mobility [cm2 V-1 ns-1].
 
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)
 Drift velocity [cm / ns].
 
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)
 Longitudinal and transverse diffusion coefficients [cm1/2].
 
virtual bool IonDissociation (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &diss)
 Dissociation coefficient.
 
virtual double IonMobility ()
 Low-field mobility [cm2 V-1 ns-1].
 
void SetFieldGrid (double emin, double emax, const size_t ne, bool logE, double bmin=0., double bmax=0., const size_t nb=1, double amin=HalfPi, double amax=HalfPi, const size_t na=1)
 Set the range of fields to be covered by the transport tables.
 
void SetFieldGrid (const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles)
 Set the fields and E-B angles to be used in the transport tables.
 
void GetFieldGrid (std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles)
 Get the fields and E-B angles used in the transport tables.
 
bool SetElectronVelocityE (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along E.
 
bool GetElectronVelocityE (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along E.
 
bool SetElectronVelocityExB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along ExB.
 
bool GetElectronVelocityExB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along ExB.
 
bool SetElectronVelocityB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along Btrans.
 
bool GetElectronVelocityB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along Btrans.
 
bool SetElectronLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetElectronLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetElectronTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetElectronTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetElectronTownsend (const size_t ie, const size_t ib, const size_t ia, const double alpha)
 Set an entry in the table of Townsend coefficients.
 
bool GetElectronTownsend (const size_t ie, const size_t ib, const size_t ia, double &alpha)
 Get an entry in the table of Townsend coefficients.
 
bool SetElectronAttachment (const size_t ie, const size_t ib, const size_t ia, const double eta)
 Set an entry in the table of attachment coefficients.
 
bool GetElectronAttachment (const size_t ie, const size_t ib, const size_t ia, double &eta)
 Get an entry in the table of attachment coefficients.
 
bool SetElectronLorentzAngle (const size_t ie, const size_t ib, const size_t ia, const double lor)
 Set an entry in the table of Lorentz angles.
 
bool GetElectronLorentzAngle (const size_t ie, const size_t ib, const size_t ia, double &lor)
 Get an entry in the table of Lorentz angles.
 
bool SetHoleVelocityE (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along E.
 
bool GetHoleVelocityE (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along E.
 
bool SetHoleVelocityExB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along ExB.
 
bool GetHoleVelocityExB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along ExB.
 
bool SetHoleVelocityB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along Btrans.
 
bool GetHoleVelocityB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along Btrans.
 
bool SetHoleLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetHoleLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetHoleTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetHoleTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetHoleTownsend (const size_t ie, const size_t ib, const size_t ia, const double alpha)
 Set an entry in the table of Townsend coefficients.
 
bool GetHoleTownsend (const size_t ie, const size_t ib, const size_t ia, double &alpha)
 Get an entry in the table of Townsend coefficients.
 
bool SetHoleAttachment (const size_t ie, const size_t ib, const size_t ia, const double eta)
 Set an entry in the table of attachment coefficients.
 
bool GetHoleAttachment (const size_t ie, const size_t ib, const size_t ia, double &eta)
 Get an entry in the table of attachment coefficients.
 
bool SetIonMobility (const std::vector< double > &fields, const std::vector< double > &mobilities)
 
bool SetIonMobility (const size_t ie, const size_t ib, const size_t ia, const double mu)
 Set an entry in the table of ion mobilities.
 
bool GetIonMobility (const size_t ie, const size_t ib, const size_t ia, double &mu)
 Get an entry in the table of ion mobilities.
 
bool SetIonLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetIonLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetIonTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetIonTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetIonDissociation (const size_t ie, const size_t ib, const size_t ia, const double diss)
 Set an entry in the table of dissociation coefficients.
 
bool GetIonDissociation (const size_t ie, const size_t ib, const size_t ia, double &diss)
 Get an entry in the table of dissociation coefficients.
 
virtual void ResetTables ()
 Reset all tables of transport parameters.
 
void ResetElectronVelocity ()
 
void ResetElectronDiffusion ()
 
void ResetElectronTownsend ()
 
void ResetElectronAttachment ()
 
void ResetElectronLorentzAngle ()
 
void ResetHoleVelocity ()
 
void ResetHoleDiffusion ()
 
void ResetHoleTownsend ()
 
void ResetHoleAttachment ()
 
void ResetIonMobility ()
 
void ResetIonDiffusion ()
 
void ResetIonDissociation ()
 
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)
 Set the degree of polynomial interpolation (usually 2).
 
void SetInterpolationMethodDiffusion (const unsigned int intrp)
 
void SetInterpolationMethodTownsend (const unsigned int intrp)
 
void SetInterpolationMethodAttachment (const unsigned int intrp)
 
void SetInterpolationMethodIonMobility (const unsigned int intrp)
 
void SetInterpolationMethodIonDissociation (const unsigned int intrp)
 
virtual double ScaleElectricField (const double e) const
 
virtual double UnScaleElectricField (const double e) const
 
virtual double ScaleVelocity (const double v) const
 
virtual double ScaleDiffusion (const double d) const
 
virtual double ScaleDiffusionTensor (const double d) const
 
virtual double ScaleTownsend (const double alpha) const
 
virtual double ScaleAttachment (const double eta) const
 
virtual double ScaleLorentzAngle (const double lor) const
 
virtual double ScaleDissociation (const double diss) const
 
virtual bool GetOpticalDataRange (double &emin, double &emax, const unsigned int i=0)
 Get the energy range [eV] of the available optical data.
 
virtual bool GetDielectricFunction (const double e, double &eps1, double &eps2, const unsigned int i=0)
 Get the complex dielectric function at a given energy.
 
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 ()
 Switch on/off debugging messages.
 
void DisableDebugging ()
 

Additional Inherited Members

- Protected Member Functions inherited from Garfield::Medium
bool Velocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &velE, const std::vector< std::vector< std::vector< double > > > &velB, const std::vector< std::vector< std::vector< double > > > &velX, const double q, double &vx, double &vy, double &vz) const
 
bool Diffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &difL, const std::vector< std::vector< std::vector< double > > > &difT, double &dl, double &dt) const
 
bool Diffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< std::vector< double > > > > &diff, double cov[3][3]) const
 
bool Alpha (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &tab, unsigned int intp, const unsigned int thr, const std::pair< unsigned int, unsigned int > &extr, double &alpha) const
 
double GetAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const double e, const double b) const
 
bool Interpolate (const double e, const double b, const double a, const std::vector< std::vector< std::vector< double > > > &table, double &y, const unsigned int intp, const std::pair< unsigned int, unsigned int > &extr) const
 
double Interpolate1D (const double e, const std::vector< double > &table, const std::vector< double > &fields, const unsigned int intpMeth, const std::pair< unsigned int, unsigned int > &extr) const
 
bool SetEntry (const size_t i, const size_t j, const size_t k, const std::string &fcn, std::vector< std::vector< std::vector< double > > > &tab, const double val)
 
bool GetEntry (const size_t i, const size_t j, const size_t k, const std::string &fcn, const std::vector< std::vector< std::vector< double > > > &tab, double &val) const
 
void SetExtrapolationMethod (const std::string &low, const std::string &high, std::pair< unsigned int, unsigned int > &extr, const std::string &fcn)
 
bool GetExtrapolationIndex (std::string str, unsigned int &nb) const
 
size_t SetThreshold (const std::vector< std::vector< std::vector< double > > > &tab) const
 
void Clone (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 std::pair< unsigned int, unsigned int > &extr, const double init, const std::string &label)
 
void Clone (std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const size_t n, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int intp, const std::pair< unsigned int, unsigned int > &extr, const double init, const std::string &label)
 
void Init (const size_t nE, const size_t nB, const size_t nA, std::vector< std::vector< std::vector< double > > > &tab, const double val)
 
void Init (const size_t nE, const size_t nB, const size_t nA, const size_t nT, std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const double val)
 
- Protected Attributes inherited from Garfield::Medium
std::string m_className = "Medium"
 
int m_id
 
unsigned int m_nComponents = 1
 
std::string m_name = ""
 
double m_temperature = 293.15
 
double m_pressure = 760.
 
double m_epsilon = 1.
 
double m_z = 1.
 
double m_a = 0.
 
double m_density = 0.
 
double m_w = 0.
 
double m_fano = 0.
 
bool m_driftable = false
 
bool m_microscopic = false
 
bool m_ionisable = false
 
bool m_isChanged = true
 
bool m_debug = false
 
bool m_tab2d = false
 
std::vector< double > m_eFields
 
std::vector< double > m_bFields
 
std::vector< double > m_bAngles
 
std::vector< std::vector< std::vector< double > > > m_eVelE
 
std::vector< std::vector< std::vector< double > > > m_eVelX
 
std::vector< std::vector< std::vector< double > > > m_eVelB
 
std::vector< std::vector< std::vector< double > > > m_eDifL
 
std::vector< std::vector< std::vector< double > > > m_eDifT
 
std::vector< std::vector< std::vector< double > > > m_eAlp
 
std::vector< std::vector< std::vector< double > > > m_eAtt
 
std::vector< std::vector< std::vector< double > > > m_eLor
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
 
std::vector< std::vector< std::vector< double > > > m_hVelE
 
std::vector< std::vector< std::vector< double > > > m_hVelX
 
std::vector< std::vector< std::vector< double > > > m_hVelB
 
std::vector< std::vector< std::vector< double > > > m_hDifL
 
std::vector< std::vector< std::vector< double > > > m_hDifT
 
std::vector< std::vector< std::vector< double > > > m_hAlp
 
std::vector< std::vector< std::vector< double > > > m_hAtt
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_hDifM
 
std::vector< std::vector< std::vector< double > > > m_iMob
 
std::vector< std::vector< std::vector< double > > > m_iDifL
 
std::vector< std::vector< std::vector< double > > > m_iDifT
 
std::vector< std::vector< std::vector< double > > > m_iDis
 
unsigned int m_eThrAlp = 0
 
unsigned int m_eThrAtt = 0
 
unsigned int m_hThrAlp = 0
 
unsigned int m_hThrAtt = 0
 
unsigned int m_iThrDis = 0
 
std::pair< unsigned int, unsigned int > m_extrVel = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrDif = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrAlp = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrAtt = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrLor = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrMob = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrDis = {0, 1}
 
unsigned int m_intpVel = 2
 
unsigned int m_intpDif = 2
 
unsigned int m_intpAlp = 2
 
unsigned int m_intpAtt = 2
 
unsigned int m_intpLor = 2
 
unsigned int m_intpMob = 2
 
unsigned int m_intpDis = 2
 
- Static Protected Attributes inherited from Garfield::Medium
static int m_idCounter = -1
 

Detailed Description

Solid crystalline silicon

Definition at line 13 of file MediumSilicon.hh.

Constructor & Destructor Documentation

◆ MediumSilicon()

Garfield::MediumSilicon::MediumSilicon ( )

Constructor.

Definition at line 15 of file MediumSilicon.cc.

16 : Medium(),
17 m_eStepXL(m_eFinalXL / nEnergyStepsXL),
18 m_eStepG(m_eFinalG / nEnergyStepsG),
19 m_eStepV(m_eFinalV / nEnergyStepsV) {
20 m_className = "MediumSilicon";
21 m_name = "Si";
22
23 SetTemperature(300.);
25 SetAtomicNumber(14.);
26 SetAtomicWeight(28.0855);
27 SetMassDensity(2.329);
28
31 m_microscopic = true;
32
33 m_w = 3.6;
34 m_fano = 0.11;
35
36 m_ieMinL = int(m_eMinL / m_eStepXL) + 1;
37 m_ieMinG = int(m_eMinG / m_eStepG) + 1;
38
39 // Load the density of states table.
40 InitialiseDensityOfStates();
41}
void SetTemperature(const double t)
Set the temperature [K].
Definition: Medium.cc:71
bool m_microscopic
Definition: Medium.hh:523
double m_fano
Definition: Medium.hh:519
virtual void EnableDrift(const bool on=true)
Switch electron/ion/hole on/off.
Definition: Medium.hh:67
std::string m_name
Definition: Medium.hh:502
virtual void EnablePrimaryIonisation(const bool on=true)
Make the medium ionisable or non-ionisable.
Definition: Medium.hh:69
virtual void SetAtomicNumber(const double z)
Set the effective atomic number.
Definition: Medium.cc:114
void SetDielectricConstant(const double eps)
Set the relative static dielectric constant.
Definition: Medium.cc:91
virtual void SetMassDensity(const double rho)
Set the mass density [g/cm3].
Definition: Medium.cc:144
Medium()
Constructor.
Definition: Medium.cc:60
virtual void SetAtomicWeight(const double a)
Set the effective atomic weight.
Definition: Medium.cc:124
std::string m_className
Definition: Medium.hh:493

◆ ~MediumSilicon()

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

Destructor.

Definition at line 18 of file MediumSilicon.hh.

18{}

Member Function Documentation

◆ ComputeSecondaries()

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

Definition at line 2921 of file MediumSilicon.cc.

2922 {
2923 const int nPointsValence = m_fbDosValence.size();
2924 const int nPointsConduction = m_fbDosConduction.size();
2925 const double widthValenceBand = m_eStepDos * nPointsValence;
2926 const double widthConductionBand = m_eStepDos * nPointsConduction;
2927
2928 bool ok = false;
2929 while (!ok) {
2930 // Sample a hole energy according to the valence band DOS.
2931 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2932 int ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2933 while (RndmUniform() > m_fbDosValence[ih] / m_fbDosMaxV) {
2934 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2935 ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2936 }
2937 // Sample an electron energy according to the conduction band DOS.
2938 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2939 int ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2940 while (RndmUniform() > m_fbDosConduction[ie] / m_fbDosMaxC) {
2941 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2942 ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2943 }
2944 // Calculate the energy of the primary electron.
2945 const double ep = e0 - m_bandGap - eh - ee;
2946 if (ep < Small) continue;
2947 if (ep > 5.) return;
2948 // Check if the primary electron energy is consistent with the DOS.
2949 int ip = std::min(int(ep / m_eStepDos), nPointsConduction - 1);
2950 if (RndmUniform() > m_fbDosConduction[ip] / m_fbDosMaxC) continue;
2951 ok = true;
2952 }
2953}
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition: Random.hh:17

Referenced by GetElectronCollision().

◆ 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 
)
overridevirtual

Attachment coefficient [cm-1].

Reimplemented from Garfield::Medium.

Definition at line 228 of file MediumSilicon.cc.

231 {
232 eta = 0.;
233 if (m_isChanged) {
234 if (!UpdateTransportParameters()) {
235 std::cerr << m_className << "::ElectronAttachment:\n"
236 << " Error calculating the transport parameters.\n";
237 return false;
238 }
239 m_isChanged = false;
240 }
241
242 if (!m_eAtt.empty()) {
243 // Interpolation in user table.
244 return Medium::ElectronAttachment(ex, ey, ez, bx, by, bz, eta);
245 }
246
247 switch (m_trappingModel) {
248 case 0:
249 eta = m_eTrapCs * m_eTrapDensity;
250 break;
251 case 1:
252 double vx, vy, vz;
253 ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
254 eta = m_eTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
255 if (eta > 0.) eta = -1. / eta;
256 break;
257 default:
258 std::cerr << m_className << "::ElectronAttachment: Unknown model. Bug!\n";
259 return false;
260 }
261
262 return true;
263}
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) override
Drift velocity [cm / ns].
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition: Medium.hh:547
virtual bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:416
bool m_isChanged
Definition: Medium.hh:527
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

◆ ElectronMobility()

double Garfield::MediumSilicon::ElectronMobility ( )
inlineoverridevirtual

Low-field mobility [cm2 V-1 ns-1].

Reimplemented from Garfield::Medium.

Definition at line 44 of file MediumSilicon.hh.

44{ return m_eMobility; }

◆ 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 
)
overridevirtual

Ionisation coefficient [cm-1].

Reimplemented from Garfield::Medium.

Definition at line 193 of file MediumSilicon.cc.

196 {
197 alpha = 0.;
198 if (m_isChanged) {
199 if (!UpdateTransportParameters()) {
200 std::cerr << m_className << "::ElectronTownsend:\n"
201 << " Error calculating the transport parameters.\n";
202 return false;
203 }
204 m_isChanged = false;
205 }
206
207 if (!m_eAlp.empty()) {
208 // Interpolation in user table.
209 return Medium::ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
210 }
211
212 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
213
214 switch (m_impactIonisationModel) {
215 case ImpactIonisation::VanOverstraeten:
216 return ElectronImpactIonisationVanOverstraetenDeMan(emag, alpha);
217 case ImpactIonisation::Grant:
218 return ElectronImpactIonisationGrant(emag, alpha);
219 case ImpactIonisation::Massey:
220 return ElectronImpactIonisationMassey(emag, alpha);
221 default:
222 std::cerr << m_className << "::ElectronTownsend: Unknown model. Bug!\n";
223 break;
224 }
225 return false;
226}
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition: Medium.hh:546
virtual bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:403

◆ 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 
)
overridevirtual

Drift velocity [cm / ns].

Reimplemented from Garfield::Medium.

Definition at line 135 of file MediumSilicon.cc.

138 {
139 vx = vy = vz = 0.;
140 if (m_isChanged) {
141 if (!UpdateTransportParameters()) {
142 std::cerr << m_className << "::ElectronVelocity:\n"
143 << " Error calculating the transport parameters.\n";
144 return false;
145 }
146 m_isChanged = false;
147 }
148
149 if (!m_eVelE.empty()) {
150 // Interpolation in user table.
151 return Medium::ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
152 }
153
154 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
155
156 // Calculate the mobility
157 double mu;
158 switch (m_highFieldMobilityModel) {
159 case HighFieldMobility::Minimos:
160 ElectronMobilityMinimos(emag, mu);
161 break;
162 case HighFieldMobility::Canali:
163 ElectronMobilityCanali(emag, mu);
164 break;
165 case HighFieldMobility::Reggiani:
166 ElectronMobilityReggiani(emag, mu);
167 break;
168 default:
169 mu = m_eMobility;
170 break;
171 }
172 mu = -mu;
173
174 const double b2 = bx * bx + by * by + bz * bz;
175 if (b2 < Small) {
176 vx = mu * ex;
177 vy = mu * ey;
178 vz = mu * ez;
179 } else {
180 // Hall mobility
181 const double muH = m_eHallFactor * mu;
182 const double mu2 = muH * muH;
183 const double f = mu / (1. + mu2 * b2);
184 const double eb = bx * ex + by * ey + bz * ez;
185 // Compute the drift velocity using the Langevin equation.
186 vx = f * (ex + muH * (ey * bz - ez * by) + mu2 * bx * eb);
187 vy = f * (ey + muH * (ez * bx - ex * bz) + mu2 * by * eb);
188 vz = f * (ez + muH * (ex * by - ey * bx) + mu2 * bz * eb);
189 }
190 return true;
191}
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)
Drift velocity [cm / ns].
Definition: Medium.cc:379
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition: Medium.hh:541

Referenced by ElectronAttachment().

◆ EnableAnisotropy()

void Garfield::MediumSilicon::EnableAnisotropy ( const bool  on = true)
inline

Definition at line 110 of file MediumSilicon.hh.

110{ m_anisotropic = on; }

◆ EnableFullBandDensityOfStates()

void Garfield::MediumSilicon::EnableFullBandDensityOfStates ( const bool  on = true)
inline

Definition at line 107 of file MediumSilicon.hh.

107 {
108 m_fullBandDos = on;
109 }

◆ EnableNonParabolicity()

void Garfield::MediumSilicon::EnableNonParabolicity ( const bool  on = true)
inline

Definition at line 106 of file MediumSilicon.hh.

106{ m_nonParabolic = on; }

◆ EnableScatteringRateOutput()

void Garfield::MediumSilicon::EnableScatteringRateOutput ( const bool  on = true)
inline

Definition at line 105 of file MediumSilicon.hh.

105{ m_cfOutput = on; }

◆ GetConductionBandDensityOfStates()

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

Definition at line 2786 of file MediumSilicon.cc.

2787 {
2788 if (band < 0) {
2789 int iE = int(e / m_eStepDos);
2790 const int nPoints = m_fbDosConduction.size();
2791 if (iE >= nPoints || iE < 0) {
2792 return 0.;
2793 } else if (iE == nPoints - 1) {
2794 return m_fbDosConduction[nPoints - 1];
2795 }
2796
2797 const double dos = m_fbDosConduction[iE] +
2798 (m_fbDosConduction[iE + 1] - m_fbDosConduction[iE]) *
2799 (e / m_eStepDos - iE);
2800 return dos * 1.e21;
2801
2802 } else if (band < m_nValleysX) {
2803 // X valleys
2804 if (e <= 0.) return 0.;
2805 // Density-of-states effective mass (cube)
2806 const double md3 = pow(ElectronMass, 3) * m_mLongX * m_mTransX * m_mTransX;
2807
2808 if (m_fullBandDos) {
2809 if (e < m_eMinL) {
2810 return GetConductionBandDensityOfStates(e, -1) / m_nValleysX;
2811 } else if (e < m_eMinG) {
2812 // Subtract the fraction of the full-band density of states
2813 // attributed to the L valleys.
2814 const double dosX =
2816 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL;
2817 return dosX / m_nValleysX;
2818 } else {
2819 // Subtract the fraction of the full-band density of states
2820 // attributed to the L valleys and the higher bands.
2821 const double dosX =
2823 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL -
2824 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
2825 if (dosX <= 0.) return 0.;
2826 return dosX / m_nValleysX;
2827 }
2828 } else if (m_nonParabolic) {
2829 const double alpha = m_alphaX;
2830 return sqrt(md3 * e * (1. + alpha * e) / 2.) * (1. + 2 * alpha * e) /
2831 (Pi2 * pow(HbarC, 3.));
2832 } else {
2833 return sqrt(md3 * e / 2.) / (Pi2 * pow(HbarC, 3.));
2834 }
2835 } else if (band < m_nValleysX + m_nValleysL) {
2836 // L valleys
2837 if (e <= m_eMinL) return 0.;
2838
2839 // Density-of-states effective mass (cube)
2840 const double md3 = pow(ElectronMass, 3) * m_mLongL * m_mTransL * m_mTransL;
2841 // Non-parabolicity parameter
2842 const double alpha = m_alphaL;
2843
2844 if (m_fullBandDos) {
2845 // Energy up to which the non-parabolic approximation is used.
2846 const double ej = m_eMinL + 0.5;
2847 if (e <= ej) {
2848 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
2849 (1. + 2 * alpha * (e - m_eMinL)) /
2850 (Sqrt2 * Pi2 * pow(HbarC, 3.));
2851 } else {
2852 // Fraction of full-band density of states attributed to L valleys
2853 double fL = sqrt(md3 * (ej - m_eMinL) * (1. + alpha * (ej - m_eMinL))) *
2854 (1. + 2 * alpha * (ej - m_eMinL)) /
2855 (Sqrt2 * Pi2 * pow(HbarC, 3.));
2856 fL = m_nValleysL * fL / GetConductionBandDensityOfStates(ej, -1);
2857
2858 double dosXL = GetConductionBandDensityOfStates(e, -1);
2859 if (e > m_eMinG) {
2860 dosXL -=
2861 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
2862 }
2863 if (dosXL <= 0.) return 0.;
2864 return fL * dosXL / 8.;
2865 }
2866 } else if (m_nonParabolic) {
2867 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
2868 (1. + 2 * alpha * (e - m_eMinL)) / (Sqrt2 * Pi2 * pow(HbarC, 3.));
2869 } else {
2870 return sqrt(md3 * (e - m_eMinL) / 2.) / (Pi2 * pow(HbarC, 3.));
2871 }
2872 } else if (band == m_nValleysX + m_nValleysL) {
2873 // Higher bands
2874 const double ej = 2.7;
2875 if (m_eMinG >= ej) {
2876 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2877 << " Cannot determine higher band density-of-states.\n"
2878 << " Program bug. Check offset energy!\n";
2879 return 0.;
2880 }
2881 if (e < m_eMinG) {
2882 return 0.;
2883 } else if (e < ej) {
2884 // Coexistence of XL and higher bands.
2885 const double dj = GetConductionBandDensityOfStates(ej, -1);
2886 // Assume linear increase of density-of-states.
2887 return dj * (e - m_eMinG) / (ej - m_eMinG);
2888 } else {
2890 }
2891 }
2892
2893 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2894 << " Band index (" << band << ") out of range.\n";
2895 return ElectronMass * sqrt(ElectronMass * e / 2.) / (Pi2 * pow(HbarC, 3.));
2896}
double GetConductionBandDensityOfStates(const double e, const int band=0)
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337

Referenced by GetConductionBandDensityOfStates(), and GetElectronMomentum().

◆ GetDielectricFunction()

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

Get the complex dielectric function at a given energy.

Reimplemented from Garfield::Medium.

Definition at line 1170 of file MediumSilicon.cc.

1171 {
1172 if (i != 0) {
1173 std::cerr << m_className + "::GetDielectricFunction: Index out of range.\n";
1174 return false;
1175 }
1176
1177 // Make sure the optical data table has been loaded.
1178 if (m_opticalDataEnergies.empty()) {
1179 if (!LoadOpticalData(m_opticalDataFile)) {
1180 std::cerr << m_className << "::GetDielectricFunction:\n";
1181 std::cerr << " Optical data table could not be loaded.\n";
1182 return false;
1183 }
1184 }
1185
1186 // Make sure the requested energy is within the range of the table.
1187 const double emin = m_opticalDataEnergies.front();
1188 const double emax = m_opticalDataEnergies.back();
1189 if (e < emin || e > emax) {
1190 std::cerr << m_className << "::GetDielectricFunction:\n"
1191 << " Requested energy (" << e << " eV) "
1192 << " is outside the range of the optical data table.\n"
1193 << " " << emin << " < E [eV] < " << emax << "\n";
1194 eps1 = eps2 = 0.;
1195 return false;
1196 }
1197
1198 // Locate the requested energy in the table.
1199 const auto begin = m_opticalDataEnergies.cbegin();
1200 const auto it1 = std::upper_bound(begin, m_opticalDataEnergies.cend(), e);
1201 if (it1 == begin) {
1202 eps1 = m_opticalDataEpsilon.front().first;
1203 eps2 = m_opticalDataEpsilon.front().second;
1204 return true;
1205 }
1206 const auto it0 = std::prev(it1);
1207
1208 // Interpolate the real part of dielectric function.
1209 const double x0 = *it0;
1210 const double x1 = *it1;
1211 const double lnx0 = log(*it0);
1212 const double lnx1 = log(*it1);
1213 const double lnx = log(e);
1214 const double y0 = m_opticalDataEpsilon[it0 - begin].first;
1215 const double y1 = m_opticalDataEpsilon[it1 - begin].first;
1216 if (y0 <= 0. || y1 <= 0.) {
1217 // Use linear interpolation if one of the values is negative.
1218 eps1 = y0 + (e - x0) * (y1 - y0) / (x1 - x0);
1219 } else {
1220 // Otherwise use log-log interpolation.
1221 const double lny0 = log(y0);
1222 const double lny1 = log(y1);
1223 eps1 = lny0 + (lnx - lnx0) * (lny1 - lny0) / (lnx1 - lnx0);
1224 eps1 = exp(eps1);
1225 }
1226
1227 // Interpolate the imaginary part of dielectric function,
1228 // using log-log interpolation.
1229 const double lnz0 = log(m_opticalDataEpsilon[it0 - begin].second);
1230 const double lnz1 = log(m_opticalDataEpsilon[it1 - begin].second);
1231 eps2 = lnz0 + (lnx - lnx0) * (lnz1 - lnz0) / (lnx1 - lnx0);
1232 eps2 = exp(eps2);
1233 return true;
1234}
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377

◆ GetDoping()

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

Retrieve doping concentration.

Definition at line 79 of file MediumSilicon.cc.

79 {
80 type = m_dopingType;
81 c = m_dopingConcentration;
82}

◆ GetElectronBandPopulation()

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

Definition at line 1136 of file MediumSilicon.cc.

1136 {
1137 const int nBands = m_nValleysX + m_nValleysL + 1;
1138 if (band < 0 || band >= nBands) {
1139 std::cerr << m_className << "::GetElectronBandPopulation:\n";
1140 std::cerr << " Band index (" << band << ") out of range.\n";
1141 return 0;
1142 }
1143 return m_nCollElectronBand[band];
1144}

◆ GetElectronCollision()

bool Garfield::MediumSilicon::GetElectronCollision ( const double  e,
int &  type,
int &  level,
double &  e1,
double &  dx,
double &  dy,
double &  dz,
std::vector< std::pair< int, double > > &  secondaries,
int &  ndxc,
int &  band 
)
overridevirtual

Sample the collision type. Update energy and direction vector.

Reimplemented from Garfield::Medium.

Definition at line 786 of file MediumSilicon.cc.

789 {
790 if (e > m_eFinalG) {
791 std::cerr << m_className << "::GetElectronCollision:\n"
792 << " Requested electron energy (" << e << " eV) exceeds the "
793 << "current energy range (" << m_eFinalG << " eV).\n"
794 << " Increasing energy range to " << 1.05 * e << " eV.\n";
795 SetMaxElectronEnergy(1.05 * e);
796 } else if (e <= 0.) {
797 std::cerr << m_className << "::GetElectronCollision:\n"
798 << " Electron energy must be greater than zero.\n";
799 return false;
800 }
801
802 if (m_isChanged) {
803 if (!UpdateTransportParameters()) {
804 std::cerr << m_className << "::GetElectronCollision:\n"
805 << " Error calculating the collision rates table.\n";
806 return false;
807 }
808 m_isChanged = false;
809 }
810
811 // Energy loss
812 double loss = 0.;
813 // Sample the scattering process.
814 if (band >= 0 && band < m_nValleysX) {
815 // X valley
816 // Get the energy interval.
817 int iE = int(e / m_eStepXL);
818 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
819 if (iE < 0) iE = 0;
820 // Select the scattering process.
821 const double r = RndmUniform();
822 if (r <= m_cfElectronsX[iE][0]) {
823 level = 0;
824 } else if (r >= m_cfElectronsX[iE][m_nLevelsX - 1]) {
825 level = m_nLevelsX - 1;
826 } else {
827 const auto begin = m_cfElectronsX[iE].cbegin();
828 level = std::lower_bound(begin, begin + m_nLevelsX, r) - begin;
829 }
830
831 // Get the collision type.
832 type = m_scatTypeElectronsX[level];
833 // Fill the collision counters.
834 ++m_nCollElectronDetailed[level];
835 ++m_nCollElectronBand[band];
836 if (type == ElectronCollisionTypeAcousticPhonon) {
837 ++m_nCollElectronAcoustic;
838 } else if (type == ElectronCollisionTypeIntervalleyG) {
839 // Intervalley scattering (g type)
840 ++m_nCollElectronIntervalley;
841 // Final valley is in opposite direction.
842 switch (band) {
843 case 0:
844 band = 1;
845 break;
846 case 1:
847 band = 0;
848 break;
849 case 2:
850 band = 3;
851 break;
852 case 3:
853 band = 2;
854 break;
855 case 4:
856 band = 5;
857 break;
858 case 5:
859 band = 4;
860 break;
861 default:
862 break;
863 }
864 } else if (type == ElectronCollisionTypeIntervalleyF) {
865 // Intervalley scattering (f type)
866 ++m_nCollElectronIntervalley;
867 // Final valley is perpendicular.
868 switch (band) {
869 case 0:
870 case 1:
871 band = int(RndmUniform() * 4) + 2;
872 break;
873 case 2:
874 case 3:
875 band = int(RndmUniform() * 4);
876 if (band > 1) band += 2;
877 break;
878 case 4:
879 case 5:
880 band = int(RndmUniform() * 4);
881 break;
882 }
883 } else if (type == ElectronCollisionTypeInterbandXL) {
884 // XL scattering
885 ++m_nCollElectronIntervalley;
886 // Final valley is in L band.
887 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
888 if (band >= m_nValleysX + m_nValleysL)
889 band = m_nValleysX + m_nValleysL - 1;
890 } else if (type == ElectronCollisionTypeInterbandXG) {
891 ++m_nCollElectronIntervalley;
892 band = m_nValleysX + m_nValleysL;
893 } else if (type == ElectronCollisionTypeImpurity) {
894 ++m_nCollElectronImpurity;
895 } else if (type == ElectronCollisionTypeIonisation) {
896 ++m_nCollElectronIonisation;
897 } else {
898 std::cerr << m_className << "::GetElectronCollision:\n";
899 std::cerr << " Unexpected collision type (" << type << ").\n";
900 }
901
902 // Get the energy loss.
903 loss = m_energyLossElectronsX[level];
904
905 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
906 // L valley
907 // Get the energy interval.
908 int iE = int(e / m_eStepXL);
909 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
910 if (iE < m_ieMinL) iE = m_ieMinL;
911 // Select the scattering process.
912 const double r = RndmUniform();
913 if (r <= m_cfElectronsL[iE][0]) {
914 level = 0;
915 } else if (r >= m_cfElectronsL[iE][m_nLevelsL - 1]) {
916 level = m_nLevelsL - 1;
917 } else {
918 const auto begin = m_cfElectronsL[iE].cbegin();
919 level = std::lower_bound(begin, begin + m_nLevelsL, r) - begin;
920 }
921
922 // Get the collision type.
923 type = m_scatTypeElectronsL[level];
924 // Fill the collision counters.
925 ++m_nCollElectronDetailed[m_nLevelsX + level];
926 ++m_nCollElectronBand[band];
927 if (type == ElectronCollisionTypeAcousticPhonon) {
928 ++m_nCollElectronAcoustic;
929 } else if (type == ElectronCollisionTypeOpticalPhonon) {
930 ++m_nCollElectronOptical;
931 } else if (type == ElectronCollisionTypeIntervalleyG ||
932 type == ElectronCollisionTypeIntervalleyF) {
933 // Equivalent intervalley scattering
934 ++m_nCollElectronIntervalley;
935 // Randomise the final valley.
936 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
937 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysX + m_nValleysL;
938 } else if (type == ElectronCollisionTypeInterbandXL) {
939 // LX scattering
940 ++m_nCollElectronIntervalley;
941 // Randomise the final valley.
942 band = int(RndmUniform() * m_nValleysX);
943 if (band >= m_nValleysX) band = m_nValleysX - 1;
944 } else if (type == ElectronCollisionTypeInterbandLG) {
945 // LG scattering
946 ++m_nCollElectronIntervalley;
947 band = m_nValleysX + m_nValleysL;
948 } else if (type == ElectronCollisionTypeImpurity) {
949 ++m_nCollElectronImpurity;
950 } else if (type == ElectronCollisionTypeIonisation) {
951 ++m_nCollElectronIonisation;
952 } else {
953 std::cerr << m_className << "::GetElectronCollision:\n";
954 std::cerr << " Unexpected collision type (" << type << ").\n";
955 }
956
957 // Get the energy loss.
958 loss = m_energyLossElectronsL[level];
959 } else if (band == m_nValleysX + m_nValleysL) {
960 // Higher bands
961 // Get the energy interval.
962 int iE = int(e / m_eStepG);
963 if (iE >= nEnergyStepsG) iE = nEnergyStepsG - 1;
964 if (iE < m_ieMinG) iE = m_ieMinG;
965 // Select the scattering process.
966 const double r = RndmUniform();
967 if (r <= m_cfElectronsG[iE][0]) {
968 level = 0;
969 } else if (r >= m_cfElectronsG[iE][m_nLevelsG - 1]) {
970 level = m_nLevelsG - 1;
971 } else {
972 const auto begin = m_cfElectronsG[iE].cbegin();
973 level = std::lower_bound(begin, begin + m_nLevelsG, r) - begin;
974 }
975
976 // Get the collision type.
977 type = m_scatTypeElectronsG[level];
978 // Fill the collision counters.
979 ++m_nCollElectronDetailed[m_nLevelsX + m_nLevelsL + level];
980 ++m_nCollElectronBand[band];
981 if (type == ElectronCollisionTypeAcousticPhonon) {
982 ++m_nCollElectronAcoustic;
983 } else if (type == ElectronCollisionTypeOpticalPhonon) {
984 ++m_nCollElectronOptical;
985 } else if (type == ElectronCollisionTypeIntervalleyG ||
986 type == ElectronCollisionTypeIntervalleyF) {
987 // Equivalent intervalley scattering
988 ++m_nCollElectronIntervalley;
989 } else if (type == ElectronCollisionTypeInterbandXG) {
990 // GX scattering
991 ++m_nCollElectronIntervalley;
992 // Randomise the final valley.
993 band = int(RndmUniform() * m_nValleysX);
994 if (band >= m_nValleysX) band = m_nValleysX - 1;
995 } else if (type == ElectronCollisionTypeInterbandLG) {
996 // GL scattering
997 ++m_nCollElectronIntervalley;
998 // Randomise the final valley.
999 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
1000 if (band >= m_nValleysX + m_nValleysL)
1001 band = m_nValleysX + m_nValleysL - 1;
1002 } else if (type == ElectronCollisionTypeIonisation) {
1003 ++m_nCollElectronIonisation;
1004 } else {
1005 std::cerr << m_className << "::GetElectronCollision:\n";
1006 std::cerr << " Unexpected collision type (" << type << ").\n";
1007 }
1008
1009 // Get the energy loss.
1010 loss = m_energyLossElectronsG[level];
1011 } else {
1012 std::cerr << m_className << "::GetElectronCollision:\n";
1013 std::cerr << " Band index (" << band << ") out of range.\n";
1014 return false;
1015 }
1016
1017 // Secondaries
1018 ndxc = 0;
1019 // Ionising collision
1020 if (type == ElectronCollisionTypeIonisation) {
1021 double ee = 0., eh = 0.;
1022 ComputeSecondaries(e, ee, eh);
1023 loss = ee + eh + m_bandGap;
1024 // Add the secondary electron.
1025 secondaries.emplace_back(std::make_pair(IonProdTypeElectron, ee));
1026 // Add the hole.
1027 secondaries.emplace_back(std::make_pair(IonProdTypeHole, eh));
1028 }
1029
1030 if (e < loss) loss = e - 0.00001;
1031 // Update the energy.
1032 e1 = e - loss;
1033 if (e1 < Small) e1 = Small;
1034
1035 // Update the momentum.
1036 if (band >= 0 && band < m_nValleysX) {
1037 // X valleys
1038 double pstar = sqrt(2. * ElectronMass * e1);
1039 if (m_nonParabolic) {
1040 const double alpha = m_alphaX;
1041 pstar *= sqrt(1. + alpha * e1);
1042 }
1043
1044 const double ctheta = 1. - 2. * RndmUniform();
1045 const double stheta = sqrt(1. - ctheta * ctheta);
1046 const double phi = TwoPi * RndmUniform();
1047
1048 if (m_anisotropic) {
1049 const double pl = pstar * sqrt(m_mLongX);
1050 const double pt = pstar * sqrt(m_mTransX);
1051 switch (band) {
1052 case 0:
1053 case 1:
1054 // 100
1055 px = pl * ctheta;
1056 py = pt * cos(phi) * stheta;
1057 pz = pt * sin(phi) * stheta;
1058 break;
1059 case 2:
1060 case 3:
1061 // 010
1062 px = pt * sin(phi) * stheta;
1063 py = pl * ctheta;
1064 pz = pt * cos(phi) * stheta;
1065 break;
1066 case 4:
1067 case 5:
1068 // 001
1069 px = pt * cos(phi) * stheta;
1070 py = pt * sin(phi) * stheta;
1071 pz = pl * ctheta;
1072 break;
1073 default:
1074 return false;
1075 }
1076 } else {
1077 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
1078 px = pstar * cos(phi) * stheta;
1079 py = pstar * sin(phi) * stheta;
1080 pz = pstar * ctheta;
1081 }
1082 return true;
1083
1084 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
1085 // L valleys
1086 double pstar = sqrt(2. * ElectronMass * (e1 - m_eMinL));
1087 if (m_nonParabolic) {
1088 const double alpha = m_alphaL;
1089 pstar *= sqrt(1. + alpha * (e1 - m_eMinL));
1090 }
1091 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
1092 RndmDirection(px, py, pz, pstar);
1093 return true;
1094 }
1095 const double pstar = sqrt(2. * ElectronMass * e1);
1096 RndmDirection(px, py, pz, pstar);
1097 return true;
1098}
void ComputeSecondaries(const double e0, double &ee, double &eh)
bool SetMaxElectronEnergy(const double e)
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:107
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384

◆ GetElectronCollisionRate()

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

Collision rate [ns-1] for given electron energy.

Reimplemented from Garfield::Medium.

Definition at line 734 of file MediumSilicon.cc.

734 {
735 if (e <= 0.) {
736 std::cerr << m_className << "::GetElectronCollisionRate:\n"
737 << " Electron energy must be positive.\n";
738 return 0.;
739 }
740
741 if (e > m_eFinalG) {
742 std::cerr << m_className << "::GetElectronCollisionRate:\n"
743 << " Collision rate at " << e << " eV (band " << band
744 << ") is not included in the current table.\n"
745 << " Increasing energy range to " << 1.05 * e << " eV.\n";
746 SetMaxElectronEnergy(1.05 * e);
747 }
748
749 if (m_isChanged) {
750 if (!UpdateTransportParameters()) {
751 std::cerr << m_className << "::GetElectronCollisionRate:\n"
752 << " Error calculating the collision rates table.\n";
753 return 0.;
754 }
755 m_isChanged = false;
756 }
757
758 if (band >= 0 && band < m_nValleysX) {
759 int iE = int(e / m_eStepXL);
760 if (iE >= nEnergyStepsXL)
761 iE = nEnergyStepsXL - 1;
762 else if (iE < 0)
763 iE = 0;
764 return m_cfTotElectronsX[iE];
765 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
766 int iE = int(e / m_eStepXL);
767 if (iE >= nEnergyStepsXL)
768 iE = nEnergyStepsXL - 1;
769 else if (iE < m_ieMinL)
770 iE = m_ieMinL;
771 return m_cfTotElectronsL[iE];
772 } else if (band == m_nValleysX + m_nValleysL) {
773 int iE = int(e / m_eStepG);
774 if (iE >= nEnergyStepsG)
775 iE = nEnergyStepsG - 1;
776 else if (iE < m_ieMinG)
777 iE = m_ieMinG;
778 return m_cfTotElectronsG[iE];
779 }
780
781 std::cerr << m_className << "::GetElectronCollisionRate:\n"
782 << " Band index (" << band << ") out of range.\n";
783 return 0.;
784}

◆ GetElectronEnergy()

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

Dispersion relation (energy vs. wave vector)

Reimplemented from Garfield::Medium.

Definition at line 520 of file MediumSilicon.cc.

522 {
523 // Effective masses
524 double mx = ElectronMass, my = ElectronMass, mz = ElectronMass;
525 // Energy offset
526 double e0 = 0.;
527 if (band >= 0 && band < m_nValleysX) {
528 // X valley
529 if (m_anisotropic) {
530 switch (band) {
531 case 0:
532 case 1:
533 // X 100, -100
534 mx *= m_mLongX;
535 my *= m_mTransX;
536 mz *= m_mTransX;
537 break;
538 case 2:
539 case 3:
540 // X 010, 0-10
541 mx *= m_mTransX;
542 my *= m_mLongX;
543 mz *= m_mTransX;
544 break;
545 case 4:
546 case 5:
547 // X 001, 00-1
548 mx *= m_mTransX;
549 my *= m_mTransX;
550 mz *= m_mLongX;
551 break;
552 default:
553 std::cerr << m_className << "::GetElectronEnergy:\n"
554 << " Unexpected band index " << band << "!\n";
555 break;
556 }
557 } else {
558 // Conduction effective mass
559 const double mc = 3. / (1. / m_mLongX + 2. / m_mTransX);
560 mx *= mc;
561 my *= mc;
562 mz *= mc;
563 }
564 } else if (band < m_nValleysX + m_nValleysL) {
565 // L valley, isotropic approximation
566 e0 = m_eMinL;
567 // Effective mass
568 const double mc = 3. / (1. / m_mLongL + 2. / m_mTransL);
569 mx *= mc;
570 my *= mc;
571 mz *= mc;
572 } else if (band == m_nValleysX + m_nValleysL) {
573 // Higher band(s)
574 }
575
576 if (m_nonParabolic) {
577 // Non-parabolicity parameter
578 double alpha = 0.;
579 if (band < m_nValleysX) {
580 // X valley
581 alpha = m_alphaX;
582 } else if (band < m_nValleysX + m_nValleysL) {
583 // L valley
584 alpha = m_alphaL;
585 }
586
587 const double p2 = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
588 if (alpha > 0.) {
589 const double e = 0.5 * (sqrt(1. + 4 * alpha * p2) - 1.) / alpha;
590 const double a = SpeedOfLight / (1. + 2 * alpha * e);
591 vx = a * px / mx;
592 vy = a * py / my;
593 vz = a * pz / mz;
594 return e0 + e;
595 }
596 }
597
598 const double e = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
599 vx = SpeedOfLight * px / mx;
600 vy = SpeedOfLight * py / my;
601 vz = SpeedOfLight * pz / mz;
602 return e0 + e;
603}

◆ GetElectronMomentum()

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

Sample the momentum vector for a given energy (only meaningful in semiconductors).

Reimplemented from Garfield::Medium.

Definition at line 605 of file MediumSilicon.cc.

606 {
607 int nBands = m_nValleysX;
608 if (e > m_eMinL) nBands += m_nValleysL;
609 if (e > m_eMinG) ++nBands;
610
611 // If the band index is out of range, choose one at random.
612 if (band < 0 || band > m_nValleysX + m_nValleysL ||
613 (e < m_eMinL || band >= m_nValleysX) ||
614 (e < m_eMinG || band == m_nValleysX + m_nValleysL)) {
615 if (e < m_eMinL) {
616 band = int(m_nValleysX * RndmUniform());
617 if (band >= m_nValleysX) band = m_nValleysX - 1;
618 } else {
619 const double dosX = GetConductionBandDensityOfStates(e, 0);
620 const double dosL = GetConductionBandDensityOfStates(e, m_nValleysX);
621 const double dosG =
622 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
623 const double dosSum = m_nValleysX * dosX + m_nValleysL * dosL + dosG;
624 if (dosSum < Small) {
625 band = m_nValleysX + m_nValleysL;
626 } else {
627 const double r = RndmUniform() * dosSum;
628 if (r < dosX) {
629 band = int(m_nValleysX * RndmUniform());
630 if (band >= m_nValleysX) band = m_nValleysX - 1;
631 } else if (r < dosX + dosL) {
632 band = m_nValleysX + int(m_nValleysL * RndmUniform());
633 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysL - 1;
634 } else {
635 band = m_nValleysX + m_nValleysL;
636 }
637 }
638 }
639 if (m_debug) {
640 std::cout << m_className << "::GetElectronMomentum:\n"
641 << " Randomised band index: " << band << "\n";
642 }
643 }
644 if (band < m_nValleysX) {
645 // X valleys
646 double pstar = sqrt(2. * ElectronMass * e);
647 if (m_nonParabolic) {
648 const double alpha = m_alphaX;
649 pstar *= sqrt(1. + alpha * e);
650 }
651
652 const double ctheta = 1. - 2. * RndmUniform();
653 const double stheta = sqrt(1. - ctheta * ctheta);
654 const double phi = TwoPi * RndmUniform();
655
656 if (m_anisotropic) {
657 const double pl = pstar * sqrt(m_mLongX);
658 const double pt = pstar * sqrt(m_mTransX);
659 switch (band) {
660 case 0:
661 case 1:
662 // 100
663 px = pl * ctheta;
664 py = pt * cos(phi) * stheta;
665 pz = pt * sin(phi) * stheta;
666 break;
667 case 2:
668 case 3:
669 // 010
670 px = pt * sin(phi) * stheta;
671 py = pl * ctheta;
672 pz = pt * cos(phi) * stheta;
673 break;
674 case 4:
675 case 5:
676 // 001
677 px = pt * cos(phi) * stheta;
678 py = pt * sin(phi) * stheta;
679 pz = pl * ctheta;
680 break;
681 default:
682 // Other band; should not occur.
683 std::cerr << m_className << "::GetElectronMomentum:\n"
684 << " Unexpected band index (" << band << ").\n";
685 px = pstar * stheta * cos(phi);
686 py = pstar * stheta * sin(phi);
687 pz = pstar * ctheta;
688 break;
689 }
690 } else {
691 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
692 px = pstar * cos(phi) * stheta;
693 py = pstar * sin(phi) * stheta;
694 pz = pstar * ctheta;
695 }
696 } else if (band < m_nValleysX + m_nValleysL) {
697 // L valleys
698 double pstar = sqrt(2. * ElectronMass * (e - m_eMinL));
699 if (m_nonParabolic) {
700 const double alpha = m_alphaL;
701 pstar *= sqrt(1. + alpha * (e - m_eMinL));
702 }
703 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
704 RndmDirection(px, py, pz, pstar);
705 } else if (band == m_nValleysX + m_nValleysL) {
706 // Higher band
707 const double pstar = sqrt(2. * ElectronMass * e);
708 RndmDirection(px, py, pz, pstar);
709 }
710}

◆ GetElectronNullCollisionRate()

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

Null-collision rate [ns-1].

Reimplemented from Garfield::Medium.

Definition at line 712 of file MediumSilicon.cc.

712 {
713 if (m_isChanged) {
714 if (!UpdateTransportParameters()) {
715 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
716 << " Error calculating the collision rates table.\n";
717 return 0.;
718 }
719 m_isChanged = false;
720 }
721
722 if (band >= 0 && band < m_nValleysX) {
723 return m_cfNullElectronsX;
724 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
725 return m_cfNullElectronsL;
726 } else if (band == m_nValleysX + m_nValleysL) {
727 return m_cfNullElectronsG;
728 }
729 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
730 << " Band index (" << band << ") out of range.\n";
731 return 0.;
732}

◆ GetMaxElectronEnergy()

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

Definition at line 99 of file MediumSilicon.hh.

99{ return m_eFinalG; }

◆ GetNumberOfElectronBands()

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

Definition at line 1132 of file MediumSilicon.cc.

1132 {
1133 return m_nValleysX + m_nValleysL + 1;
1134}

◆ GetNumberOfElectronCollisions() [1/2]

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

Definition at line 1111 of file MediumSilicon.cc.

1111 {
1112 return m_nCollElectronAcoustic + m_nCollElectronOptical +
1113 m_nCollElectronIntervalley + m_nCollElectronImpurity +
1114 m_nCollElectronIonisation;
1115}

◆ GetNumberOfElectronCollisions() [2/2]

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

Definition at line 1121 of file MediumSilicon.cc.

1122 {
1123 const unsigned int nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1124 if (level >= nLevels) {
1125 std::cerr << m_className << "::GetNumberOfElectronCollisions:\n"
1126 << " Scattering rate term (" << level << ") does not exist.\n";
1127 return 0;
1128 }
1129 return m_nCollElectronDetailed[level];
1130}

◆ GetNumberOfLevels()

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

Definition at line 1117 of file MediumSilicon.cc.

1117 {
1118 return m_nLevelsX + m_nLevelsL + m_nLevelsG;
1119}

◆ GetOpticalDataRange()

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

Get the energy range [eV] of the available optical data.

Reimplemented from Garfield::Medium.

Definition at line 1146 of file MediumSilicon.cc.

1147 {
1148 if (i != 0) {
1149 std::cerr << m_className << "::GetOpticalDataRange: Index out of range.\n";
1150 }
1151
1152 // Make sure the optical data table has been loaded.
1153 if (m_opticalDataEnergies.empty()) {
1154 if (!LoadOpticalData(m_opticalDataFile)) {
1155 std::cerr << m_className << "::GetOpticalDataRange:\n"
1156 << " Optical data table could not be loaded.\n";
1157 return false;
1158 }
1159 }
1160
1161 emin = m_opticalDataEnergies.front();
1162 emax = m_opticalDataEnergies.back();
1163 if (m_debug) {
1164 std::cout << m_className << "::GetOpticalDataRange:\n"
1165 << " " << emin << " < E [eV] < " << emax << "\n";
1166 }
1167 return true;
1168}

◆ GetValenceBandDensityOfStates()

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

Definition at line 2898 of file MediumSilicon.cc.

2899 {
2900 if (band <= 0) {
2901 // Total (full-band) density of states.
2902 const int nPoints = m_fbDosValence.size();
2903 int iE = int(e / m_eStepDos);
2904 if (iE >= nPoints || iE < 0) {
2905 return 0.;
2906 } else if (iE == nPoints - 1) {
2907 return m_fbDosValence[nPoints - 1];
2908 }
2909
2910 const double dos =
2911 m_fbDosValence[iE] +
2912 (m_fbDosValence[iE + 1] - m_fbDosValence[iE]) * (e / m_eStepDos - iE);
2913 return dos * 1.e21;
2914 }
2915
2916 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2917 << " Band index (" << band << ") out of range.\n";
2918 return 0.;
2919}

◆ 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 
)
overridevirtual

Attachment coefficient [cm-1].

Reimplemented from Garfield::Medium.

Definition at line 356 of file MediumSilicon.cc.

359 {
360 eta = 0.;
361 if (m_isChanged) {
362 if (!UpdateTransportParameters()) {
363 std::cerr << m_className << "::HoleAttachment:\n"
364 << " Error calculating the transport parameters.\n";
365 return false;
366 }
367 m_isChanged = false;
368 }
369
370 if (!m_hAtt.empty()) {
371 // Interpolation in user table.
372 return Medium::HoleAttachment(ex, ey, ez, bx, by, bz, eta);
373 }
374
375 switch (m_trappingModel) {
376 case 0:
377 eta = m_hTrapCs * m_hTrapDensity;
378 break;
379 case 1:
380 double vx, vy, vz;
381 HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
382 eta = m_hTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
383 if (eta > 0.) eta = -1. / eta;
384 break;
385 default:
386 std::cerr << m_className << "::HoleAttachment: Unknown model. Bug!\n";
387 return false;
388 }
389
390 return true;
391}
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) override
Drift velocity [cm / ns].
std::vector< std::vector< std::vector< double > > > m_hAtt
Definition: Medium.hh:559
virtual bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Attachment coefficient [cm-1].
Definition: Medium.cc:547

◆ HoleMobility()

double Garfield::MediumSilicon::HoleMobility ( )
inlineoverridevirtual

Low-field mobility [cm2 V-1 ns-1].

Reimplemented from Garfield::Medium.

Definition at line 55 of file MediumSilicon.hh.

55{ return m_hMobility; }

◆ 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 
)
overridevirtual

Ionisation coefficient [cm-1].

Reimplemented from Garfield::Medium.

Definition at line 321 of file MediumSilicon.cc.

324 {
325 alpha = 0.;
326 if (m_isChanged) {
327 if (!UpdateTransportParameters()) {
328 std::cerr << m_className << "::HoleTownsend:\n"
329 << " Error calculating the transport parameters.\n";
330 return false;
331 }
332 m_isChanged = false;
333 }
334
335 if (!m_hAlp.empty()) {
336 // Interpolation in user table.
337 return Medium::HoleTownsend(ex, ey, ez, bx, by, bz, alpha);
338 }
339
340 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
341
342 switch (m_impactIonisationModel) {
343 case ImpactIonisation::VanOverstraeten:
344 return HoleImpactIonisationVanOverstraetenDeMan(emag, alpha);
345 case ImpactIonisation::Grant:
346 return HoleImpactIonisationGrant(emag, alpha);
347 case ImpactIonisation::Massey:
348 return HoleImpactIonisationMassey(emag, alpha);
349 default:
350 std::cerr << m_className << "::HoleTownsend: Unknown model. Bug!\n";
351 break;
352 }
353 return false;
354}
virtual bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Ionisation coefficient [cm-1].
Definition: Medium.cc:534
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition: Medium.hh:558

◆ 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 
)
overridevirtual

Drift velocity [cm / ns].

Reimplemented from Garfield::Medium.

Definition at line 265 of file MediumSilicon.cc.

268 {
269 vx = vy = vz = 0.;
270 if (m_isChanged) {
271 if (!UpdateTransportParameters()) {
272 std::cerr << m_className << "::HoleVelocity:\n"
273 << " Error calculating the transport parameters.\n";
274 return false;
275 }
276 m_isChanged = false;
277 }
278
279 if (!m_hVelE.empty()) {
280 // Interpolation in user table.
281 return Medium::HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
282 }
283
284 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
285
286 // Calculate the mobility
287 double mu;
288 switch (m_highFieldMobilityModel) {
289 case HighFieldMobility::Minimos:
290 HoleMobilityMinimos(emag, mu);
291 break;
292 case HighFieldMobility::Canali:
293 HoleMobilityCanali(emag, mu);
294 break;
295 case HighFieldMobility::Reggiani:
296 HoleMobilityReggiani(emag, mu);
297 break;
298 default:
299 mu = m_hMobility;
300 }
301
302 const double b2 = bx * bx + by * by + bz * bz;
303 if (b2 < Small) {
304 vx = mu * ex;
305 vy = mu * ey;
306 vz = mu * ez;
307 } else {
308 // Hall mobility
309 const double muH = m_hHallFactor * mu;
310 const double mu2 = muH * muH;
311 const double f = mu / (1. + mu2 * b2);
312 const double eb = bx * ex + by * ey + bz * ez;
313 // Compute the drift velocity using the Langevin equation.
314 vx = f * (ex + muH * (ey * bz - ez * by) + mu2 * bx * eb);
315 vy = f * (ey + muH * (ez * bx - ex * bz) + mu2 * by * eb);
316 vz = f * (ez + muH * (ex * by - ey * bx) + mu2 * bz * eb);
317 }
318 return true;
319}
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)
Drift velocity [cm / ns].
Definition: Medium.cc:513
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition: Medium.hh:553

Referenced by HoleAttachment().

◆ Initialise()

bool Garfield::MediumSilicon::Initialise ( )

Definition at line 1236 of file MediumSilicon.cc.

1236 {
1237 if (!m_isChanged) {
1238 if (m_debug) {
1239 std::cerr << m_className << "::Initialise: Nothing changed.\n";
1240 }
1241 return true;
1242 }
1243 if (!UpdateTransportParameters()) {
1244 std::cerr << m_className << "::Initialise:\n Error preparing "
1245 << "transport parameters/calculating collision rates.\n";
1246 return false;
1247 }
1248 return true;
1249}

◆ IsSemiconductor()

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

Is this medium a semiconductor?

Reimplemented from Garfield::Medium.

Definition at line 20 of file MediumSilicon.hh.

20{ return true; }

◆ ResetCollisionCounters()

void Garfield::MediumSilicon::ResetCollisionCounters ( )

Definition at line 1100 of file MediumSilicon.cc.

1100 {
1101 m_nCollElectronAcoustic = m_nCollElectronOptical = 0;
1102 m_nCollElectronIntervalley = 0;
1103 m_nCollElectronImpurity = 0;
1104 m_nCollElectronIonisation = 0;
1105 const auto nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1106 m_nCollElectronDetailed.assign(nLevels, 0);
1107 const auto nBands = m_nValleysX + m_nValleysL + 1;
1108 m_nCollElectronBand.assign(nBands, 0);
1109}

◆ SetDiffusionScaling()

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

Apply a scaling factor to the diffusion coefficients.

Definition at line 95 of file MediumSilicon.hh.

95{ m_diffScale = d; }

◆ SetDoping()

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

Set doping concentration [cm-3] and type ('i', 'n', 'p').

Definition at line 43 of file MediumSilicon.cc.

43 {
44 if (toupper(type) == 'N') {
45 m_dopingType = 'n';
46 if (c > Small) {
47 m_dopingConcentration = c;
48 } else {
49 std::cerr << m_className << "::SetDoping:\n"
50 << " Doping concentration must be greater than zero.\n"
51 << " Using default value for n-type silicon "
52 << "(10^12 cm-3) instead.\n";
53 m_dopingConcentration = 1.e12;
54 }
55 } else if (toupper(type) == 'P') {
56 m_dopingType = 'p';
57 if (c > Small) {
58 m_dopingConcentration = c;
59 } else {
60 std::cerr << m_className << "::SetDoping:\n"
61 << " Doping concentration must be greater than zero.\n"
62 << " Using default value for p-type silicon "
63 << "(10^18 cm-3) instead.\n";
64 m_dopingConcentration = 1.e18;
65 }
66 } else if (toupper(type) == 'I') {
67 m_dopingType = 'i';
68 m_dopingConcentration = 0.;
69 } else {
70 std::cerr << m_className << "::SetDoping:\n"
71 << " Unknown dopant type (" << type << ").\n"
72 << " Available types are n, p and i (intrinsic).\n";
73 return;
74 }
75
76 m_isChanged = true;
77}

◆ SetDopingMobilityModelMasetti()

void Garfield::MediumSilicon::SetDopingMobilityModelMasetti ( )

Use the Masetti model for the doping-dependence of the mobility (default).

Definition at line 430 of file MediumSilicon.cc.

430 {
431 m_dopingMobilityModel = DopingMobility::Masetti;
432 m_hasUserMobility = false;
433 m_isChanged = true;
434}

◆ SetDopingMobilityModelMinimos()

void Garfield::MediumSilicon::SetDopingMobilityModelMinimos ( )

Use the Minimos model for the doping-dependence of the mobility.

Definition at line 424 of file MediumSilicon.cc.

424 {
425 m_dopingMobilityModel = DopingMobility::Minimos;
426 m_hasUserMobility = false;
427 m_isChanged = true;
428}

◆ SetHighFieldMobilityModelCanali()

void Garfield::MediumSilicon::SetHighFieldMobilityModelCanali ( )

Parameterize the high-field mobility using the Canali model (default).

Definition at line 474 of file MediumSilicon.cc.

474 {
475 m_highFieldMobilityModel = HighFieldMobility::Canali;
476 m_isChanged = true;
477}

◆ SetHighFieldMobilityModelConstant()

void Garfield::MediumSilicon::SetHighFieldMobilityModelConstant ( )

Make the velocity proportional to the electric field (no saturation).

Definition at line 484 of file MediumSilicon.cc.

484 {
485 m_highFieldMobilityModel = HighFieldMobility::Constant;
486}

◆ SetHighFieldMobilityModelMinimos()

void Garfield::MediumSilicon::SetHighFieldMobilityModelMinimos ( )

Parameterize the high-field mobility using the Minimos model.

Definition at line 469 of file MediumSilicon.cc.

469 {
470 m_highFieldMobilityModel = HighFieldMobility::Minimos;
471 m_isChanged = true;
472}

◆ SetHighFieldMobilityModelReggiani()

void Garfield::MediumSilicon::SetHighFieldMobilityModelReggiani ( )

Parameterize the high-field mobility using the Reggiani model.

Definition at line 479 of file MediumSilicon.cc.

479 {
480 m_highFieldMobilityModel = HighFieldMobility::Reggiani;
481 m_isChanged = true;
482}

◆ SetImpactIonisationModelGrant()

void Garfield::MediumSilicon::SetImpactIonisationModelGrant ( )

Calculate α using the Grant model.

Definition at line 493 of file MediumSilicon.cc.

493 {
494 m_impactIonisationModel = ImpactIonisation::Grant;
495 m_isChanged = true;
496}

◆ SetImpactIonisationModelMassey()

void Garfield::MediumSilicon::SetImpactIonisationModelMassey ( )

Calculate α using the Massey model.

Definition at line 498 of file MediumSilicon.cc.

498 {
499 m_impactIonisationModel = ImpactIonisation::Massey;
500 m_isChanged = true;
501}

◆ SetImpactIonisationModelVanOverstraetenDeMan()

void Garfield::MediumSilicon::SetImpactIonisationModelVanOverstraetenDeMan ( )

Calculate α using the van Overstraeten-de Man model (default).

Definition at line 488 of file MediumSilicon.cc.

488 {
489 m_impactIonisationModel = ImpactIonisation::VanOverstraeten;
490 m_isChanged = true;
491}

◆ SetLatticeMobilityModelMinimos()

void Garfield::MediumSilicon::SetLatticeMobilityModelMinimos ( )

Calculate the lattice mobility using the Minimos model.

Definition at line 406 of file MediumSilicon.cc.

406 {
407 m_latticeMobilityModel = LatticeMobility::Minimos;
408 m_hasUserMobility = false;
409 m_isChanged = true;
410}

◆ SetLatticeMobilityModelReggiani()

void Garfield::MediumSilicon::SetLatticeMobilityModelReggiani ( )

Calculate the lattice mobility using the Reggiani model.

Definition at line 418 of file MediumSilicon.cc.

418 {
419 m_latticeMobilityModel = LatticeMobility::Reggiani;
420 m_hasUserMobility = false;
421 m_isChanged = true;
422}

◆ SetLatticeMobilityModelSentaurus()

void Garfield::MediumSilicon::SetLatticeMobilityModelSentaurus ( )

Calculate the lattice mobility using the Sentaurus model (default).

Definition at line 412 of file MediumSilicon.cc.

412 {
413 m_latticeMobilityModel = LatticeMobility::Sentaurus;
414 m_hasUserMobility = false;
415 m_isChanged = true;
416}

◆ SetLowFieldMobility()

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

Specify the low field values of the electron and hole mobilities.

Definition at line 393 of file MediumSilicon.cc.

393 {
394 if (mue <= 0. || muh <= 0.) {
395 std::cerr << m_className << "::SetLowFieldMobility:\n"
396 << " Mobility must be greater than zero.\n";
397 return;
398 }
399
400 m_eMobility = mue;
401 m_hMobility = muh;
402 m_hasUserMobility = true;
403 m_isChanged = true;
404}

◆ SetMaxElectronEnergy()

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

Definition at line 503 of file MediumSilicon.cc.

503 {
504 if (e <= m_eMinG + Small) {
505 std::cerr << m_className << "::SetMaxElectronEnergy:\n"
506 << " Requested upper electron energy limit (" << e
507 << " eV) is too small.\n";
508 return false;
509 }
510
511 m_eFinalG = e;
512 // Determine the energy interval size.
513 m_eStepG = m_eFinalG / nEnergyStepsG;
514
515 m_isChanged = true;
516
517 return true;
518}

Referenced by GetElectronCollision(), and GetElectronCollisionRate().

◆ SetSaturationVelocity()

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

Specify the saturation velocities of electrons and holes.

Definition at line 436 of file MediumSilicon.cc.

437 {
438 if (vsate <= 0. || vsath <= 0.) {
439 std::cout << m_className << "::SetSaturationVelocity:\n"
440 << " Restoring default values.\n";
441 m_hasUserSaturationVelocity = false;
442 } else {
443 m_eSatVel = vsate;
444 m_hSatVel = vsath;
445 m_hasUserSaturationVelocity = true;
446 }
447
448 m_isChanged = true;
449}

◆ SetSaturationVelocityModelCanali()

void Garfield::MediumSilicon::SetSaturationVelocityModelCanali ( )

Calculate the saturation velocities using the Canali model (default).

Definition at line 457 of file MediumSilicon.cc.

457 {
458 m_saturationVelocityModel = SaturationVelocity::Canali;
459 m_hasUserSaturationVelocity = false;
460 m_isChanged = true;
461}

◆ SetSaturationVelocityModelMinimos()

void Garfield::MediumSilicon::SetSaturationVelocityModelMinimos ( )

Calculate the saturation velocities using the Minimos model.

Definition at line 451 of file MediumSilicon.cc.

451 {
452 m_saturationVelocityModel = SaturationVelocity::Minimos;
453 m_hasUserSaturationVelocity = false;
454 m_isChanged = true;
455}

◆ SetSaturationVelocityModelReggiani()

void Garfield::MediumSilicon::SetSaturationVelocityModelReggiani ( )

Calculate the saturation velocities using the Reggiani model.

Definition at line 463 of file MediumSilicon.cc.

463 {
464 m_saturationVelocityModel = SaturationVelocity::Reggiani;
465 m_hasUserSaturationVelocity = false;
466 m_isChanged = true;
467}

◆ SetTrapCrossSection()

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

Trapping cross-sections for electrons and holes.

Definition at line 84 of file MediumSilicon.cc.

84 {
85 if (ecs < 0.) {
86 std::cerr << m_className << "::SetTrapCrossSection:\n"
87 << " Capture cross-section [cm2] must non-negative.\n";
88 } else {
89 m_eTrapCs = ecs;
90 }
91
92 if (hcs < 0.) {
93 std::cerr << m_className << "::SetTrapCrossSection:\n"
94 << " Capture cross-section [cm2] must be non-negative.n";
95 } else {
96 m_hTrapCs = hcs;
97 }
98
99 m_trappingModel = 0;
100 m_isChanged = true;
101}

◆ SetTrapDensity()

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

Trap density [cm-3], by default set to zero.

Definition at line 103 of file MediumSilicon.cc.

103 {
104 if (n < 0.) {
105 std::cerr << m_className << "::SetTrapDensity:\n"
106 << " Trap density [cm-3] must be non-negative.\n";
107 } else {
108 m_eTrapDensity = n;
109 m_hTrapDensity = n;
110 }
111
112 m_trappingModel = 0;
113 m_isChanged = true;
114}

◆ SetTrappingTime()

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

Set time constant for trapping of electrons and holes [ns].

Definition at line 116 of file MediumSilicon.cc.

116 {
117 if (etau <= 0.) {
118 std::cerr << m_className << "::SetTrappingTime:\n"
119 << " Trapping time [ns-1] must be positive.\n";
120 } else {
121 m_eTrapTime = etau;
122 }
123
124 if (htau <= 0.) {
125 std::cerr << m_className << "::SetTrappingTime:\n"
126 << " Trapping time [ns-1] must be positive.\n";
127 } else {
128 m_hTrapTime = htau;
129 }
130
131 m_trappingModel = 1;
132 m_isChanged = true;
133}

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