Garfield++ 3.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 unsigned int ie, const unsigned int ib, const unsigned int ia, const double v)
 Set an entry in the table of drift speeds along E.
 
bool GetElectronVelocityE (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 Get an entry in the table of drift speeds along E.
 
bool SetElectronVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double v)
 Set an entry in the table of drift speeds along ExB.
 
bool GetElectronVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 Get an entry in the table of drift speeds along ExB.
 
bool SetElectronVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double v)
 Set an entry in the table of drift speeds along Btrans.
 
bool GetElectronVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 Get an entry in the table of drift speeds along Btrans.
 
bool SetElectronLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetElectronLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetElectronTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetElectronTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetElectronTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double alpha)
 Set an entry in the table of Townsend coefficients.
 
bool GetElectronTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &alpha)
 Get an entry in the table of Townsend coefficients.
 
bool SetElectronAttachment (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double eta)
 Set an entry in the table of attachment coefficients.
 
bool GetElectronAttachment (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &eta)
 Get an entry in the table of attachment coefficients.
 
bool SetElectronLorentzAngle (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double lor)
 Set an entry in the table of Lorentz angles.
 
bool GetElectronLorentzAngle (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &lor)
 Get an entry in the table of Lorentz angles.
 
bool SetHoleVelocityE (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double v)
 Set an entry in the table of drift speeds along E.
 
bool GetHoleVelocityE (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 Get an entry in the table of drift speeds along E.
 
bool SetHoleVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double v)
 Set an entry in the table of drift speeds along ExB.
 
bool GetHoleVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 Get an entry in the table of drift speeds along ExB.
 
bool SetHoleVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double v)
 Set an entry in the table of drift speeds along Btrans.
 
bool GetHoleVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 Get an entry in the table of drift speeds along Btrans.
 
bool SetHoleLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetHoleLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetHoleTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetHoleTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetHoleTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double alpha)
 Set an entry in the table of Townsend coefficients.
 
bool GetHoleTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &alpha)
 Get an entry in the table of Townsend coefficients.
 
bool SetHoleAttachment (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double eta)
 Set an entry in the table of attachment coefficients.
 
bool GetHoleAttachment (const unsigned int ie, const unsigned int ib, const unsigned int 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 unsigned int ie, const unsigned int ib, const unsigned int ia, const double mu)
 Set an entry in the table of ion mobilities.
 
bool GetIonMobility (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &mu)
 Get an entry in the table of ion mobilities.
 
bool SetIonLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetIonLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetIonTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetIonTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetIonDissociation (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double diss)
 Set an entry in the table of dissociation coefficients.
 
bool GetIonDissociation (const unsigned int ie, const unsigned int ib, const unsigned int 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 unsigned int i, const unsigned int j, const unsigned int k, const std::string &fcn, std::vector< std::vector< std::vector< double > > > &tab, const double val)
 
bool GetEntry (const unsigned int i, const unsigned int j, const unsigned int 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
 
unsigned int 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 unsigned int 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
 
std::string m_name = ""
 
double m_temperature = 293.15
 
double m_pressure = 760.
 
double m_epsilon = 1.
 
unsigned int m_nComponents = 1
 
double m_z = 1.
 
double m_a = 0.
 
double m_density = 0.
 
bool m_driftable = false
 
bool m_microscopic = false
 
bool m_ionisable = false
 
double m_w = 0.
 
double m_fano = 0.
 
bool m_isChanged = true
 
bool m_debug = false
 
std::vector< double > m_eFields
 
std::vector< double > m_bFields
 
std::vector< double > m_bAngles
 
bool m_tab2d = false
 
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 9 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:531
double m_fano
Definition: Medium.hh:537
virtual void EnableDrift(const bool on=true)
Switch electron/ion/hole on/off.
Definition: Medium.hh:67
std::string m_name
Definition: Medium.hh:513
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:506

◆ ~MediumSilicon()

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

Destructor.

Definition at line 14 of file MediumSilicon.hh.

14{}

Member Function Documentation

◆ ComputeSecondaries()

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

Definition at line 2937 of file MediumSilicon.cc.

2938 {
2939 const int nPointsValence = m_fbDosValence.size();
2940 const int nPointsConduction = m_fbDosConduction.size();
2941 const double widthValenceBand = m_eStepDos * nPointsValence;
2942 const double widthConductionBand = m_eStepDos * nPointsConduction;
2943
2944 bool ok = false;
2945 while (!ok) {
2946 // Sample a hole energy according to the valence band DOS.
2947 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2948 int ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2949 while (RndmUniform() > m_fbDosValence[ih] / m_fbDosMaxV) {
2950 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2951 ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2952 }
2953 // Sample an electron energy according to the conduction band DOS.
2954 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2955 int ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2956 while (RndmUniform() > m_fbDosConduction[ie] / m_fbDosMaxC) {
2957 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2958 ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2959 }
2960 // Calculate the energy of the primary electron.
2961 const double ep = e0 - m_bandGap - eh - ee;
2962 if (ep < Small) continue;
2963 if (ep > 5.) return;
2964 // Check if the primary electron energy is consistent with the DOS.
2965 int ip = std::min(int(ep / m_eStepDos), nPointsConduction - 1);
2966 if (RndmUniform() > m_fbDosConduction[ip] / m_fbDosMaxC) continue;
2967 ok = true;
2968 }
2969}
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 232 of file MediumSilicon.cc.

235 {
236 eta = 0.;
237 if (m_isChanged) {
238 if (!UpdateTransportParameters()) {
239 std::cerr << m_className << "::ElectronAttachment:\n"
240 << " Error calculating the transport parameters.\n";
241 return false;
242 }
243 m_isChanged = false;
244 }
245
246 if (!m_eAtt.empty()) {
247 // Interpolation in user table.
248 return Medium::ElectronAttachment(ex, ey, ez, bx, by, bz, eta);
249 }
250
251 switch (m_trappingModel) {
252 case 0:
253 eta = m_eTrapCs * m_eTrapDensity;
254 break;
255 case 1:
256 double vx, vy, vz;
257 ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
258 eta = m_eTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
259 if (eta > 0.) eta = 1. / eta;
260 break;
261 default:
262 std::cerr << m_className << "::ElectronAttachment: Unknown model. Bug!\n";
263 return false;
264 break;
265 }
266
267 return true;
268}
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:559
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:540
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 40 of file MediumSilicon.hh.

40{ 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 194 of file MediumSilicon.cc.

197 {
198 alpha = 0.;
199 if (m_isChanged) {
200 if (!UpdateTransportParameters()) {
201 std::cerr << m_className << "::ElectronTownsend:\n"
202 << " Error calculating the transport parameters.\n";
203 return false;
204 }
205 m_isChanged = false;
206 }
207
208 if (!m_eAlp.empty()) {
209 // Interpolation in user table.
210 return Medium::ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
211 }
212
213 const double e = sqrt(ex * ex + ey * ey + ez * ez);
214
215 switch (m_impactIonisationModel) {
216 case ImpactIonisation::VanOverstraeten:
217 return ElectronImpactIonisationVanOverstraetenDeMan(e, alpha);
218 break;
219 case ImpactIonisation::Grant:
220 return ElectronImpactIonisationGrant(e, alpha);
221 break;
222 case ImpactIonisation::Massey:
223 return ElectronImpactIonisationMassey(e, alpha);
224 break;
225 default:
226 std::cerr << m_className << "::ElectronTownsend: Unknown model. Bug!\n";
227 break;
228 }
229 return false;
230}
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition: Medium.hh:558
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 e = 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(e, mu);
161 break;
162 case HighFieldMobility::Canali:
163 ElectronMobilityCanali(e, mu);
164 break;
165 case HighFieldMobility::Reggiani:
166 ElectronMobilityReggiani(e, mu);
167 break;
168 default:
169 mu = m_eMobility;
170 break;
171 }
172 mu = -mu;
173
174 const double b = sqrt(bx * bx + by * by + bz * bz);
175 if (b < 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 eb = bx * ex + by * ey + bz * ez;
183 const double mub = muH * b;
184 const double f = mu / (1. + mub * mub);
185 const double muH2 = muH * muH;
186 // Compute the drift velocity using the Langevin equation.
187 vx = f * (ex + muH * (ey * bz - ez * by) + muH2 * bx * eb);
188 vy = f * (ey + muH * (ez * bx - ex * bz) + muH2 * by * eb);
189 vz = f * (ez + muH * (ex * by - ey * bx) + muH2 * bz * eb);
190 }
191 return true;
192}
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:553

Referenced by ElectronAttachment().

◆ EnableAnisotropy()

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

Definition at line 106 of file MediumSilicon.hh.

106{ m_anisotropic = on; }

◆ EnableFullBandDensityOfStates()

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

Definition at line 103 of file MediumSilicon.hh.

103 {
104 m_fullBandDos = on;
105 }

◆ EnableNonParabolicity()

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

Definition at line 102 of file MediumSilicon.hh.

102{ m_nonParabolic = on; }

◆ EnableScatteringRateOutput()

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

Definition at line 101 of file MediumSilicon.hh.

101{ m_cfOutput = on; }

◆ GetConductionBandDensityOfStates()

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

Definition at line 2802 of file MediumSilicon.cc.

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

1189 {
1190 if (i != 0) {
1191 std::cerr << m_className + "::GetDielectricFunction: Index out of range.\n";
1192 return false;
1193 }
1194
1195 // Make sure the optical data table has been loaded.
1196 if (m_opticalDataEnergies.empty()) {
1197 if (!LoadOpticalData(m_opticalDataFile)) {
1198 std::cerr << m_className << "::GetDielectricFunction:\n";
1199 std::cerr << " Optical data table could not be loaded.\n";
1200 return false;
1201 }
1202 }
1203
1204 // Make sure the requested energy is within the range of the table.
1205 const double emin = m_opticalDataEnergies.front();
1206 const double emax = m_opticalDataEnergies.back();
1207 if (e < emin || e > emax) {
1208 std::cerr << m_className << "::GetDielectricFunction:\n"
1209 << " Requested energy (" << e << " eV) "
1210 << " is outside the range of the optical data table.\n"
1211 << " " << emin << " < E [eV] < " << emax << "\n";
1212 eps1 = eps2 = 0.;
1213 return false;
1214 }
1215
1216 // Locate the requested energy in the table.
1217 const auto begin = m_opticalDataEnergies.cbegin();
1218 const auto it1 = std::upper_bound(begin, m_opticalDataEnergies.cend(), e);
1219 if (it1 == begin) {
1220 eps1 = m_opticalDataEpsilon.front().first;
1221 eps2 = m_opticalDataEpsilon.front().second;
1222 return true;
1223 }
1224 const auto it0 = std::prev(it1);
1225
1226 // Interpolate the real part of dielectric function.
1227 const double x0 = *it0;
1228 const double x1 = *it1;
1229 const double lnx0 = log(*it0);
1230 const double lnx1 = log(*it1);
1231 const double lnx = log(e);
1232 const double y0 = m_opticalDataEpsilon[it0 - begin].first;
1233 const double y1 = m_opticalDataEpsilon[it1 - begin].first;
1234 if (y0 <= 0. || y1 <= 0.) {
1235 // Use linear interpolation if one of the values is negative.
1236 eps1 = y0 + (e - x0) * (y1 - y0) / (x1 - x0);
1237 } else {
1238 // Otherwise use log-log interpolation.
1239 const double lny0 = log(y0);
1240 const double lny1 = log(y1);
1241 eps1 = lny0 + (lnx - lnx0) * (lny1 - lny0) / (lnx1 - lnx0);
1242 eps1 = exp(eps1);
1243 }
1244
1245 // Interpolate the imaginary part of dielectric function,
1246 // using log-log interpolation.
1247 const double lnz0 = log(m_opticalDataEpsilon[it0 - begin].second);
1248 const double lnz1 = log(m_opticalDataEpsilon[it1 - begin].second);
1249 eps2 = lnz0 + (lnx - lnx0) * (lnz1 - lnz0) / (lnx1 - lnx0);
1250 eps2 = exp(eps2);
1251 return true;
1252}
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 1154 of file MediumSilicon.cc.

1154 {
1155 const int nBands = m_nValleysX + m_nValleysL + 1;
1156 if (band < 0 || band >= nBands) {
1157 std::cerr << m_className << "::GetElectronBandPopulation:\n";
1158 std::cerr << " Band index (" << band << ") out of range.\n";
1159 return 0;
1160 }
1161 return m_nCollElectronBand[band];
1162}

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

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

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

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

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

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

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

◆ GetElectronNullCollisionRate()

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

Null-collision rate [ns-1].

Reimplemented from Garfield::Medium.

Definition at line 722 of file MediumSilicon.cc.

722 {
723 if (m_isChanged) {
724 if (!UpdateTransportParameters()) {
725 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
726 << " Error calculating the collision rates table.\n";
727 return 0.;
728 }
729 m_isChanged = false;
730 }
731
732 if (band >= 0 && band < m_nValleysX) {
733 return m_cfNullElectronsX;
734 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
735 return m_cfNullElectronsL;
736 } else if (band == m_nValleysX + m_nValleysL) {
737 return m_cfNullElectronsG;
738 }
739 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
740 << " Band index (" << band << ") out of range.\n";
741 return 0.;
742}

◆ GetMaxElectronEnergy()

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

Definition at line 95 of file MediumSilicon.hh.

95{ return m_eFinalG; }

◆ GetNumberOfElectronBands()

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

Definition at line 1150 of file MediumSilicon.cc.

1150 {
1151 return m_nValleysX + m_nValleysL + 1;
1152}

◆ GetNumberOfElectronCollisions() [1/2]

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

Definition at line 1129 of file MediumSilicon.cc.

1129 {
1130 return m_nCollElectronAcoustic + m_nCollElectronOptical +
1131 m_nCollElectronIntervalley + m_nCollElectronImpurity +
1132 m_nCollElectronIonisation;
1133}

◆ GetNumberOfElectronCollisions() [2/2]

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

Definition at line 1139 of file MediumSilicon.cc.

1140 {
1141 const unsigned int nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1142 if (level >= nLevels) {
1143 std::cerr << m_className << "::GetNumberOfElectronCollisions:\n"
1144 << " Scattering rate term (" << level << ") does not exist.\n";
1145 return 0;
1146 }
1147 return m_nCollElectronDetailed[level];
1148}

◆ GetNumberOfLevels()

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

Definition at line 1135 of file MediumSilicon.cc.

1135 {
1136 return m_nLevelsX + m_nLevelsL + m_nLevelsG;
1137}

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

1165 {
1166 if (i != 0) {
1167 std::cerr << m_className << "::GetOpticalDataRange: Index out of range.\n";
1168 }
1169
1170 // Make sure the optical data table has been loaded.
1171 if (m_opticalDataEnergies.empty()) {
1172 if (!LoadOpticalData(m_opticalDataFile)) {
1173 std::cerr << m_className << "::GetOpticalDataRange:\n"
1174 << " Optical data table could not be loaded.\n";
1175 return false;
1176 }
1177 }
1178
1179 emin = m_opticalDataEnergies.front();
1180 emax = m_opticalDataEnergies.back();
1181 if (m_debug) {
1182 std::cout << m_className << "::GetOpticalDataRange:\n"
1183 << " " << emin << " < E [eV] < " << emax << "\n";
1184 }
1185 return true;
1186}

◆ GetValenceBandDensityOfStates()

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

Definition at line 2914 of file MediumSilicon.cc.

2915 {
2916 if (band <= 0) {
2917 // Total (full-band) density of states.
2918 const int nPoints = m_fbDosValence.size();
2919 int iE = int(e / m_eStepDos);
2920 if (iE >= nPoints || iE < 0) {
2921 return 0.;
2922 } else if (iE == nPoints - 1) {
2923 return m_fbDosValence[nPoints - 1];
2924 }
2925
2926 const double dos =
2927 m_fbDosValence[iE] +
2928 (m_fbDosValence[iE + 1] - m_fbDosValence[iE]) * (e / m_eStepDos - iE);
2929 return dos * 1.e21;
2930 }
2931
2932 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2933 << " Band index (" << band << ") out of range.\n";
2934 return 0.;
2935}

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

368 {
369 eta = 0.;
370 if (m_isChanged) {
371 if (!UpdateTransportParameters()) {
372 std::cerr << m_className << "::HoleAttachment:\n"
373 << " Error calculating the transport parameters.\n";
374 return false;
375 }
376 m_isChanged = false;
377 }
378
379 if (!m_hAtt.empty()) {
380 // Interpolation in user table.
381 return Medium::HoleAttachment(ex, ey, ez, bx, by, bz, eta);
382 }
383
384 switch (m_trappingModel) {
385 case 0:
386 eta = m_hTrapCs * m_hTrapDensity;
387 break;
388 case 1:
389 double vx, vy, vz;
390 HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
391 eta = m_hTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
392 if (eta > 0.) eta = 1. / eta;
393 break;
394 default:
395 std::cerr << m_className << "::HoleAttachment: Unknown model. Bug!\n";
396 return false;
397 break;
398 }
399
400 return true;
401}
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:571
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 51 of file MediumSilicon.hh.

51{ 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 327 of file MediumSilicon.cc.

330 {
331 alpha = 0.;
332 if (m_isChanged) {
333 if (!UpdateTransportParameters()) {
334 std::cerr << m_className << "::HoleTownsend:\n"
335 << " Error calculating the transport parameters.\n";
336 return false;
337 }
338 m_isChanged = false;
339 }
340
341 if (!m_hAlp.empty()) {
342 // Interpolation in user table.
343 return Medium::HoleTownsend(ex, ey, ez, bx, by, bz, alpha);
344 }
345
346 const double e = sqrt(ex * ex + ey * ey + ez * ez);
347
348 switch (m_impactIonisationModel) {
349 case ImpactIonisation::VanOverstraeten:
350 return HoleImpactIonisationVanOverstraetenDeMan(e, alpha);
351 break;
352 case ImpactIonisation::Grant:
353 return HoleImpactIonisationGrant(e, alpha);
354 break;
355 case ImpactIonisation::Massey:
356 return HoleImpactIonisationMassey(e, alpha);
357 break;
358 default:
359 std::cerr << m_className << "::HoleTownsend: Unknown model. Bug!\n";
360 break;
361 }
362 return false;
363}
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:570

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

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

Referenced by HoleAttachment().

◆ Initialise()

bool Garfield::MediumSilicon::Initialise ( )

Definition at line 1254 of file MediumSilicon.cc.

1254 {
1255 if (!m_isChanged) {
1256 if (m_debug) {
1257 std::cerr << m_className << "::Initialise: Nothing changed.\n";
1258 }
1259 return true;
1260 }
1261 if (!UpdateTransportParameters()) {
1262 std::cerr << m_className << "::Initialise:\n Error preparing "
1263 << "transport parameters/calculating collision rates.\n";
1264 return false;
1265 }
1266 return true;
1267}

◆ IsSemiconductor()

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

Is this medium a semiconductor?

Reimplemented from Garfield::Medium.

Definition at line 16 of file MediumSilicon.hh.

16{ return true; }

◆ ResetCollisionCounters()

void Garfield::MediumSilicon::ResetCollisionCounters ( )

Definition at line 1118 of file MediumSilicon.cc.

1118 {
1119 m_nCollElectronAcoustic = m_nCollElectronOptical = 0;
1120 m_nCollElectronIntervalley = 0;
1121 m_nCollElectronImpurity = 0;
1122 m_nCollElectronIonisation = 0;
1123 const auto nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1124 m_nCollElectronDetailed.assign(nLevels, 0);
1125 const auto nBands = m_nValleysX + m_nValleysL + 1;
1126 m_nCollElectronBand.assign(nBands, 0);
1127}

◆ SetDiffusionScaling()

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

Apply a scaling factor to the diffusion coefficients.

Definition at line 91 of file MediumSilicon.hh.

91{ 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 440 of file MediumSilicon.cc.

440 {
441 m_dopingMobilityModel = DopingMobility::Masetti;
442 m_hasUserMobility = false;
443 m_isChanged = true;
444}

◆ SetDopingMobilityModelMinimos()

void Garfield::MediumSilicon::SetDopingMobilityModelMinimos ( )

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

Definition at line 434 of file MediumSilicon.cc.

434 {
435 m_dopingMobilityModel = DopingMobility::Minimos;
436 m_hasUserMobility = false;
437 m_isChanged = true;
438}

◆ SetHighFieldMobilityModelCanali()

void Garfield::MediumSilicon::SetHighFieldMobilityModelCanali ( )

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

Definition at line 484 of file MediumSilicon.cc.

484 {
485 m_highFieldMobilityModel = HighFieldMobility::Canali;
486 m_isChanged = true;
487}

◆ SetHighFieldMobilityModelConstant()

void Garfield::MediumSilicon::SetHighFieldMobilityModelConstant ( )

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

Definition at line 494 of file MediumSilicon.cc.

494 {
495 m_highFieldMobilityModel = HighFieldMobility::Constant;
496}

◆ SetHighFieldMobilityModelMinimos()

void Garfield::MediumSilicon::SetHighFieldMobilityModelMinimos ( )

Parameterize the high-field mobility using the Minimos model.

Definition at line 479 of file MediumSilicon.cc.

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

◆ SetHighFieldMobilityModelReggiani()

void Garfield::MediumSilicon::SetHighFieldMobilityModelReggiani ( )

Parameterize the high-field mobility using the Reggiani model.

Definition at line 489 of file MediumSilicon.cc.

489 {
490 m_highFieldMobilityModel = HighFieldMobility::Reggiani;
491 m_isChanged = true;
492}

◆ SetImpactIonisationModelGrant()

void Garfield::MediumSilicon::SetImpactIonisationModelGrant ( )

Calculate α using the Grant model.

Definition at line 503 of file MediumSilicon.cc.

503 {
504 m_impactIonisationModel = ImpactIonisation::Grant;
505 m_isChanged = true;
506}

◆ SetImpactIonisationModelMassey()

void Garfield::MediumSilicon::SetImpactIonisationModelMassey ( )

Calculate α using the Massey model.

Definition at line 508 of file MediumSilicon.cc.

508 {
509 m_impactIonisationModel = ImpactIonisation::Massey;
510 m_isChanged = true;
511}

◆ SetImpactIonisationModelVanOverstraetenDeMan()

void Garfield::MediumSilicon::SetImpactIonisationModelVanOverstraetenDeMan ( )

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

Definition at line 498 of file MediumSilicon.cc.

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

◆ SetLatticeMobilityModelMinimos()

void Garfield::MediumSilicon::SetLatticeMobilityModelMinimos ( )

Calculate the lattice mobility using the Minimos model.

Definition at line 416 of file MediumSilicon.cc.

416 {
417 m_latticeMobilityModel = LatticeMobility::Minimos;
418 m_hasUserMobility = false;
419 m_isChanged = true;
420}

◆ SetLatticeMobilityModelReggiani()

void Garfield::MediumSilicon::SetLatticeMobilityModelReggiani ( )

Calculate the lattice mobility using the Reggiani model.

Definition at line 428 of file MediumSilicon.cc.

428 {
429 m_latticeMobilityModel = LatticeMobility::Reggiani;
430 m_hasUserMobility = false;
431 m_isChanged = true;
432}

◆ SetLatticeMobilityModelSentaurus()

void Garfield::MediumSilicon::SetLatticeMobilityModelSentaurus ( )

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

Definition at line 422 of file MediumSilicon.cc.

422 {
423 m_latticeMobilityModel = LatticeMobility::Sentaurus;
424 m_hasUserMobility = false;
425 m_isChanged = true;
426}

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

403 {
404 if (mue <= 0. || muh <= 0.) {
405 std::cerr << m_className << "::SetLowFieldMobility:\n"
406 << " Mobility must be greater than zero.\n";
407 return;
408 }
409
410 m_eMobility = mue;
411 m_hMobility = muh;
412 m_hasUserMobility = true;
413 m_isChanged = true;
414}

◆ SetMaxElectronEnergy()

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

Definition at line 513 of file MediumSilicon.cc.

513 {
514 if (e <= m_eMinG + Small) {
515 std::cerr << m_className << "::SetMaxElectronEnergy:\n"
516 << " Requested upper electron energy limit (" << e
517 << " eV) is too small.\n";
518 return false;
519 }
520
521 m_eFinalG = e;
522 // Determine the energy interval size.
523 m_eStepG = m_eFinalG / nEnergyStepsG;
524
525 m_isChanged = true;
526
527 return true;
528}

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

447 {
448 if (vsate <= 0. || vsath <= 0.) {
449 std::cout << m_className << "::SetSaturationVelocity:\n"
450 << " Restoring default values.\n";
451 m_hasUserSaturationVelocity = false;
452 } else {
453 m_eSatVel = vsate;
454 m_hSatVel = vsath;
455 m_hasUserSaturationVelocity = true;
456 }
457
458 m_isChanged = true;
459}

◆ SetSaturationVelocityModelCanali()

void Garfield::MediumSilicon::SetSaturationVelocityModelCanali ( )

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

Definition at line 467 of file MediumSilicon.cc.

467 {
468 m_saturationVelocityModel = SaturationVelocity::Canali;
469 m_hasUserSaturationVelocity = false;
470 m_isChanged = true;
471}

◆ SetSaturationVelocityModelMinimos()

void Garfield::MediumSilicon::SetSaturationVelocityModelMinimos ( )

Calculate the saturation velocities using the Minimos model.

Definition at line 461 of file MediumSilicon.cc.

461 {
462 m_saturationVelocityModel = SaturationVelocity::Minimos;
463 m_hasUserSaturationVelocity = false;
464 m_isChanged = true;
465}

◆ SetSaturationVelocityModelReggiani()

void Garfield::MediumSilicon::SetSaturationVelocityModelReggiani ( )

Calculate the saturation velocities using the Reggiani model.

Definition at line 473 of file MediumSilicon.cc.

473 {
474 m_saturationVelocityModel = SaturationVelocity::Reggiani;
475 m_hasUserSaturationVelocity = false;
476 m_isChanged = true;
477}

◆ 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: