Garfield++ 5.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 SetImpactIonisationModelOkutoCrowell ()
 Calculate α using the Okuto-Crowell 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 ElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< std::pair< Particle, 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 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 transport 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 () const
 Get the W value.
 
void SetFanoFactor (const double f)
 Set the Fano factor.
 
double GetFanoFactor () const
 Get the Fano factor.
 
void PlotVelocity (const std::string &carriers, TPad *pad)
 Plot the drift velocity as function of the electric field.
 
void PlotDiffusion (const std::string &carriers, TPad *pad)
 Plot the diffusion coefficients as function of the electric field.
 
void PlotTownsend (const std::string &carriers, TPad *pad)
 Plot the Townsend coefficient(s) as function of the electric field.
 
void PlotAttachment (const std::string &carriers, TPad *pad)
 Plot the attachment coefficient(s) as function of the electric field.
 
void PlotAlphaEta (const std::string &carriers, TPad *pad)
 Plot Townsend and attachment coefficients.
 
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 ElectronLorentzAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
 Lorentz angle.
 
virtual unsigned int GetNumberOfDeexcitationProducts () const
 
virtual bool GetDeexcitationProduct (const unsigned int i, double &t, double &s, int &type, double &energy) const
 
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 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)
 Ion 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 ion mobility [cm2 V-1 ns-1].
 
virtual bool NegativeIonVelocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
 Negative ion drift velocity [cm / ns].
 
virtual double NegativeIonMobility ()
 Low-field negative ion 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, const bool negativeIons=false)
 
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.
 
bool SetNegativeIonMobility (const size_t ie, const size_t ib, const size_t ia, const double mu)
 Set an entry in the table of negative ion mobilities.
 
bool GetNegativeIonMobility (const size_t ie, const size_t ib, const size_t ia, double &mu)
 Get an entry in the table of negative ion mobilities.
 
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 ResetNegativeIonMobility ()
 
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 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)
 
- Static Protected Member Functions inherited from Garfield::Medium
static void Langevin (const double ex, const double ey, const double ez, double bx, double by, double bz, const double mu, double &vx, double &vy, double &vz)
 
static void Langevin (const double ex, const double ey, const double ez, double bx, double by, double bz, const double mu, const double muH, double &vx, double &vy, double &vz)
 
- 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
 
std::vector< std::vector< std::vector< double > > > m_nMob
 
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.);
28
29 m_driftable = true;
30 m_ionisable = true;
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:72
std::string m_name
Definition Medium.hh:538
virtual void SetAtomicNumber(const double z)
Set the effective atomic number.
Definition Medium.cc:115
void SetDielectricConstant(const double eps)
Set the relative static dielectric constant.
Definition Medium.cc:92
virtual void SetMassDensity(const double rho)
Set the mass density [g/cm3].
Definition Medium.cc:145
Medium()
Constructor.
Definition Medium.cc:61
virtual void SetAtomicWeight(const double a)
Set the effective atomic weight.
Definition Medium.cc:125
std::string m_className
Definition Medium.hh:529

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

2727 {
2728 const int nPointsValence = m_fbDosValence.size();
2729 const int nPointsConduction = m_fbDosConduction.size();
2730 const double widthValenceBand = m_eStepDos * nPointsValence;
2731 const double widthConductionBand = m_eStepDos * nPointsConduction;
2732
2733 bool ok = false;
2734 while (!ok) {
2735 // Sample a hole energy according to the valence band DOS.
2736 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2737 int ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2738 while (RndmUniform() > m_fbDosValence[ih] / m_fbDosMaxV) {
2739 eh = RndmUniformPos() * std::min(widthValenceBand, e0);
2740 ih = std::min(int(eh / m_eStepDos), nPointsValence - 1);
2741 }
2742 // Sample an electron energy according to the conduction band DOS.
2743 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2744 int ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2745 while (RndmUniform() > m_fbDosConduction[ie] / m_fbDosMaxC) {
2746 ee = RndmUniformPos() * std::min(widthConductionBand, e0);
2747 ie = std::min(int(ee / m_eStepDos), nPointsConduction - 1);
2748 }
2749 // Calculate the energy of the primary electron.
2750 const double ep = e0 - m_bandGap - eh - ee;
2751 if (ep < Small) continue;
2752 if (ep > 5.) return;
2753 // Check if the primary electron energy is consistent with the DOS.
2754 int ip = std::min(int(ep / m_eStepDos), nPointsConduction - 1);
2755 if (RndmUniform() > m_fbDosConduction[ip] / m_fbDosMaxC) continue;
2756 ok = true;
2757 }
2758}
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 ElectronCollision().

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

195 {
196 eta = 0.;
197 if (m_isChanged) {
198 if (!UpdateTransportParameters()) {
199 std::cerr << m_className << "::ElectronAttachment:\n"
200 << " Error calculating the transport parameters.\n";
201 return false;
202 }
203 m_isChanged = false;
204 }
205
206 if (!m_eAtt.empty()) {
207 // Interpolation in user table.
208 return Medium::ElectronAttachment(ex, ey, ez, bx, by, bz, eta);
209 }
210
211 switch (m_trappingModel) {
212 case 0:
213 eta = m_eTrapCs * m_eTrapDensity;
214 break;
215 case 1:
216 double vx, vy, vz;
217 ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
218 eta = m_eTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
219 if (eta > 0.) eta = -1. / eta;
220 break;
221 default:
222 std::cerr << m_className << "::ElectronAttachment: Unknown model. Bug!\n";
223 return false;
224 }
225
226 return true;
227}
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:583
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:481
DoubleAc sqrt(const DoubleAc &f)
Definition DoubleAc.cpp:314

◆ ElectronCollision()

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

Sample the collision type. Update energy and direction vector.

Reimplemented from Garfield::Medium.

Definition at line 717 of file MediumSilicon.cc.

720 {
721 if (e > m_eFinalG) {
722 std::cerr << m_className << "::ElectronCollision:\n"
723 << " Requested electron energy (" << e << " eV) exceeds the "
724 << "current energy range (" << m_eFinalG << " eV).\n"
725 << " Increasing energy range to " << 1.05 * e << " eV.\n";
726 SetMaxElectronEnergy(1.05 * e);
727 } else if (e <= 0.) {
728 std::cerr << m_className << "::ElectronCollision:\n"
729 << " Electron energy must be greater than zero.\n";
730 return false;
731 }
732
733 if (m_isChanged) {
734 if (!UpdateTransportParameters()) {
735 std::cerr << m_className << "::ElectronCollision:\n"
736 << " Error calculating the collision rates table.\n";
737 return false;
738 }
739 m_isChanged = false;
740 }
741
742 // Energy loss
743 double loss = 0.;
744 // Sample the scattering process.
745 if (band >= 0 && band < m_nValleysX) {
746 // X valley
747 // Get the energy interval.
748 int iE = int(e / m_eStepXL);
749 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
750 if (iE < 0) iE = 0;
751 // Select the scattering process.
752 const double r = RndmUniform();
753 if (r <= m_cfElectronsX[iE][0]) {
754 level = 0;
755 } else if (r >= m_cfElectronsX[iE][m_nLevelsX - 1]) {
756 level = m_nLevelsX - 1;
757 } else {
758 const auto begin = m_cfElectronsX[iE].cbegin();
759 level = std::lower_bound(begin, begin + m_nLevelsX, r) - begin;
760 }
761
762 // Get the collision type.
763 type = m_scatTypeElectronsX[level];
764 // Fill the collision counters.
765 ++m_nCollElectronDetailed[level];
766 ++m_nCollElectronBand[band];
767 if (type == ElectronCollisionTypeAcousticPhonon) {
768 ++m_nCollElectronAcoustic;
769 } else if (type == ElectronCollisionTypeIntervalleyG) {
770 // Intervalley scattering (g type)
771 ++m_nCollElectronIntervalley;
772 // Final valley is in opposite direction.
773 switch (band) {
774 case 0:
775 band = 1;
776 break;
777 case 1:
778 band = 0;
779 break;
780 case 2:
781 band = 3;
782 break;
783 case 3:
784 band = 2;
785 break;
786 case 4:
787 band = 5;
788 break;
789 case 5:
790 band = 4;
791 break;
792 default:
793 break;
794 }
795 } else if (type == ElectronCollisionTypeIntervalleyF) {
796 // Intervalley scattering (f type)
797 ++m_nCollElectronIntervalley;
798 // Final valley is perpendicular.
799 switch (band) {
800 case 0:
801 case 1:
802 band = int(RndmUniform() * 4) + 2;
803 break;
804 case 2:
805 case 3:
806 band = int(RndmUniform() * 4);
807 if (band > 1) band += 2;
808 break;
809 case 4:
810 case 5:
811 band = int(RndmUniform() * 4);
812 break;
813 }
814 } else if (type == ElectronCollisionTypeInterbandXL) {
815 // XL scattering
816 ++m_nCollElectronIntervalley;
817 // Final valley is in L band.
818 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
819 if (band >= m_nValleysX + m_nValleysL)
820 band = m_nValleysX + m_nValleysL - 1;
821 } else if (type == ElectronCollisionTypeInterbandXG) {
822 ++m_nCollElectronIntervalley;
823 band = m_nValleysX + m_nValleysL;
824 } else if (type == ElectronCollisionTypeImpurity) {
825 ++m_nCollElectronImpurity;
826 } else if (type == ElectronCollisionTypeIonisation) {
827 ++m_nCollElectronIonisation;
828 } else {
829 std::cerr << m_className << "::ElectronCollision:\n";
830 std::cerr << " Unexpected collision type (" << type << ").\n";
831 }
832
833 // Get the energy loss.
834 loss = m_energyLossElectronsX[level];
835
836 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
837 // L valley
838 // Get the energy interval.
839 int iE = int(e / m_eStepXL);
840 if (iE >= nEnergyStepsXL) iE = nEnergyStepsXL - 1;
841 if (iE < m_ieMinL) iE = m_ieMinL;
842 // Select the scattering process.
843 const double r = RndmUniform();
844 if (r <= m_cfElectronsL[iE][0]) {
845 level = 0;
846 } else if (r >= m_cfElectronsL[iE][m_nLevelsL - 1]) {
847 level = m_nLevelsL - 1;
848 } else {
849 const auto begin = m_cfElectronsL[iE].cbegin();
850 level = std::lower_bound(begin, begin + m_nLevelsL, r) - begin;
851 }
852
853 // Get the collision type.
854 type = m_scatTypeElectronsL[level];
855 // Fill the collision counters.
856 ++m_nCollElectronDetailed[m_nLevelsX + level];
857 ++m_nCollElectronBand[band];
858 if (type == ElectronCollisionTypeAcousticPhonon) {
859 ++m_nCollElectronAcoustic;
860 } else if (type == ElectronCollisionTypeOpticalPhonon) {
861 ++m_nCollElectronOptical;
862 } else if (type == ElectronCollisionTypeIntervalleyG ||
863 type == ElectronCollisionTypeIntervalleyF) {
864 // Equivalent intervalley scattering
865 ++m_nCollElectronIntervalley;
866 // Randomise the final valley.
867 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
868 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysX + m_nValleysL;
869 } else if (type == ElectronCollisionTypeInterbandXL) {
870 // LX scattering
871 ++m_nCollElectronIntervalley;
872 // Randomise the final valley.
873 band = int(RndmUniform() * m_nValleysX);
874 if (band >= m_nValleysX) band = m_nValleysX - 1;
875 } else if (type == ElectronCollisionTypeInterbandLG) {
876 // LG scattering
877 ++m_nCollElectronIntervalley;
878 band = m_nValleysX + m_nValleysL;
879 } else if (type == ElectronCollisionTypeImpurity) {
880 ++m_nCollElectronImpurity;
881 } else if (type == ElectronCollisionTypeIonisation) {
882 ++m_nCollElectronIonisation;
883 } else {
884 std::cerr << m_className << "::ElectronCollision:\n"
885 << " Unexpected collision type (" << type << ").\n";
886 }
887
888 // Get the energy loss.
889 loss = m_energyLossElectronsL[level];
890 } else if (band == m_nValleysX + m_nValleysL) {
891 // Higher bands
892 // Get the energy interval.
893 int iE = int(e / m_eStepG);
894 if (iE >= nEnergyStepsG) iE = nEnergyStepsG - 1;
895 if (iE < m_ieMinG) iE = m_ieMinG;
896 // Select the scattering process.
897 const double r = RndmUniform();
898 if (r <= m_cfElectronsG[iE][0]) {
899 level = 0;
900 } else if (r >= m_cfElectronsG[iE][m_nLevelsG - 1]) {
901 level = m_nLevelsG - 1;
902 } else {
903 const auto begin = m_cfElectronsG[iE].cbegin();
904 level = std::lower_bound(begin, begin + m_nLevelsG, r) - begin;
905 }
906
907 // Get the collision type.
908 type = m_scatTypeElectronsG[level];
909 // Fill the collision counters.
910 ++m_nCollElectronDetailed[m_nLevelsX + m_nLevelsL + level];
911 ++m_nCollElectronBand[band];
912 if (type == ElectronCollisionTypeAcousticPhonon) {
913 ++m_nCollElectronAcoustic;
914 } else if (type == ElectronCollisionTypeOpticalPhonon) {
915 ++m_nCollElectronOptical;
916 } else if (type == ElectronCollisionTypeIntervalleyG ||
917 type == ElectronCollisionTypeIntervalleyF) {
918 // Equivalent intervalley scattering
919 ++m_nCollElectronIntervalley;
920 } else if (type == ElectronCollisionTypeInterbandXG) {
921 // GX scattering
922 ++m_nCollElectronIntervalley;
923 // Randomise the final valley.
924 band = int(RndmUniform() * m_nValleysX);
925 if (band >= m_nValleysX) band = m_nValleysX - 1;
926 } else if (type == ElectronCollisionTypeInterbandLG) {
927 // GL scattering
928 ++m_nCollElectronIntervalley;
929 // Randomise the final valley.
930 band = m_nValleysX + int(RndmUniform() * m_nValleysL);
931 if (band >= m_nValleysX + m_nValleysL)
932 band = m_nValleysX + m_nValleysL - 1;
933 } else if (type == ElectronCollisionTypeIonisation) {
934 ++m_nCollElectronIonisation;
935 } else {
936 std::cerr << m_className << "::ElectronCollision:\n"
937 << " Unexpected collision type (" << type << ").\n";
938 }
939
940 // Get the energy loss.
941 loss = m_energyLossElectronsG[level];
942 } else {
943 std::cerr << m_className << "::ElectronCollision:\n"
944 << " Band index (" << band << ") out of range.\n";
945 return false;
946 }
947
948 // Secondaries
949 ndxc = 0;
950 // Ionising collision
951 if (type == ElectronCollisionTypeIonisation) {
952 double ee = 0., eh = 0.;
953 ComputeSecondaries(e, ee, eh);
954 loss = ee + eh + m_bandGap;
955 // Add the secondary electron.
956 secondaries.emplace_back(std::make_pair(Particle::Electron, ee));
957 // Add the hole.
958 secondaries.emplace_back(std::make_pair(Particle::Hole, eh));
959 }
960
961 if (e < loss) loss = e - 0.00001;
962 // Update the energy.
963 e1 = e - loss;
964 if (e1 < Small) e1 = Small;
965
966 // Update the momentum.
967 if (band >= 0 && band < m_nValleysX) {
968 // X valleys
969 double pstar = sqrt(2. * ElectronMass * e1);
970 if (m_nonParabolic) {
971 const double alpha = m_alphaX;
972 pstar *= sqrt(1. + alpha * e1);
973 }
974
975 const double ctheta = 1. - 2. * RndmUniform();
976 const double stheta = sqrt(1. - ctheta * ctheta);
977 const double phi = TwoPi * RndmUniform();
978
979 if (m_anisotropic) {
980 const double pl = pstar * sqrt(m_mLongX);
981 const double pt = pstar * sqrt(m_mTransX);
982 switch (band) {
983 case 0:
984 case 1:
985 // 100
986 px = pl * ctheta;
987 py = pt * cos(phi) * stheta;
988 pz = pt * sin(phi) * stheta;
989 break;
990 case 2:
991 case 3:
992 // 010
993 px = pt * sin(phi) * stheta;
994 py = pl * ctheta;
995 pz = pt * cos(phi) * stheta;
996 break;
997 case 4:
998 case 5:
999 // 001
1000 px = pt * cos(phi) * stheta;
1001 py = pt * sin(phi) * stheta;
1002 pz = pl * ctheta;
1003 break;
1004 default:
1005 return false;
1006 }
1007 } else {
1008 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
1009 px = pstar * cos(phi) * stheta;
1010 py = pstar * sin(phi) * stheta;
1011 pz = pstar * ctheta;
1012 }
1013 return true;
1014
1015 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
1016 // L valleys
1017 double pstar = sqrt(2. * ElectronMass * (e1 - m_eMinL));
1018 if (m_nonParabolic) {
1019 const double alpha = m_alphaL;
1020 pstar *= sqrt(1. + alpha * (e1 - m_eMinL));
1021 }
1022 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
1023 RndmDirection(px, py, pz, pstar);
1024 return true;
1025 }
1026 const double pstar = sqrt(2. * ElectronMass * e1);
1027 RndmDirection(px, py, pz, pstar);
1028 return true;
1029}
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:128
DoubleAc cos(const DoubleAc &f)
Definition DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition DoubleAc.cpp:384

◆ 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; }

Referenced by ElectronVelocity().

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

171 {
172 alpha = 0.;
173 if (m_isChanged) {
174 if (!UpdateTransportParameters()) {
175 std::cerr << m_className << "::ElectronTownsend:\n"
176 << " Error calculating the transport parameters.\n";
177 return false;
178 }
179 m_isChanged = false;
180 }
181
182 if (!m_eAlp.empty()) {
183 // Interpolation in user table.
184 return Medium::ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
185 }
186
187 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
188 alpha = ElectronAlpha(emag);
189 return true;
190}
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition Medium.hh:582
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:468

◆ 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 // Calculate the mobility.
155 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
156 const double mu = -ElectronMobility(emag);
157
158 if (fabs(bx) < Small && fabs(by) < Small && fabs(bz) < Small) {
159 vx = mu * ex;
160 vy = mu * ey;
161 vz = mu * ez;
162 } else {
163 Langevin(ex, ey, ez, bx, by, bz, mu, m_eHallFactor * mu, vx, vy, vz);
164 }
165 return true;
166}
double ElectronMobility() override
Low-field mobility [cm2 V-1 ns-1].
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:444
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition Medium.hh:577
static void Langevin(const double ex, const double ey, const double ez, double bx, double by, double bz, const double mu, double &vx, double &vy, double &vz)
Definition Medium.cc:301
DoubleAc fabs(const DoubleAc &f)
Definition DoubleAc.h:615

Referenced by ElectronAttachment().

◆ EnableAnisotropy()

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

Definition at line 112 of file MediumSilicon.hh.

112{ m_anisotropic = on; }

◆ EnableFullBandDensityOfStates()

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

Definition at line 109 of file MediumSilicon.hh.

109 {
110 m_fullBandDos = on;
111 }

◆ EnableNonParabolicity()

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

Definition at line 108 of file MediumSilicon.hh.

108{ m_nonParabolic = on; }

◆ EnableScatteringRateOutput()

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

Definition at line 107 of file MediumSilicon.hh.

107{ m_cfOutput = on; }

◆ GetConductionBandDensityOfStates()

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

Definition at line 2591 of file MediumSilicon.cc.

2592 {
2593 if (band < 0) {
2594 int iE = int(e / m_eStepDos);
2595 const int nPoints = m_fbDosConduction.size();
2596 if (iE >= nPoints || iE < 0) {
2597 return 0.;
2598 } else if (iE == nPoints - 1) {
2599 return m_fbDosConduction[nPoints - 1];
2600 }
2601
2602 const double dos = m_fbDosConduction[iE] +
2603 (m_fbDosConduction[iE + 1] - m_fbDosConduction[iE]) *
2604 (e / m_eStepDos - iE);
2605 return dos * 1.e21;
2606
2607 } else if (band < m_nValleysX) {
2608 // X valleys
2609 if (e <= 0.) return 0.;
2610 // Density-of-states effective mass (cube)
2611 const double md3 = pow(ElectronMass, 3) * m_mLongX * m_mTransX * m_mTransX;
2612
2613 if (m_fullBandDos) {
2614 if (e < m_eMinL) {
2615 return GetConductionBandDensityOfStates(e, -1) / m_nValleysX;
2616 } else if (e < m_eMinG) {
2617 // Subtract the fraction of the full-band density of states
2618 // attributed to the L valleys.
2619 const double dosX =
2621 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL;
2622 return dosX / m_nValleysX;
2623 } else {
2624 // Subtract the fraction of the full-band density of states
2625 // attributed to the L valleys and the higher bands.
2626 const double dosX =
2628 GetConductionBandDensityOfStates(e, m_nValleysX) * m_nValleysL -
2629 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
2630 if (dosX <= 0.) return 0.;
2631 return dosX / m_nValleysX;
2632 }
2633 } else if (m_nonParabolic) {
2634 const double alpha = m_alphaX;
2635 return sqrt(md3 * e * (1. + alpha * e) / 2.) * (1. + 2 * alpha * e) /
2636 (Pi2 * pow(HbarC, 3.));
2637 } else {
2638 return sqrt(md3 * e / 2.) / (Pi2 * pow(HbarC, 3.));
2639 }
2640 } else if (band < m_nValleysX + m_nValleysL) {
2641 // L valleys
2642 if (e <= m_eMinL) return 0.;
2643
2644 // Density-of-states effective mass (cube)
2645 const double md3 = pow(ElectronMass, 3) * m_mLongL * m_mTransL * m_mTransL;
2646 // Non-parabolicity parameter
2647 const double alpha = m_alphaL;
2648
2649 if (m_fullBandDos) {
2650 // Energy up to which the non-parabolic approximation is used.
2651 const double ej = m_eMinL + 0.5;
2652 if (e <= ej) {
2653 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
2654 (1. + 2 * alpha * (e - m_eMinL)) /
2655 (Sqrt2 * Pi2 * pow(HbarC, 3.));
2656 } else {
2657 // Fraction of full-band density of states attributed to L valleys
2658 double fL = sqrt(md3 * (ej - m_eMinL) * (1. + alpha * (ej - m_eMinL))) *
2659 (1. + 2 * alpha * (ej - m_eMinL)) /
2660 (Sqrt2 * Pi2 * pow(HbarC, 3.));
2661 fL = m_nValleysL * fL / GetConductionBandDensityOfStates(ej, -1);
2662
2663 double dosXL = GetConductionBandDensityOfStates(e, -1);
2664 if (e > m_eMinG) {
2665 dosXL -=
2666 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
2667 }
2668 if (dosXL <= 0.) return 0.;
2669 return fL * dosXL / 8.;
2670 }
2671 } else if (m_nonParabolic) {
2672 return sqrt(md3 * (e - m_eMinL) * (1. + alpha * (e - m_eMinL))) *
2673 (1. + 2 * alpha * (e - m_eMinL)) / (Sqrt2 * Pi2 * pow(HbarC, 3.));
2674 } else {
2675 return sqrt(md3 * (e - m_eMinL) / 2.) / (Pi2 * pow(HbarC, 3.));
2676 }
2677 } else if (band == m_nValleysX + m_nValleysL) {
2678 // Higher bands
2679 const double ej = 2.7;
2680 if (m_eMinG >= ej) {
2681 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2682 << " Cannot determine higher band density-of-states.\n"
2683 << " Program bug. Check offset energy!\n";
2684 return 0.;
2685 }
2686 if (e < m_eMinG) {
2687 return 0.;
2688 } else if (e < ej) {
2689 // Coexistence of XL and higher bands.
2690 const double dj = GetConductionBandDensityOfStates(ej, -1);
2691 // Assume linear increase of density-of-states.
2692 return dj * (e - m_eMinG) / (ej - m_eMinG);
2693 } else {
2695 }
2696 }
2697
2698 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2699 << " Band index (" << band << ") out of range.\n";
2700 return ElectronMass * sqrt(ElectronMass * e / 2.) / (Pi2 * pow(HbarC, 3.));
2701}
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 1101 of file MediumSilicon.cc.

1102 {
1103 if (i != 0) {
1104 std::cerr << m_className + "::GetDielectricFunction: Index out of range.\n";
1105 return false;
1106 }
1107
1108 // Make sure the optical data table has been loaded.
1109 if (m_opticalDataEnergies.empty()) {
1110 if (!LoadOpticalData(m_opticalDataFile)) {
1111 std::cerr << m_className << "::GetDielectricFunction:\n";
1112 std::cerr << " Optical data table could not be loaded.\n";
1113 return false;
1114 }
1115 }
1116
1117 // Make sure the requested energy is within the range of the table.
1118 const double emin = m_opticalDataEnergies.front();
1119 const double emax = m_opticalDataEnergies.back();
1120 if (e < emin || e > emax) {
1121 std::cerr << m_className << "::GetDielectricFunction:\n"
1122 << " Requested energy (" << e << " eV) "
1123 << " is outside the range of the optical data table.\n"
1124 << " " << emin << " < E [eV] < " << emax << "\n";
1125 eps1 = eps2 = 0.;
1126 return false;
1127 }
1128
1129 // Locate the requested energy in the table.
1130 const auto begin = m_opticalDataEnergies.cbegin();
1131 const auto it1 = std::upper_bound(begin, m_opticalDataEnergies.cend(), e);
1132 if (it1 == begin) {
1133 eps1 = m_opticalDataEpsilon.front().first;
1134 eps2 = m_opticalDataEpsilon.front().second;
1135 return true;
1136 }
1137 const auto it0 = std::prev(it1);
1138
1139 // Interpolate the real part of dielectric function.
1140 const double x0 = *it0;
1141 const double x1 = *it1;
1142 const double lnx0 = log(*it0);
1143 const double lnx1 = log(*it1);
1144 const double lnx = log(e);
1145 const double y0 = m_opticalDataEpsilon[it0 - begin].first;
1146 const double y1 = m_opticalDataEpsilon[it1 - begin].first;
1147 if (y0 <= 0. || y1 <= 0.) {
1148 // Use linear interpolation if one of the values is negative.
1149 eps1 = y0 + (e - x0) * (y1 - y0) / (x1 - x0);
1150 } else {
1151 // Otherwise use log-log interpolation.
1152 const double lny0 = log(y0);
1153 const double lny1 = log(y1);
1154 eps1 = lny0 + (lnx - lnx0) * (lny1 - lny0) / (lnx1 - lnx0);
1155 eps1 = exp(eps1);
1156 }
1157
1158 // Interpolate the imaginary part of dielectric function,
1159 // using log-log interpolation.
1160 const double lnz0 = log(m_opticalDataEpsilon[it0 - begin].second);
1161 const double lnz1 = log(m_opticalDataEpsilon[it1 - begin].second);
1162 eps2 = lnz0 + (lnx - lnx0) * (lnz1 - lnz0) / (lnx1 - lnx0);
1163 eps2 = exp(eps2);
1164 return true;
1165}
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 1067 of file MediumSilicon.cc.

1067 {
1068 const int nBands = m_nValleysX + m_nValleysL + 1;
1069 if (band < 0 || band >= nBands) {
1070 std::cerr << m_className << "::GetElectronBandPopulation:\n";
1071 std::cerr << " Band index (" << band << ") out of range.\n";
1072 return 0;
1073 }
1074 return m_nCollElectronBand[band];
1075}

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

665 {
666 if (e <= 0.) {
667 std::cerr << m_className << "::GetElectronCollisionRate:\n"
668 << " Electron energy must be positive.\n";
669 return 0.;
670 }
671
672 if (e > m_eFinalG) {
673 std::cerr << m_className << "::GetElectronCollisionRate:\n"
674 << " Collision rate at " << e << " eV (band " << band
675 << ") is not included in the current table.\n"
676 << " Increasing energy range to " << 1.05 * e << " eV.\n";
677 SetMaxElectronEnergy(1.05 * e);
678 }
679
680 if (m_isChanged) {
681 if (!UpdateTransportParameters()) {
682 std::cerr << m_className << "::GetElectronCollisionRate:\n"
683 << " Error calculating the collision rates table.\n";
684 return 0.;
685 }
686 m_isChanged = false;
687 }
688
689 if (band >= 0 && band < m_nValleysX) {
690 int iE = int(e / m_eStepXL);
691 if (iE >= nEnergyStepsXL)
692 iE = nEnergyStepsXL - 1;
693 else if (iE < 0)
694 iE = 0;
695 return m_cfTotElectronsX[iE];
696 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
697 int iE = int(e / m_eStepXL);
698 if (iE >= nEnergyStepsXL)
699 iE = nEnergyStepsXL - 1;
700 else if (iE < m_ieMinL)
701 iE = m_ieMinL;
702 return m_cfTotElectronsL[iE];
703 } else if (band == m_nValleysX + m_nValleysL) {
704 int iE = int(e / m_eStepG);
705 if (iE >= nEnergyStepsG)
706 iE = nEnergyStepsG - 1;
707 else if (iE < m_ieMinG)
708 iE = m_ieMinG;
709 return m_cfTotElectronsG[iE];
710 }
711
712 std::cerr << m_className << "::GetElectronCollisionRate:\n"
713 << " Band index (" << band << ") out of range.\n";
714 return 0.;
715}

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

457 {
458 // Effective masses
459 double mx = ElectronMass, my = ElectronMass, mz = ElectronMass;
460 // Energy offset
461 double e0 = 0.;
462 if (band >= 0 && band < m_nValleysX) {
463 // X valley
464 if (m_anisotropic) {
465 switch (band) {
466 case 0:
467 case 1:
468 // X 100, -100
469 mx *= m_mLongX;
470 my *= m_mTransX;
471 mz *= m_mTransX;
472 break;
473 case 2:
474 case 3:
475 // X 010, 0-10
476 mx *= m_mTransX;
477 my *= m_mLongX;
478 mz *= m_mTransX;
479 break;
480 case 4:
481 case 5:
482 // X 001, 00-1
483 mx *= m_mTransX;
484 my *= m_mTransX;
485 mz *= m_mLongX;
486 break;
487 default:
488 std::cerr << m_className << "::GetElectronEnergy:\n"
489 << " Unexpected band index " << band << "!\n";
490 break;
491 }
492 } else {
493 // Conduction effective mass
494 const double mc = 3. / (1. / m_mLongX + 2. / m_mTransX);
495 mx *= mc;
496 my *= mc;
497 mz *= mc;
498 }
499 } else if (band < m_nValleysX + m_nValleysL) {
500 // L valley, isotropic approximation
501 e0 = m_eMinL;
502 // Effective mass
503 const double mc = 3. / (1. / m_mLongL + 2. / m_mTransL);
504 mx *= mc;
505 my *= mc;
506 mz *= mc;
507 } else if (band == m_nValleysX + m_nValleysL) {
508 // Higher band(s)
509 }
510
511 if (m_nonParabolic) {
512 // Non-parabolicity parameter
513 double alpha = 0.;
514 if (band < m_nValleysX) {
515 // X valley
516 alpha = m_alphaX;
517 } else if (band < m_nValleysX + m_nValleysL) {
518 // L valley
519 alpha = m_alphaL;
520 }
521
522 const double p2 = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
523 if (alpha > 0.) {
524 const double e = 0.5 * (sqrt(1. + 4 * alpha * p2) - 1.) / alpha;
525 const double a = SpeedOfLight / (1. + 2 * alpha * e);
526 vx = a * px / mx;
527 vy = a * py / my;
528 vz = a * pz / mz;
529 return e0 + e;
530 }
531 }
532
533 const double e = 0.5 * (px * px / mx + py * py / my + pz * pz / mz);
534 vx = SpeedOfLight * px / mx;
535 vy = SpeedOfLight * py / my;
536 vz = SpeedOfLight * pz / mz;
537 return e0 + e;
538}

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

541 {
542 // If the band index is out of range, choose one at random.
543 if (band < 0 || band > m_nValleysX + m_nValleysL ||
544 (e < m_eMinL || band >= m_nValleysX) ||
545 (e < m_eMinG || band == m_nValleysX + m_nValleysL)) {
546 if (e < m_eMinL) {
547 band = int(m_nValleysX * RndmUniform());
548 if (band >= m_nValleysX) band = m_nValleysX - 1;
549 } else {
550 const double dosX = GetConductionBandDensityOfStates(e, 0);
551 const double dosL = GetConductionBandDensityOfStates(e, m_nValleysX);
552 const double dosG =
553 GetConductionBandDensityOfStates(e, m_nValleysX + m_nValleysL);
554 const double dosSum = m_nValleysX * dosX + m_nValleysL * dosL + dosG;
555 if (dosSum < Small) {
556 band = m_nValleysX + m_nValleysL;
557 } else {
558 const double r = RndmUniform() * dosSum;
559 if (r < dosX) {
560 band = int(m_nValleysX * RndmUniform());
561 if (band >= m_nValleysX) band = m_nValleysX - 1;
562 } else if (r < dosX + dosL) {
563 band = m_nValleysX + int(m_nValleysL * RndmUniform());
564 if (band >= m_nValleysX + m_nValleysL) band = m_nValleysL - 1;
565 } else {
566 band = m_nValleysX + m_nValleysL;
567 }
568 }
569 }
570 if (m_debug) {
571 std::cout << m_className << "::GetElectronMomentum:\n"
572 << " Randomised band index: " << band << "\n";
573 }
574 }
575 if (band < m_nValleysX) {
576 // X valleys
577 double pstar = sqrt(2. * ElectronMass * e);
578 if (m_nonParabolic) {
579 const double alpha = m_alphaX;
580 pstar *= sqrt(1. + alpha * e);
581 }
582
583 const double ctheta = 1. - 2. * RndmUniform();
584 const double stheta = sqrt(1. - ctheta * ctheta);
585 const double phi = TwoPi * RndmUniform();
586
587 if (m_anisotropic) {
588 const double pl = pstar * sqrt(m_mLongX);
589 const double pt = pstar * sqrt(m_mTransX);
590 switch (band) {
591 case 0:
592 case 1:
593 // 100
594 px = pl * ctheta;
595 py = pt * cos(phi) * stheta;
596 pz = pt * sin(phi) * stheta;
597 break;
598 case 2:
599 case 3:
600 // 010
601 px = pt * sin(phi) * stheta;
602 py = pl * ctheta;
603 pz = pt * cos(phi) * stheta;
604 break;
605 case 4:
606 case 5:
607 // 001
608 px = pt * cos(phi) * stheta;
609 py = pt * sin(phi) * stheta;
610 pz = pl * ctheta;
611 break;
612 default:
613 // Other band; should not occur.
614 std::cerr << m_className << "::GetElectronMomentum:\n"
615 << " Unexpected band index (" << band << ").\n";
616 px = pstar * stheta * cos(phi);
617 py = pstar * stheta * sin(phi);
618 pz = pstar * ctheta;
619 break;
620 }
621 } else {
622 pstar *= sqrt(3. / (1. / m_mLongX + 2. / m_mTransX));
623 px = pstar * cos(phi) * stheta;
624 py = pstar * sin(phi) * stheta;
625 pz = pstar * ctheta;
626 }
627 } else if (band < m_nValleysX + m_nValleysL) {
628 // L valleys
629 double pstar = sqrt(2. * ElectronMass * (e - m_eMinL));
630 if (m_nonParabolic) {
631 const double alpha = m_alphaL;
632 pstar *= sqrt(1. + alpha * (e - m_eMinL));
633 }
634 pstar *= sqrt(3. / (1. / m_mLongL + 2. / m_mTransL));
635 RndmDirection(px, py, pz, pstar);
636 } else if (band == m_nValleysX + m_nValleysL) {
637 // Higher band
638 const double pstar = sqrt(2. * ElectronMass * e);
639 RndmDirection(px, py, pz, pstar);
640 }
641}

◆ GetElectronNullCollisionRate()

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

Null-collision rate [ns-1].

Reimplemented from Garfield::Medium.

Definition at line 643 of file MediumSilicon.cc.

643 {
644 if (m_isChanged) {
645 if (!UpdateTransportParameters()) {
646 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
647 << " Error calculating the collision rates table.\n";
648 return 0.;
649 }
650 m_isChanged = false;
651 }
652
653 if (band >= 0 && band < m_nValleysX) {
654 return m_cfNullElectronsX;
655 } else if (band >= m_nValleysX && band < m_nValleysX + m_nValleysL) {
656 return m_cfNullElectronsL;
657 } else if (band == m_nValleysX + m_nValleysL) {
658 return m_cfNullElectronsG;
659 }
660 std::cerr << m_className << "::GetElectronNullCollisionRate:\n"
661 << " Band index (" << band << ") out of range.\n";
662 return 0.;
663}

◆ GetMaxElectronEnergy()

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

Definition at line 101 of file MediumSilicon.hh.

101{ return m_eFinalG; }

◆ GetNumberOfElectronBands()

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

Definition at line 1063 of file MediumSilicon.cc.

1063 {
1064 return m_nValleysX + m_nValleysL + 1;
1065}

◆ GetNumberOfElectronCollisions() [1/2]

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

Definition at line 1042 of file MediumSilicon.cc.

1042 {
1043 return m_nCollElectronAcoustic + m_nCollElectronOptical +
1044 m_nCollElectronIntervalley + m_nCollElectronImpurity +
1045 m_nCollElectronIonisation;
1046}

◆ GetNumberOfElectronCollisions() [2/2]

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

Definition at line 1052 of file MediumSilicon.cc.

1053 {
1054 const unsigned int nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1055 if (level >= nLevels) {
1056 std::cerr << m_className << "::GetNumberOfElectronCollisions:\n"
1057 << " Scattering rate term (" << level << ") does not exist.\n";
1058 return 0;
1059 }
1060 return m_nCollElectronDetailed[level];
1061}

◆ GetNumberOfLevels()

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

Definition at line 1048 of file MediumSilicon.cc.

1048 {
1049 return m_nLevelsX + m_nLevelsL + m_nLevelsG;
1050}

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

1078 {
1079 if (i != 0) {
1080 std::cerr << m_className << "::GetOpticalDataRange: Index out of range.\n";
1081 }
1082
1083 // Make sure the optical data table has been loaded.
1084 if (m_opticalDataEnergies.empty()) {
1085 if (!LoadOpticalData(m_opticalDataFile)) {
1086 std::cerr << m_className << "::GetOpticalDataRange:\n"
1087 << " Optical data table could not be loaded.\n";
1088 return false;
1089 }
1090 }
1091
1092 emin = m_opticalDataEnergies.front();
1093 emax = m_opticalDataEnergies.back();
1094 if (m_debug) {
1095 std::cout << m_className << "::GetOpticalDataRange:\n"
1096 << " " << emin << " < E [eV] < " << emax << "\n";
1097 }
1098 return true;
1099}

◆ GetValenceBandDensityOfStates()

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

Definition at line 2703 of file MediumSilicon.cc.

2704 {
2705 if (band <= 0) {
2706 // Total (full-band) density of states.
2707 const int nPoints = m_fbDosValence.size();
2708 int iE = int(e / m_eStepDos);
2709 if (iE >= nPoints || iE < 0) {
2710 return 0.;
2711 } else if (iE == nPoints - 1) {
2712 return m_fbDosValence[nPoints - 1];
2713 }
2714
2715 const double dos =
2716 m_fbDosValence[iE] +
2717 (m_fbDosValence[iE + 1] - m_fbDosValence[iE]) * (e / m_eStepDos - iE);
2718 return dos * 1.e21;
2719 }
2720
2721 std::cerr << m_className << "::GetConductionBandDensityOfStates:\n"
2722 << " Band index (" << band << ") out of range.\n";
2723 return 0.;
2724}

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

289 {
290 eta = 0.;
291 if (m_isChanged) {
292 if (!UpdateTransportParameters()) {
293 std::cerr << m_className << "::HoleAttachment:\n"
294 << " Error calculating the transport parameters.\n";
295 return false;
296 }
297 m_isChanged = false;
298 }
299
300 if (!m_hAtt.empty()) {
301 // Interpolation in user table.
302 return Medium::HoleAttachment(ex, ey, ez, bx, by, bz, eta);
303 }
304
305 switch (m_trappingModel) {
306 case 0:
307 eta = m_hTrapCs * m_hTrapDensity;
308 break;
309 case 1:
310 double vx, vy, vz;
311 HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
312 eta = m_hTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
313 if (eta > 0.) eta = -1. / eta;
314 break;
315 default:
316 std::cerr << m_className << "::HoleAttachment: Unknown model. Bug!\n";
317 return false;
318 }
319
320 return true;
321}
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:595
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:612

◆ 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; }

Referenced by HoleVelocity().

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

265 {
266 alpha = 0.;
267 if (m_isChanged) {
268 if (!UpdateTransportParameters()) {
269 std::cerr << m_className << "::HoleTownsend:\n"
270 << " Error calculating the transport parameters.\n";
271 return false;
272 }
273 m_isChanged = false;
274 }
275
276 if (!m_hAlp.empty()) {
277 // Interpolation in user table.
278 return Medium::HoleTownsend(ex, ey, ez, bx, by, bz, alpha);
279 }
280
281 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
282 alpha = HoleAlpha(emag);
283 return true;
284}
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:599
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition Medium.hh:594

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

232 {
233 vx = vy = vz = 0.;
234 if (m_isChanged) {
235 if (!UpdateTransportParameters()) {
236 std::cerr << m_className << "::HoleVelocity:\n"
237 << " Error calculating the transport parameters.\n";
238 return false;
239 }
240 m_isChanged = false;
241 }
242
243 if (!m_hVelE.empty()) {
244 // Interpolation in user table.
245 return Medium::HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
246 }
247
248 // Calculate the mobility.
249 const double emag = sqrt(ex * ex + ey * ey + ez * ez);
250 const double mu = HoleMobility(emag);
251
252 if (fabs(bx) < Small && fabs(by) < Small && fabs(bz) < Small) {
253 vx = mu * ex;
254 vy = mu * ey;
255 vz = mu * ez;
256 } else {
257 Langevin(ex, ey, ez, bx, by, bz, mu, m_hHallFactor * mu, vx, vy, vz);
258 }
259 return true;
260}
double HoleMobility() override
Low-field mobility [cm2 V-1 ns-1].
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:578
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition Medium.hh:589

Referenced by HoleAttachment().

◆ Initialise()

bool Garfield::MediumSilicon::Initialise ( )

Definition at line 1167 of file MediumSilicon.cc.

1167 {
1168 if (!m_isChanged) {
1169 if (m_debug) {
1170 std::cerr << m_className << "::Initialise: Nothing changed.\n";
1171 }
1172 return true;
1173 }
1174 if (!UpdateTransportParameters()) {
1175 std::cerr << m_className << "::Initialise:\n Error preparing "
1176 << "transport parameters/calculating collision rates.\n";
1177 return false;
1178 }
1179 return true;
1180}

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

1031 {
1032 m_nCollElectronAcoustic = m_nCollElectronOptical = 0;
1033 m_nCollElectronIntervalley = 0;
1034 m_nCollElectronImpurity = 0;
1035 m_nCollElectronIonisation = 0;
1036 const auto nLevels = m_nLevelsX + m_nLevelsL + m_nLevelsG;
1037 m_nCollElectronDetailed.assign(nLevels, 0);
1038 const auto nBands = m_nValleysX + m_nValleysL + 1;
1039 m_nCollElectronBand.assign(nBands, 0);
1040}

◆ SetDiffusionScaling()

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

Apply a scaling factor to the diffusion coefficients.

Definition at line 97 of file MediumSilicon.hh.

97{ 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 360 of file MediumSilicon.cc.

360 {
361 m_dopingMobilityModel = DopingMobility::Masetti;
362 m_hasUserMobility = false;
363 m_isChanged = true;
364}

◆ SetDopingMobilityModelMinimos()

void Garfield::MediumSilicon::SetDopingMobilityModelMinimos ( )

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

Definition at line 354 of file MediumSilicon.cc.

354 {
355 m_dopingMobilityModel = DopingMobility::Minimos;
356 m_hasUserMobility = false;
357 m_isChanged = true;
358}

◆ SetHighFieldMobilityModelCanali()

void Garfield::MediumSilicon::SetHighFieldMobilityModelCanali ( )

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

Definition at line 404 of file MediumSilicon.cc.

404 {
405 m_highFieldMobilityModel = HighFieldMobility::Canali;
406 m_isChanged = true;
407}

◆ SetHighFieldMobilityModelConstant()

void Garfield::MediumSilicon::SetHighFieldMobilityModelConstant ( )

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

Definition at line 414 of file MediumSilicon.cc.

414 {
415 m_highFieldMobilityModel = HighFieldMobility::Constant;
416}

◆ SetHighFieldMobilityModelMinimos()

void Garfield::MediumSilicon::SetHighFieldMobilityModelMinimos ( )

Parameterize the high-field mobility using the Minimos model.

Definition at line 399 of file MediumSilicon.cc.

399 {
400 m_highFieldMobilityModel = HighFieldMobility::Minimos;
401 m_isChanged = true;
402}

◆ SetHighFieldMobilityModelReggiani()

void Garfield::MediumSilicon::SetHighFieldMobilityModelReggiani ( )

Parameterize the high-field mobility using the Reggiani model.

Definition at line 409 of file MediumSilicon.cc.

409 {
410 m_highFieldMobilityModel = HighFieldMobility::Reggiani;
411 m_isChanged = true;
412}

◆ SetImpactIonisationModelGrant()

void Garfield::MediumSilicon::SetImpactIonisationModelGrant ( )

Calculate α using the Grant model.

Definition at line 423 of file MediumSilicon.cc.

423 {
424 m_impactIonisationModel = ImpactIonisation::Grant;
425 m_isChanged = true;
426}

◆ SetImpactIonisationModelMassey()

void Garfield::MediumSilicon::SetImpactIonisationModelMassey ( )

Calculate α using the Massey model.

Definition at line 428 of file MediumSilicon.cc.

428 {
429 m_impactIonisationModel = ImpactIonisation::Massey;
430 m_isChanged = true;
431}

◆ SetImpactIonisationModelOkutoCrowell()

void Garfield::MediumSilicon::SetImpactIonisationModelOkutoCrowell ( )

Calculate α using the Okuto-Crowell model.

Definition at line 433 of file MediumSilicon.cc.

433 {
434 m_impactIonisationModel = ImpactIonisation::Okuto;
435 m_isChanged = true;
436}

◆ SetImpactIonisationModelVanOverstraetenDeMan()

void Garfield::MediumSilicon::SetImpactIonisationModelVanOverstraetenDeMan ( )

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

Definition at line 418 of file MediumSilicon.cc.

418 {
419 m_impactIonisationModel = ImpactIonisation::VanOverstraeten;
420 m_isChanged = true;
421}

◆ SetLatticeMobilityModelMinimos()

void Garfield::MediumSilicon::SetLatticeMobilityModelMinimos ( )

Calculate the lattice mobility using the Minimos model.

Definition at line 336 of file MediumSilicon.cc.

336 {
337 m_latticeMobilityModel = LatticeMobility::Minimos;
338 m_hasUserMobility = false;
339 m_isChanged = true;
340}

◆ SetLatticeMobilityModelReggiani()

void Garfield::MediumSilicon::SetLatticeMobilityModelReggiani ( )

Calculate the lattice mobility using the Reggiani model.

Definition at line 348 of file MediumSilicon.cc.

348 {
349 m_latticeMobilityModel = LatticeMobility::Reggiani;
350 m_hasUserMobility = false;
351 m_isChanged = true;
352}

◆ SetLatticeMobilityModelSentaurus()

void Garfield::MediumSilicon::SetLatticeMobilityModelSentaurus ( )

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

Definition at line 342 of file MediumSilicon.cc.

342 {
343 m_latticeMobilityModel = LatticeMobility::Sentaurus;
344 m_hasUserMobility = false;
345 m_isChanged = true;
346}

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

323 {
324 if (mue <= 0. || muh <= 0.) {
325 std::cerr << m_className << "::SetLowFieldMobility:\n"
326 << " Mobility must be greater than zero.\n";
327 return;
328 }
329
330 m_eMobility = mue;
331 m_hMobility = muh;
332 m_hasUserMobility = true;
333 m_isChanged = true;
334}

◆ SetMaxElectronEnergy()

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

Definition at line 438 of file MediumSilicon.cc.

438 {
439 if (e <= m_eMinG + Small) {
440 std::cerr << m_className << "::SetMaxElectronEnergy:\n"
441 << " Requested upper electron energy limit (" << e
442 << " eV) is too small.\n";
443 return false;
444 }
445
446 m_eFinalG = e;
447 // Determine the energy interval size.
448 m_eStepG = m_eFinalG / nEnergyStepsG;
449
450 m_isChanged = true;
451
452 return true;
453}

Referenced by ElectronCollision(), and GetElectronCollisionRate().

◆ SetSaturationVelocity()

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

Specify the saturation velocities of electrons and holes.

Definition at line 366 of file MediumSilicon.cc.

367 {
368 if (vsate <= 0. || vsath <= 0.) {
369 std::cout << m_className << "::SetSaturationVelocity:\n"
370 << " Restoring default values.\n";
371 m_hasUserSaturationVelocity = false;
372 } else {
373 m_eSatVel = vsate;
374 m_hSatVel = vsath;
375 m_hasUserSaturationVelocity = true;
376 }
377
378 m_isChanged = true;
379}

◆ SetSaturationVelocityModelCanali()

void Garfield::MediumSilicon::SetSaturationVelocityModelCanali ( )

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

Definition at line 387 of file MediumSilicon.cc.

387 {
388 m_saturationVelocityModel = SaturationVelocity::Canali;
389 m_hasUserSaturationVelocity = false;
390 m_isChanged = true;
391}

◆ SetSaturationVelocityModelMinimos()

void Garfield::MediumSilicon::SetSaturationVelocityModelMinimos ( )

Calculate the saturation velocities using the Minimos model.

Definition at line 381 of file MediumSilicon.cc.

381 {
382 m_saturationVelocityModel = SaturationVelocity::Minimos;
383 m_hasUserSaturationVelocity = false;
384 m_isChanged = true;
385}

◆ SetSaturationVelocityModelReggiani()

void Garfield::MediumSilicon::SetSaturationVelocityModelReggiani ( )

Calculate the saturation velocities using the Reggiani model.

Definition at line 393 of file MediumSilicon.cc.

393 {
394 m_saturationVelocityModel = SaturationVelocity::Reggiani;
395 m_hasUserSaturationVelocity = false;
396 m_isChanged = true;
397}

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