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::Medium Class Reference

Abstract base class for media. More...

#include <Medium.hh>

+ Inheritance diagram for Garfield::Medium:

Public Member Functions

 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 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 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 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)
 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)
 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 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 ()
 

Protected Member Functions

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

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

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

static int m_idCounter = -1
 

Detailed Description

Abstract base class for media.

Definition at line 16 of file Medium.hh.

Constructor & Destructor Documentation

◆ Medium()

Garfield::Medium::Medium ( )

Constructor.

Definition at line 61 of file Medium.cc.

61 : m_id(++m_idCounter) {
62 // Initialise the tables.
63 m_bFields.assign(1, 0.);
64 m_bAngles.assign(1, HalfPi);
65
66 // Set the default grid.
67 SetFieldGrid(100., 100000., 20, true, 0., 0., 1, HalfPi, HalfPi, 1);
68}
std::vector< double > m_bFields
Definition Medium.hh:573
static int m_idCounter
Definition Medium.hh:531
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.
Definition Medium.cc:791
std::vector< double > m_bAngles
Definition Medium.hh:574

Referenced by Garfield::MediumCdTe::MediumCdTe(), Garfield::MediumConductor::MediumConductor(), Garfield::MediumDiamond::MediumDiamond(), Garfield::MediumGaAs::MediumGaAs(), Garfield::MediumGaN::MediumGaN(), Garfield::MediumGas::MediumGas(), Garfield::MediumPlastic::MediumPlastic(), and Garfield::MediumSilicon::MediumSilicon().

◆ ~Medium()

Garfield::Medium::~Medium ( )
virtual

Destructor.

Definition at line 70 of file Medium.cc.

70{}

Member Function Documentation

◆ Alpha()

bool Garfield::Medium::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
protected

Definition at line 413 of file Medium.cc.

418 {
419
420 alpha = 0.;
421 if (tab.empty()) return false;
422
423 // Compute the magnitude of the electric field.
424 const double e = sqrt(ex * ex + ey * ey + ez * ez);
425 const double e0 = ScaleElectricField(e);
426 if (e < Small || e0 < Small) return true;
427
428 // Compute the magnitude of the magnetic field.
429 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
430 // Compute the angle between B field and E field.
431 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
432
433 // Interpolate.
434 if (e0 < m_eFields[thr]) intp = 1;
435 if (!Interpolate(e0, b, ebang, tab, alpha, intp, extr)) alpha = -30.;
436 if (alpha < -20.) {
437 alpha = 0.;
438 } else {
439 alpha = exp(alpha);
440 }
441 return true;
442}
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
Definition Medium.cc:1314
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
Definition Medium.cc:1300
virtual double ScaleElectricField(const double e) const
Definition Medium.hh:499
std::vector< double > m_eFields
Definition Medium.hh:572
DoubleAc exp(const DoubleAc &f)
Definition DoubleAc.cpp:377
DoubleAc sqrt(const DoubleAc &f)
Definition DoubleAc.cpp:314

Referenced by ElectronAttachment(), ElectronTownsend(), HoleAttachment(), HoleTownsend(), and IonDissociation().

◆ Clone() [1/2]

void Garfield::Medium::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 )
protected

Definition at line 1018 of file Medium.cc.

1024 {
1025 if (m_debug) {
1026 std::cout << m_className << "::Clone: Copying " << lbl << " to new grid.\n";
1027 }
1028
1029 if (tab.empty()) {
1030 if (m_debug) std::cout << m_className << "::Clone: Table is empty.\n";
1031 return;
1032 }
1033 // Get the dimensions of the new grid.
1034 const auto nE = efields.size();
1035 const auto nB = bfields.size();
1036 const auto nA = angles.size();
1037
1038 // Create a temporary table to store the values at the new grid points.
1039 std::vector<std::vector<std::vector<double> > > tabClone;
1040 Init(nE, nB, nA, tabClone, init);
1041
1042 // Fill the temporary table.
1043 for (size_t i = 0; i < nE; ++i) {
1044 const double e = efields[i];
1045 for (size_t j = 0; j < nB; ++j) {
1046 const double b = bfields[j];
1047 for (size_t k = 0; k < nA; ++k) {
1048 const double a = angles[k];
1049 double val = 0.;
1050 if (!Interpolate(e, b, a, tab, val, intp, extr)) {
1051 std::cerr << m_className << "::Clone:\n"
1052 << " Interpolation of " << lbl << " failed.\n"
1053 << " Cannot copy value to new grid at E = " << e
1054 << ", B = " << b << ", angle: " << a << "\n";
1055 continue;
1056 }
1057 tabClone[k][j][i] = val;
1058 }
1059 }
1060 }
1061 // Copy the values to the original table.
1062 tab.swap(tabClone);
1063}
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)
Definition Medium.cc:1408
std::string m_className
Definition Medium.hh:529

Referenced by SetFieldGrid().

◆ Clone() [2/2]

void Garfield::Medium::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 )
protected

Definition at line 1065 of file Medium.cc.

1070 {
1071 // If the table does not exist, do nothing.
1072 if (tab.empty()) return;
1073
1074 // Get the dimensions of the new grid.
1075 const auto nE = efields.size();
1076 const auto nB = bfields.size();
1077 const auto nA = angles.size();
1078
1079 // Create a temporary table to store the values at the new grid points.
1080 std::vector<std::vector<std::vector<std::vector<double> > > > tabClone;
1081 Init(nE, nB, nA, n, tabClone, init);
1082
1083 // Fill the temporary table.
1084 for (size_t l = 0; l < n; ++l) {
1085 for (size_t i = 0; i < nE; ++i) {
1086 const double e = efields[i];
1087 for (size_t j = 0; j < nB; ++j) {
1088 const double b = bfields[j];
1089 for (size_t k = 0; k < nA; ++k) {
1090 const double a = angles[k];
1091 double val = 0.;
1092 if (!Interpolate(e, b, a, tab[l], val, intp, extr)) {
1093 std::cerr << m_className << "::Clone:\n"
1094 << " Interpolation of " << lbl << " failed.\n"
1095 << " Cannot copy value to new grid at index " << l
1096 << ", E = " << e << ", B = " << b << ", angle: " << a
1097 << "\n";
1098 continue;
1099 }
1100 tabClone[l][k][j][i] = val;
1101 }
1102 }
1103 }
1104 }
1105 // Copy the values to the original table.
1106 tab.swap(tabClone);
1107}

◆ Diffusion() [1/2]

bool Garfield::Medium::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
protected

Definition at line 334 of file Medium.cc.

338 {
339
340 dl = dt = 0.;
341 // Compute the magnitude of the electric field.
342 const double e = sqrt(ex * ex + ey * ey + ez * ez);
343 const double e0 = ScaleElectricField(e);
344 if (e < Small || e0 < Small) return true;
345
346 // Compute the magnitude of the magnetic field.
347 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
348 // Compute the angle between B field and E field.
349 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
350
351 // Interpolate.
352 if (!difL.empty()) {
353 if (!Interpolate(e0, b, ebang, difL, dl, m_intpDif, m_extrDif)) dl = 0.;
354 }
355 if (!difT.empty()) {
356 if (!Interpolate(e0, b, ebang, difT, dt, m_intpDif, m_extrDif)) dt = 0.;
357 }
358
359 // If no data available, calculate
360 // the diffusion coefficients using the Einstein relation
361 if (difL.empty() || difT.empty()) {
362 const double d = sqrt(2. * BoltzmannConstant * m_temperature / e);
363 if (difL.empty()) dl = d;
364 if (difT.empty()) dt = d;
365 }
366 // Verify values and apply scaling.
367 dl = ScaleDiffusion(std::max(dl, 0.));
368 dt = ScaleDiffusion(std::max(dt, 0.));
369 return true;
370}
virtual double ScaleDiffusion(const double d) const
Definition Medium.hh:502
unsigned int m_intpDif
Definition Medium.hh:625
std::pair< unsigned int, unsigned int > m_extrDif
Definition Medium.hh:616
double m_temperature
Definition Medium.hh:540

Referenced by ElectronDiffusion(), ElectronDiffusion(), HoleDiffusion(), HoleDiffusion(), and IonDiffusion().

◆ Diffusion() [2/2]

bool Garfield::Medium::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
protected

Definition at line 372 of file Medium.cc.

375 {
376
377 // Initialise the tensor.
378 cov[0][0] = cov[0][1] = cov[0][2] = 0.;
379 cov[1][0] = cov[1][1] = cov[1][2] = 0.;
380 cov[2][0] = cov[2][1] = cov[2][2] = 0.;
381
382 if (diff.empty()) return false;
383
384 // Compute the magnitude of the electric field.
385 const double e = sqrt(ex * ex + ey * ey + ez * ez);
386 const double e0 = ScaleElectricField(e);
387 if (e < Small || e0 < Small) return true;
388
389 // Compute the magnitude of the magnetic field.
390 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
391 // Compute the angle between B field and E field.
392 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
393
394 for (int j = 0; j < 6; ++j) {
395 // Interpolate.
396 double y = 0.;
397 if (!Interpolate(e0, b, ebang, diff[j], y, m_intpDif, m_extrDif)) y = 0.;
398 // Apply scaling.
400 if (j < 3) {
401 cov[j][j] = y;
402 } else if (j == 3) {
403 cov[0][1] = cov[1][0] = y;
404 } else if (j == 4) {
405 cov[0][2] = cov[2][0] = y;
406 } else if (j == 5) {
407 cov[1][2] = cov[2][1] = y;
408 }
409 }
410 return true;
411}
virtual double ScaleDiffusionTensor(const double d) const
Definition Medium.hh:503

◆ DisableDebugging()

void Garfield::Medium::DisableDebugging ( )
inline

Definition at line 526 of file Medium.hh.

526{ m_debug = false; }

◆ ElectronAttachment()

bool Garfield::Medium::ElectronAttachment ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & eta )
virtual

Attachment coefficient [cm-1].

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 481 of file Medium.cc.

483 {
484
485 if (!Alpha(ex, ey, ez, bx, by, bz, m_eAtt, m_intpAtt, m_eThrAtt, m_extrAtt,
486 eta)) {
487 return false;
488 }
489 // Apply scaling.
490 eta = ScaleAttachment(eta);
491 return true;
492}
virtual double ScaleAttachment(const double eta) const
Definition Medium.hh:505
std::pair< unsigned int, unsigned int > m_extrAtt
Definition Medium.hh:618
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition Medium.hh:583
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
Definition Medium.cc:413
unsigned int m_intpAtt
Definition Medium.hh:627
unsigned int m_eThrAtt
Definition Medium.hh:609

Referenced by Garfield::MediumCdTe::ElectronAttachment(), Garfield::MediumDiamond::ElectronAttachment(), Garfield::MediumGaAs::ElectronAttachment(), Garfield::MediumGaN::ElectronAttachment(), and Garfield::MediumSilicon::ElectronAttachment().

◆ ElectronCollision()

bool Garfield::Medium::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 )
virtual

Sample the collision type. Update energy and direction vector.

Reimplemented in Garfield::MediumMagboltz, and Garfield::MediumSilicon.

Definition at line 556 of file Medium.cc.

559 {
560 type = level = -1;
561 e1 = e;
562 ndxc = band = 0;
563 RndmDirection(dx, dy, dz);
564
565 if (m_debug) PrintNotImplemented(m_className, "GetElectronCollision");
566 return false;
567}
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition Random.hh:128

◆ ElectronDiffusion() [1/2]

bool Garfield::Medium::ElectronDiffusion ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & dl,
double & dt )
virtual

Longitudinal and transverse diffusion coefficients [cm1/2].

Definition at line 452 of file Medium.cc.

455 {
456
457 return Diffusion(ex, ey, ez, bx, by, bz, m_eDifL, m_eDifT, dl, dt);
458}
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
Definition Medium.cc:334
std::vector< std::vector< std::vector< double > > > m_eDifL
Definition Medium.hh:580
std::vector< std::vector< std::vector< double > > > m_eDifT
Definition Medium.hh:581

◆ ElectronDiffusion() [2/2]

bool Garfield::Medium::ElectronDiffusion ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double cov[3][3] )
virtual

Diffusion tensor: diagonal elements are the diffusion coefficients [cm] along e, btrans, e x b, off-diagonal elements are the covariances

Definition at line 460 of file Medium.cc.

463 {
464
465 return Diffusion(ex, ey, ez, bx, by, bz, m_eDifM, cov);
466}
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
Definition Medium.hh:586

◆ ElectronLorentzAngle()

bool Garfield::Medium::ElectronLorentzAngle ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & lor )
virtual

Lorentz angle.

Definition at line 494 of file Medium.cc.

497 {
498 lor = 0.;
499 if (m_eLor.empty()) return false;
500
501 // Compute the magnitude of the electric field.
502 const double e = sqrt(ex * ex + ey * ey + ez * ez);
503 const double e0 = ScaleElectricField(e);
504 if (e < Small || e0 < Small) return true;
505
506 // Compute the magnitude of the magnetic field.
507 const double b = m_tab2d ? sqrt(bx * bx + by * by + bz * bz) : 0.;
508 // Compute the angle between B field and E field.
509 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
510
511 // Interpolate.
512 if (!Interpolate(e0, b, ebang, m_eLor, lor, m_intpLor, m_extrLor)) lor = 0.;
513 // Apply scaling.
514 lor = ScaleLorentzAngle(lor);
515 return true;
516}
virtual double ScaleLorentzAngle(const double lor) const
Definition Medium.hh:506
std::vector< std::vector< std::vector< double > > > m_eLor
Definition Medium.hh:584
std::pair< unsigned int, unsigned int > m_extrLor
Definition Medium.hh:619
unsigned int m_intpLor
Definition Medium.hh:628

◆ ElectronMobility()

double Garfield::Medium::ElectronMobility ( )
virtual

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

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 518 of file Medium.cc.

518 {
519 if (m_eVelE.empty()) return -1.;
520 return m_eVelE[0][0][0] / UnScaleElectricField(m_eFields[0]);
521}
virtual double UnScaleElectricField(const double e) const
Definition Medium.hh:500
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition Medium.hh:577

◆ ElectronTownsend()

bool Garfield::Medium::ElectronTownsend ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & alpha )
virtual

Ionisation coefficient [cm-1].

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 468 of file Medium.cc.

470 {
471
472 if (!Alpha(ex, ey, ez, bx, by, bz, m_eAlp, m_intpAlp, m_eThrAlp, m_extrAlp,
473 alpha)) {
474 return false;
475 }
476 // Apply scaling.
477 alpha = ScaleTownsend(alpha);
478 return true;
479}
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition Medium.hh:582
virtual double ScaleTownsend(const double alpha) const
Definition Medium.hh:504
std::pair< unsigned int, unsigned int > m_extrAlp
Definition Medium.hh:617
unsigned int m_eThrAlp
Definition Medium.hh:608
unsigned int m_intpAlp
Definition Medium.hh:626

Referenced by Garfield::MediumCdTe::ElectronTownsend(), Garfield::MediumDiamond::ElectronTownsend(), Garfield::MediumGaAs::ElectronTownsend(), Garfield::MediumGaN::ElectronTownsend(), and Garfield::MediumSilicon::ElectronTownsend().

◆ ElectronVelocity()

bool Garfield::Medium::ElectronVelocity ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & vx,
double & vy,
double & vz )
virtual

Drift velocity [cm / ns].

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 444 of file Medium.cc.

446 {
447
448 return Velocity(ex, ey, ez, bx, by, bz, m_eVelE, m_eVelB, m_eVelX, -1.,
449 vx, vy, vz);
450}
std::vector< std::vector< std::vector< double > > > m_eVelX
Definition Medium.hh:578
std::vector< std::vector< std::vector< double > > > m_eVelB
Definition Medium.hh:579
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
Definition Medium.cc:196

Referenced by Garfield::MediumCdTe::ElectronVelocity(), Garfield::MediumDiamond::ElectronVelocity(), Garfield::MediumGaAs::ElectronVelocity(), Garfield::MediumGaN::ElectronVelocity(), and Garfield::MediumSilicon::ElectronVelocity().

◆ EnableDebugging()

void Garfield::Medium::EnableDebugging ( )
inline

Switch on/off debugging messages.

Definition at line 525 of file Medium.hh.

525{ m_debug = true; }

◆ EnableDrift()

virtual void Garfield::Medium::EnableDrift ( const bool on = true)
inlinevirtual

Switch electron/ion/hole transport on/off.

Reimplemented in Garfield::MediumConductor, and Garfield::MediumPlastic.

Definition at line 70 of file Medium.hh.

70{ m_driftable = on; }

◆ EnablePrimaryIonisation()

virtual void Garfield::Medium::EnablePrimaryIonisation ( const bool on = true)
inlinevirtual

Make the medium ionisable or non-ionisable.

Reimplemented in Garfield::MediumConductor, and Garfield::MediumPlastic.

Definition at line 72 of file Medium.hh.

72 {
73 m_ionisable = on;
74 }

◆ GetAngle()

double Garfield::Medium::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
protected

Definition at line 1300 of file Medium.cc.

1302 {
1303 const double eb = emag * bmag;
1304 if (eb <= 0.) return m_bAngles[0];
1305 const double einb = fabs(ex * bx + ey * by + ez * bz);
1306 if (einb > 0.2 * eb) {
1307 double exb[3] = {ex * by - ey * bx, ex * bz - ez * bx, ez * by - ey * bz};
1308 return asin(
1309 std::min(1., sqrt(exb[0] * exb[0] + exb[1] * exb[1] + exb[2] * exb[2]) / eb));
1310 }
1311 return acos(std::min(1., einb / eb));
1312}
DoubleAc acos(const DoubleAc &f)
Definition DoubleAc.cpp:490
DoubleAc fabs(const DoubleAc &f)
Definition DoubleAc.h:615
DoubleAc asin(const DoubleAc &f)
Definition DoubleAc.cpp:470

Referenced by Alpha(), Diffusion(), Diffusion(), ElectronLorentzAngle(), IonVelocity(), NegativeIonVelocity(), and Velocity().

◆ GetAtomicNumber()

virtual double Garfield::Medium::GetAtomicNumber ( ) const
inlinevirtual

Get the effective atomic number.

Reimplemented in Garfield::MediumGas.

Definition at line 55 of file Medium.hh.

55{ return m_z; }

Referenced by Garfield::TrackSrim::NewTrack().

◆ GetAtomicWeight()

virtual double Garfield::Medium::GetAtomicWeight ( ) const
inlinevirtual

Get the effective atomic weight.

Reimplemented in Garfield::MediumGas.

Definition at line 59 of file Medium.hh.

59{ return m_a; }

◆ GetComponent()

void Garfield::Medium::GetComponent ( const unsigned int i,
std::string & label,
double & f )
virtual

Get the name and fraction of a given component.

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumGas.

Definition at line 106 of file Medium.cc.

106 {
107 if (i >= m_nComponents) {
108 std::cerr << m_className << "::GetComponent: Index out of range.\n";
109 }
110
111 label = m_name;
112 f = 1.;
113}
std::string m_name
Definition Medium.hh:538
unsigned int m_nComponents
Definition Medium.hh:536

Referenced by Garfield::TrackDegrade::Initialise(), and Garfield::TrackDegrade::SetupPenning().

◆ GetDeexcitationProduct()

bool Garfield::Medium::GetDeexcitationProduct ( const unsigned int i,
double & t,
double & s,
int & type,
double & energy ) const
virtual

Reimplemented in Garfield::MediumMagboltz.

Definition at line 569 of file Medium.cc.

571 {
572 if (m_debug) PrintNotImplemented(m_className, "GetDeexcitationProduct");
573 t = s = energy = 0.;
574 type = 0;
575 return false;
576}

◆ GetDielectricConstant()

double Garfield::Medium::GetDielectricConstant ( ) const
inline

Get the relative static dielectric constant.

Definition at line 45 of file Medium.hh.

45{ return m_epsilon; }
double m_epsilon
Definition Medium.hh:544

Referenced by Garfield::ComponentNeBem3d::GetVolume(), and Garfield::ComponentNeBem3d::Initialise().

◆ GetDielectricFunction()

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

Get the complex dielectric function at a given energy.

Reimplemented in Garfield::MediumSilicon.

Definition at line 734 of file Medium.cc.

735 {
736 if (i >= m_nComponents) {
737 std::cerr << m_className << "::GetDielectricFunction: Index out of range.\n";
738 return false;
739 }
740
741 if (e < 0.) {
742 std::cerr << m_className << "::GetDielectricFunction: Energy must be > 0.\n";
743 return false;
744 }
745
746 if (m_debug) PrintNotImplemented(m_className, "GetDielectricFunction");
747 eps1 = 1.;
748 eps2 = 0.;
749 return false;
750}

◆ GetElectronAttachment()

bool Garfield::Medium::GetElectronAttachment ( const size_t ie,
const size_t ib,
const size_t ia,
double & eta )
inline

Get an entry in the table of attachment coefficients.

Definition at line 291 of file Medium.hh.

292 {
293 return GetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
294 }
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
Definition Medium.cc:982

◆ GetElectronCollisionRate()

double Garfield::Medium::GetElectronCollisionRate ( const double e,
const int band = 0 )
virtual

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

Reimplemented in Garfield::MediumMagboltz, and Garfield::MediumSilicon.

Definition at line 550 of file Medium.cc.

551 {
552 if (m_debug) PrintNotImplemented(m_className, "GetElectronCollisionRate");
553 return 0.;
554}

◆ GetElectronEnergy()

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

Dispersion relation (energy vs. wave vector)

Reimplemented in Garfield::MediumSilicon.

Definition at line 523 of file Medium.cc.

525 {
526 if (band > 0) {
527 std::cerr << m_className << "::GetElectronEnergy:\n";
528 std::cerr << " Unknown band index.\n";
529 }
530
531 vx = SpeedOfLight * px / ElectronMass;
532 vy = SpeedOfLight * py / ElectronMass;
533 vz = SpeedOfLight * pz / ElectronMass;
534
535 return 0.5 * (px * px + py * py + pz * pz) / ElectronMass;
536}

◆ GetElectronLongitudinalDiffusion()

bool Garfield::Medium::GetElectronLongitudinalDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
double & dl )
inline

Get an entry in the table of longitudinal diffusion coefficients.

Definition at line 261 of file Medium.hh.

262 {
263 return GetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
264 }

◆ GetElectronLorentzAngle()

bool Garfield::Medium::GetElectronLorentzAngle ( const size_t ie,
const size_t ib,
const size_t ia,
double & lor )
inline

Get an entry in the table of Lorentz angles.

Definition at line 302 of file Medium.hh.

303 {
304 return GetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
305 }

◆ GetElectronMomentum()

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

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

Reimplemented in Garfield::MediumSilicon.

Definition at line 538 of file Medium.cc.

539 {
540 const double p = sqrt(2. * ElectronMass * e) / SpeedOfLight;
541 RndmDirection(px, py, pz, p);
542 band = -1;
543}

◆ GetElectronNullCollisionRate()

double Garfield::Medium::GetElectronNullCollisionRate ( const int band = 0)
virtual

Null-collision rate [ns-1].

Reimplemented in Garfield::MediumMagboltz, and Garfield::MediumSilicon.

Definition at line 545 of file Medium.cc.

545 {
546 if (m_debug) PrintNotImplemented(m_className, "GetElectronNullCollisionRate");
547 return 0.;
548}

◆ GetElectronTownsend()

bool Garfield::Medium::GetElectronTownsend ( const size_t ie,
const size_t ib,
const size_t ia,
double & alpha )
inline

Get an entry in the table of Townsend coefficients.

Definition at line 281 of file Medium.hh.

282 {
283 return GetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
284 }

◆ GetElectronTransverseDiffusion()

bool Garfield::Medium::GetElectronTransverseDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
double & dt )
inline

Get an entry in the table of transverse diffusion coefficients.

Definition at line 271 of file Medium.hh.

272 {
273 return GetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
274 }

◆ GetElectronVelocityB()

bool Garfield::Medium::GetElectronVelocityB ( const size_t ie,
const size_t ib,
const size_t ia,
double & v )
inline

Get an entry in the table of drift speeds along Btrans.

Definition at line 251 of file Medium.hh.

252 {
253 return GetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
254 }

◆ GetElectronVelocityE()

bool Garfield::Medium::GetElectronVelocityE ( const size_t ie,
const size_t ib,
const size_t ia,
double & v )
inline

Get an entry in the table of drift speeds along E.

Definition at line 231 of file Medium.hh.

232 {
233 return GetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
234 }

◆ GetElectronVelocityExB()

bool Garfield::Medium::GetElectronVelocityExB ( const size_t ie,
const size_t ib,
const size_t ia,
double & v )
inline

Get an entry in the table of drift speeds along ExB.

Definition at line 241 of file Medium.hh.

242 {
243 return GetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
244 }

◆ GetEntry()

bool Garfield::Medium::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
protected

◆ GetExtrapolationIndex()

bool Garfield::Medium::GetExtrapolationIndex ( std::string str,
unsigned int & nb ) const
protected

Definition at line 1235 of file Medium.cc.

1235 {
1236 // Convert to upper-case.
1237 std::transform(str.begin(), str.end(), str.begin(), toupper);
1238
1239 if (str == "CONST" || str == "CONSTANT") {
1240 nb = 0;
1241 } else if (str == "LIN" || str == "LINEAR") {
1242 nb = 1;
1243 } else if (str == "EXP" || str == "EXPONENTIAL") {
1244 nb = 2;
1245 } else {
1246 return false;
1247 }
1248
1249 return true;
1250}

Referenced by SetExtrapolationMethod().

◆ GetFanoFactor()

double Garfield::Medium::GetFanoFactor ( ) const
inline

Get the Fano factor.

Definition at line 90 of file Medium.hh.

90{ return m_fano; }

Referenced by Garfield::TrackSrim::NewTrack(), and Garfield::TrackTrim::NewTrack().

◆ GetFieldGrid()

void Garfield::Medium::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.

Definition at line 958 of file Medium.cc.

960 {
961 efields = m_eFields;
962 bfields = m_bFields;
963 angles = m_bAngles;
964}

◆ GetHoleAttachment()

bool Garfield::Medium::GetHoleAttachment ( const size_t ie,
const size_t ib,
const size_t ia,
double & eta )
inline

Get an entry in the table of attachment coefficients.

Definition at line 374 of file Medium.hh.

375 {
376 return GetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
377 }
std::vector< std::vector< std::vector< double > > > m_hAtt
Definition Medium.hh:595

◆ GetHoleLongitudinalDiffusion()

bool Garfield::Medium::GetHoleLongitudinalDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
double & dl )
inline

Get an entry in the table of longitudinal diffusion coefficients.

Definition at line 344 of file Medium.hh.

345 {
346 return GetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
347 }
std::vector< std::vector< std::vector< double > > > m_hDifL
Definition Medium.hh:592

◆ GetHoleTownsend()

bool Garfield::Medium::GetHoleTownsend ( const size_t ie,
const size_t ib,
const size_t ia,
double & alpha )
inline

Get an entry in the table of Townsend coefficients.

Definition at line 364 of file Medium.hh.

365 {
366 return GetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
367 }
std::vector< std::vector< std::vector< double > > > m_hAlp
Definition Medium.hh:594

◆ GetHoleTransverseDiffusion()

bool Garfield::Medium::GetHoleTransverseDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
double & dt )
inline

Get an entry in the table of transverse diffusion coefficients.

Definition at line 354 of file Medium.hh.

355 {
356 return GetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
357 }
std::vector< std::vector< std::vector< double > > > m_hDifT
Definition Medium.hh:593

◆ GetHoleVelocityB()

bool Garfield::Medium::GetHoleVelocityB ( const size_t ie,
const size_t ib,
const size_t ia,
double & v )
inline

Get an entry in the table of drift speeds along Btrans.

Definition at line 333 of file Medium.hh.

334 {
335 return GetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
336 }
std::vector< std::vector< std::vector< double > > > m_hVelB
Definition Medium.hh:591

◆ GetHoleVelocityE()

bool Garfield::Medium::GetHoleVelocityE ( const size_t ie,
const size_t ib,
const size_t ia,
double & v )
inline

Get an entry in the table of drift speeds along E.

Definition at line 313 of file Medium.hh.

314 {
315 return GetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
316 }
std::vector< std::vector< std::vector< double > > > m_hVelE
Definition Medium.hh:589

◆ GetHoleVelocityExB()

bool Garfield::Medium::GetHoleVelocityExB ( const size_t ie,
const size_t ib,
const size_t ia,
double & v )
inline

Get an entry in the table of drift speeds along ExB.

Definition at line 323 of file Medium.hh.

324 {
325 return GetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
326 }
std::vector< std::vector< std::vector< double > > > m_hVelX
Definition Medium.hh:590

◆ GetId()

int Garfield::Medium::GetId ( ) const
inline

Return the id number of the class instance.

Definition at line 24 of file Medium.hh.

24{ return m_id; }

Referenced by Garfield::ComponentNeBem3d::GetVolume(), and Garfield::ViewGeometry::Plot3d().

◆ GetIonDissociation()

bool Garfield::Medium::GetIonDissociation ( const size_t ie,
const size_t ib,
const size_t ia,
double & diss )
inline

Get an entry in the table of dissociation coefficients.

Definition at line 421 of file Medium.hh.

422 {
423 return GetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
424 }
std::vector< std::vector< std::vector< double > > > m_iDis
Definition Medium.hh:603

◆ GetIonLongitudinalDiffusion()

bool Garfield::Medium::GetIonLongitudinalDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
double & dl )
inline

Get an entry in the table of longitudinal diffusion coefficients.

Definition at line 401 of file Medium.hh.

402 {
403 return GetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
404 }
std::vector< std::vector< std::vector< double > > > m_iDifL
Definition Medium.hh:601

◆ GetIonMobility()

bool Garfield::Medium::GetIonMobility ( const size_t ie,
const size_t ib,
const size_t ia,
double & mu )
inline

Get an entry in the table of ion mobilities.

Definition at line 390 of file Medium.hh.

391 {
392 return GetEntry(ie, ib, ia, "IonMobility", m_iMob, mu);
393 }
std::vector< std::vector< std::vector< double > > > m_iMob
Definition Medium.hh:600

◆ GetIonTransverseDiffusion()

bool Garfield::Medium::GetIonTransverseDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
double & dt )
inline

Get an entry in the table of transverse diffusion coefficients.

Definition at line 411 of file Medium.hh.

412 {
413 return GetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
414 }
std::vector< std::vector< std::vector< double > > > m_iDifT
Definition Medium.hh:602

◆ GetMassDensity()

double Garfield::Medium::GetMassDensity ( ) const
virtual

Get the mass density [g/cm3].

Reimplemented in Garfield::MediumGas.

Definition at line 102 of file Medium.cc.

102 {
103 return m_density * AtomicMassUnit * m_a;
104}
double m_density
Definition Medium.hh:550

Referenced by Garfield::TrackHeed::NewTrack(), Garfield::TrackSrim::NewTrack(), Garfield::GeometryRoot::SetMedium(), Garfield::TrackHeed::TransportDeltaElectron(), and Garfield::TrackHeed::TransportPhoton().

◆ GetName()

◆ GetNegativeIonMobility()

bool Garfield::Medium::GetNegativeIonMobility ( const size_t ie,
const size_t ib,
const size_t ia,
double & mu )
inline

Get an entry in the table of negative ion mobilities.

Definition at line 432 of file Medium.hh.

433 {
434 return GetEntry(ie, ib, ia, "NegativeIonMobility", m_nMob, mu);
435 }
std::vector< std::vector< std::vector< double > > > m_nMob
Definition Medium.hh:605

◆ GetNumberDensity()

virtual double Garfield::Medium::GetNumberDensity ( ) const
inlinevirtual

Get the number density [cm-3].

Reimplemented in Garfield::MediumGas.

Definition at line 63 of file Medium.hh.

63{ return m_density; }

Referenced by Garfield::TrackBichsel::Initialise(), Garfield::TrackElectron::NewTrack(), Garfield::TrackPAI::NewTrack(), and Garfield::TrackSrim::NewTrack().

◆ GetNumberOfComponents()

unsigned int Garfield::Medium::GetNumberOfComponents ( ) const
inline

Get number of components of the medium.

Definition at line 48 of file Medium.hh.

48{ return m_nComponents; }

Referenced by Garfield::TrackDegrade::Initialise(), and Garfield::TrackDegrade::SetupPenning().

◆ GetNumberOfDeexcitationProducts()

virtual unsigned int Garfield::Medium::GetNumberOfDeexcitationProducts ( ) const
inlinevirtual

Reimplemented in Garfield::MediumMagboltz.

Definition at line 159 of file Medium.hh.

159{ return 0; }

◆ GetOpticalDataRange()

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

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

Reimplemented in Garfield::MediumSilicon.

Definition at line 722 of file Medium.cc.

723 {
724 if (i >= m_nComponents) {
725 std::cerr << m_className << "::GetOpticalDataRange: Index out of range.\n";
726 return false;
727 }
728
729 if (m_debug) PrintNotImplemented(m_className, "GetOpticalDataRange");
730 emin = emax = 0.;
731 return false;
732}

◆ GetPhotoAbsorptionCrossSection()

bool Garfield::Medium::GetPhotoAbsorptionCrossSection ( const double e,
double & sigma,
const unsigned int i = 0 )
virtual

Reimplemented in Garfield::MediumGas.

Definition at line 752 of file Medium.cc.

753 {
754 if (i >= m_nComponents) {
755 std::cerr << m_className << "::GetPhotoAbsorptionCrossSection:\n";
756 std::cerr << " Component " << i << " does not exist.\n";
757 return false;
758 }
759
760 if (e < 0.) {
761 std::cerr << m_className << "::GetPhotoAbsorptionCrossSection:\n";
762 std::cerr << " Energy must be > 0.\n";
763 return false;
764 }
765
766 if (m_debug) {
767 PrintNotImplemented(m_className, "GetPhotoAbsorptionCrossSection");
768 }
769 sigma = 0.;
770 return false;
771}

Referenced by GetPhotonCollisionRate().

◆ GetPhotonCollision()

bool Garfield::Medium::GetPhotonCollision ( const double e,
int & type,
int & level,
double & e1,
double & ctheta,
int & nsec,
double & esec )
virtual

Reimplemented in Garfield::MediumMagboltz.

Definition at line 780 of file Medium.cc.

782 {
783 type = level = -1;
784 e1 = e;
785 ctheta = 1.;
786 nsec = 0;
787 esec = 0.;
788 return false;
789}

◆ GetPhotonCollisionRate()

double Garfield::Medium::GetPhotonCollisionRate ( const double e)
virtual

Reimplemented in Garfield::MediumMagboltz.

Definition at line 773 of file Medium.cc.

773 {
774 double sigma = 0.;
775 if (!GetPhotoAbsorptionCrossSection(e, sigma)) return 0.;
776
777 return sigma * m_density * SpeedOfLight;
778}
virtual bool GetPhotoAbsorptionCrossSection(const double e, double &sigma, const unsigned int i=0)
Definition Medium.cc:752

◆ GetPressure()

double Garfield::Medium::GetPressure ( ) const
inline

◆ GetTemperature()

double Garfield::Medium::GetTemperature ( ) const
inline

Get the temperature [K].

Definition at line 37 of file Medium.hh.

37{ return m_temperature; }

Referenced by Garfield::TrackDegrade::Initialise(), Garfield::TrackDegrade::IsInside(), and Garfield::TrackDegrade::NewTrack().

◆ GetW()

double Garfield::Medium::GetW ( ) const
inline

Get the W value.

Definition at line 86 of file Medium.hh.

86{ return m_w; }

Referenced by Garfield::TrackSrim::NewTrack(), and Garfield::TrackTrim::NewTrack().

◆ HoleAttachment()

bool Garfield::Medium::HoleAttachment ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & eta )
virtual

Attachment coefficient [cm-1].

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 612 of file Medium.cc.

614 {
615
616 if (!Alpha(ex, ey, ez, bx, by, bz, m_hAtt, m_intpAtt, m_hThrAtt, m_extrAtt,
617 eta)) {
618 return false;
619 }
620 // Apply scaling.
621 eta = ScaleAttachment(eta);
622 return true;
623}
unsigned int m_hThrAtt
Definition Medium.hh:611

Referenced by Garfield::MediumCdTe::HoleAttachment(), Garfield::MediumDiamond::HoleAttachment(), Garfield::MediumGaAs::HoleAttachment(), Garfield::MediumGaN::HoleAttachment(), and Garfield::MediumSilicon::HoleAttachment().

◆ HoleDiffusion() [1/2]

bool Garfield::Medium::HoleDiffusion ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & dl,
double & dt )
virtual

Longitudinal and transverse diffusion coefficients [cm1/2].

Definition at line 586 of file Medium.cc.

588 {
589 return Diffusion(ex, ey, ez, bx, by, bz, m_hDifL, m_hDifT, dl, dt);
590}

◆ HoleDiffusion() [2/2]

bool Garfield::Medium::HoleDiffusion ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double cov[3][3] )
virtual

Diffusion tensor.

Definition at line 592 of file Medium.cc.

594 {
595
596 return Diffusion(ex, ey, ez, bx, by, bz, m_hDifM, cov);
597}
std::vector< std::vector< std::vector< std::vector< double > > > > m_hDifM
Definition Medium.hh:597

◆ HoleMobility()

double Garfield::Medium::HoleMobility ( )
virtual

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

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 625 of file Medium.cc.

625 {
626 if (m_hVelE.empty()) return -1.;
627 return m_hVelE[0][0][0] / UnScaleElectricField(m_eFields[0]);
628}

◆ HoleTownsend()

bool Garfield::Medium::HoleTownsend ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & alpha )
virtual

Ionisation coefficient [cm-1].

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 599 of file Medium.cc.

601 {
602
603 if (!Alpha(ex, ey, ez, bx, by, bz, m_hAlp, m_intpAlp, m_hThrAlp, m_extrAlp,
604 alpha)) {
605 return false;
606 }
607 // Apply scaling.
608 alpha = ScaleTownsend(alpha);
609 return true;
610}
unsigned int m_hThrAlp
Definition Medium.hh:610

Referenced by Garfield::MediumCdTe::HoleTownsend(), Garfield::MediumDiamond::HoleTownsend(), Garfield::MediumGaAs::HoleTownsend(), Garfield::MediumGaN::HoleTownsend(), and Garfield::MediumSilicon::HoleTownsend().

◆ HoleVelocity()

bool Garfield::Medium::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 )
virtual

◆ Init() [1/2]

void Garfield::Medium::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

Definition at line 1419 of file Medium.cc.

1422 {
1423 if (nE == 0 || nB == 0 || nA == 0 || nT == 0) {
1424 std::cerr << m_className << "::Init: Invalid grid.\n";
1425 return;
1426 }
1427
1428 tab.assign(nT, std::vector<std::vector<std::vector<double> > >(
1429 nA, std::vector<std::vector<double> >(
1430 nB, std::vector<double>(nE, val))));
1431}

◆ Init() [2/2]

void Garfield::Medium::Init ( const size_t nE,
const size_t nB,
const size_t nA,
std::vector< std::vector< std::vector< double > > > & tab,
const double val )
protected

Definition at line 1408 of file Medium.cc.

1410 {
1411 if (nE == 0 || nB == 0 || nA == 0) {
1412 std::cerr << m_className << "::Init: Invalid grid.\n";
1413 return;
1414 }
1415 tab.assign(
1416 nA, std::vector<std::vector<double> >(nB, std::vector<double>(nE, val)));
1417}

Referenced by Clone(), Clone(), Garfield::MediumMagboltz::GenerateGasTable(), Garfield::MediumGas::LoadGasFile(), Garfield::MediumGas::MergeGasFile(), SetEntry(), and SetIonMobility().

◆ Interpolate()

bool Garfield::Medium::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
protected

Definition at line 1314 of file Medium.cc.

1318 {
1319 if (table.empty()) {
1320 y = 0.;
1321 return false; // TODO: true!
1322 }
1323
1324 if (m_tab2d) {
1326 m_bAngles.size(), m_bFields.size(), m_eFields.size(),
1327 a, b, e, y, intp);
1328 } else {
1329 y = Interpolate1D(e, table[0][0], m_eFields, intp, extr);
1330 }
1331 return true;
1332}
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
Definition Medium.cc:1334
bool Boxin3(const std::vector< std::vector< std::vector< double > > > &value, const std::vector< double > &xAxis, const std::vector< double > &yAxis, const std::vector< double > &zAxis, const int nx, const int ny, const int nz, const double xx, const double yy, const double zz, double &f, const int iOrder)
Definition Numerics.cc:1545

Referenced by Alpha(), Clone(), Clone(), Diffusion(), Diffusion(), ElectronLorentzAngle(), IonVelocity(), NegativeIonVelocity(), and Velocity().

◆ Interpolate1D()

double Garfield::Medium::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
protected

Definition at line 1334 of file Medium.cc.

1337 {
1338 // This function is a generalized version of the Fortran functions
1339 // GASVEL, GASVT1, GASVT2, GASLOR, GASMOB, GASDFT, and GASDFL
1340 // for the case of a 1D table. All variables are generic.
1341
1342 const auto nSizeTable = fields.size();
1343
1344 if (e < 0. || nSizeTable < 1) return 0.;
1345
1346 double result = 0.;
1347
1348 if (nSizeTable == 1) {
1349 // Only one point
1350 result = table[0];
1351 } else if (e < fields[0]) {
1352 // Extrapolation towards small fields
1353 if (fields[0] >= fields[1]) {
1354 if (m_debug) {
1355 std::cerr << m_className << "::Interpolate1D:\n";
1356 std::cerr << " First two field values coincide.\n";
1357 std::cerr << " No extrapolation to lower fields.\n";
1358 }
1359 result = table[0];
1360 } else if (extr.first == 1) {
1361 // Linear extrapolation
1362 const double extr4 = (table[1] - table[0]) / (fields[1] - fields[0]);
1363 const double extr3 = table[0] - extr4 * fields[0];
1364 result = extr3 + extr4 * e;
1365 } else if (extr.first == 2) {
1366 // Logarithmic extrapolation
1367 const double extr4 = log(table[1] / table[0]) / (fields[1] - fields[0]);
1368 const double extr3 = log(table[0] - extr4 * fields[0]);
1369 result = std::exp(std::min(50., extr3 + extr4 * e));
1370 } else {
1371 result = table[0];
1372 }
1373 } else if (e > fields[nSizeTable - 1]) {
1374 // Extrapolation towards large fields
1375 if (fields[nSizeTable - 1] <= fields[nSizeTable - 2]) {
1376 if (m_debug) {
1377 std::cerr << m_className << "::Interpolate1D:\n";
1378 std::cerr << " Last two field values coincide.\n";
1379 std::cerr << " No extrapolation to higher fields.\n";
1380 }
1381 result = table[nSizeTable - 1];
1382 } else if (extr.second == 1) {
1383 // Linear extrapolation
1384 const double extr2 = (table[nSizeTable - 1] - table[nSizeTable - 2]) /
1385 (fields[nSizeTable - 1] - fields[nSizeTable - 2]);
1386 const double extr1 =
1387 table[nSizeTable - 1] - extr2 * fields[nSizeTable - 1];
1388 result = extr1 + extr2 * e;
1389 } else if (extr.second == 2) {
1390 // Logarithmic extrapolation
1391 const double extr2 = log(table[nSizeTable - 1] / table[nSizeTable - 2]) /
1392 (fields[nSizeTable - 1] - fields[nSizeTable - 2]);
1393 const double extr1 =
1394 log(table[nSizeTable - 1]) - extr2 * fields[nSizeTable - 1];
1395 result = exp(std::min(50., extr1 + extr2 * e));
1396 } else {
1397 result = table[nSizeTable - 1];
1398 }
1399 } else {
1400 // Intermediate points, spline interpolation (not implemented).
1401 // Intermediate points, Newtonian interpolation
1402 result = Numerics::Divdif(table, fields, nSizeTable, e, intpMeth);
1403 }
1404
1405 return result;
1406}
double Divdif(const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
Definition Numerics.cc:1255

Referenced by Interpolate(), and SetIonMobility().

◆ IonDiffusion()

bool Garfield::Medium::IonDiffusion ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & dl,
double & dt )
virtual

Longitudinal and transverse diffusion coefficients [cm1/2].

Definition at line 660 of file Medium.cc.

662 {
663
664 return Diffusion(ex, ey, ez, bx, by, bz, m_iDifL, m_iDifT, dl, dt);
665}

◆ IonDissociation()

bool Garfield::Medium::IonDissociation ( const double ex,
const double ey,
const double ez,
const double bx,
const double by,
const double bz,
double & diss )
virtual

Dissociation coefficient.

Definition at line 667 of file Medium.cc.

669 {
670
671 if (!Alpha(ex, ey, ez, bx, by, bz, m_iDis, m_intpDis, m_iThrDis, m_extrDis,
672 diss)) {
673 return false;
674 }
675 // Apply scaling.
676 diss = ScaleDissociation(diss);
677 return true;
678}
unsigned int m_intpDis
Definition Medium.hh:630
virtual double ScaleDissociation(const double diss) const
Definition Medium.hh:507
std::pair< unsigned int, unsigned int > m_extrDis
Definition Medium.hh:621
unsigned int m_iThrDis
Definition Medium.hh:612

◆ IonMobility()

double Garfield::Medium::IonMobility ( )
virtual

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

Definition at line 680 of file Medium.cc.

680 {
681 return m_iMob.empty() ? -1. : m_iMob[0][0][0];
682}

Referenced by NegativeIonMobility().

◆ IonVelocity()

bool Garfield::Medium::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 )
virtual

Ion drift velocity [cm / ns].

Definition at line 630 of file Medium.cc.

632 {
633 vx = vy = vz = 0.;
634 if (m_iMob.empty()) return false;
635 // Compute the magnitude of the electric field.
636 const double e = sqrt(ex * ex + ey * ey + ez * ez);
637 const double e0 = ScaleElectricField(e);
638 if (e < Small || e0 < Small) return true;
639 // Compute the magnitude of the electric field.
640 const double b = sqrt(bx * bx + by * by + bz * bz);
641
642 // Compute the angle between B field and E field.
643 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
644 double mu = 0.;
645 if (!Interpolate(e0, b, ebang, m_iMob, mu, m_intpMob, m_extrMob)) mu = 0.;
646
647 constexpr double q = 1.;
648 mu *= q;
649 if (b < Small) {
650 vx = mu * ex;
651 vy = mu * ey;
652 vz = mu * ez;
653 } else {
654 Langevin(ex, ey, ez, bx, by, bz, mu, vx, vy, vz);
655 }
656
657 return true;
658}
unsigned int m_intpMob
Definition Medium.hh:629
std::pair< unsigned int, unsigned int > m_extrMob
Definition Medium.hh:620
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

◆ IsConductor()

virtual bool Garfield::Medium::IsConductor ( ) const
inlinevirtual

Is this medium a conductor?

Reimplemented in Garfield::MediumConductor.

Definition at line 32 of file Medium.hh.

32{ return false; }

Referenced by Garfield::ComponentNeBem3d::Initialise().

◆ IsDriftable()

◆ IsGas()

virtual bool Garfield::Medium::IsGas ( ) const
inlinevirtual

◆ IsIonisable()

◆ IsMicroscopic()

bool Garfield::Medium::IsMicroscopic ( ) const
inline

Does the medium have electron scattering rates?

Definition at line 79 of file Medium.hh.

79{ return m_microscopic; }

◆ IsSemiconductor()

virtual bool Garfield::Medium::IsSemiconductor ( ) const
inlinevirtual

Is this medium a semiconductor?

Reimplemented in Garfield::MediumCdTe, Garfield::MediumDiamond, Garfield::MediumGaAs, Garfield::MediumGaN, and Garfield::MediumSilicon.

Definition at line 30 of file Medium.hh.

30{ return false; }

Referenced by Garfield::ViewGeometry::Plot3d().

◆ Langevin() [1/2]

void Garfield::Medium::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 )
staticprotected

Definition at line 317 of file Medium.cc.

320 {
321
322 bx *= Tesla2Internal;
323 by *= Tesla2Internal;
324 bz *= Tesla2Internal;
325 const double b2 = bx * bx + by * by + bz * bz;
326 const double mu2 = muH * muH;
327 const double f = mu / (1. + mu2 * b2);
328 const double eb = bx * ex + by * ey + bz * ez;
329 vx = f * (ex + muH * (ey * bz - ez * by) + mu2 * bx * eb);
330 vy = f * (ey + muH * (ez * bx - ex * bz) + mu2 * by * eb);
331 vz = f * (ez + muH * (ex * by - ey * bx) + mu2 * bz * eb);
332}

◆ Langevin() [2/2]

void Garfield::Medium::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 )
staticprotected

Definition at line 301 of file Medium.cc.

303 {
304
305 bx *= Tesla2Internal;
306 by *= Tesla2Internal;
307 bz *= Tesla2Internal;
308 const double b2 = bx * bx + by * by + bz * bz;
309 const double mu2 = mu * mu;
310 const double eb = bx * ex + by * ey + bz * ez;
311 const double f = mu / (1. + mu2 * b2);
312 vx = f * (ex + mu * (ey * bz - ez * by) + mu2 * bx * eb);
313 vy = f * (ey + mu * (ez * bx - ex * bz) + mu2 * by * eb);
314 vz = f * (ez + mu * (ex * by - ey * bx) + mu2 * bz * eb);
315}

Referenced by Garfield::MediumCdTe::ElectronVelocity(), Garfield::MediumDiamond::ElectronVelocity(), Garfield::MediumGaAs::ElectronVelocity(), Garfield::MediumGaN::ElectronVelocity(), Garfield::MediumSilicon::ElectronVelocity(), Garfield::MediumCdTe::HoleVelocity(), Garfield::MediumDiamond::HoleVelocity(), Garfield::MediumGaAs::HoleVelocity(), Garfield::MediumGaN::HoleVelocity(), Garfield::MediumSilicon::HoleVelocity(), IonVelocity(), NegativeIonVelocity(), and Velocity().

◆ NegativeIonMobility()

double Garfield::Medium::NegativeIonMobility ( )
virtual

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

Definition at line 718 of file Medium.cc.

718 {
719 return m_nMob.empty() ? IonMobility() : m_nMob[0][0][0];
720}
virtual double IonMobility()
Low-field ion mobility [cm2 V-1 ns-1].
Definition Medium.cc:680

◆ NegativeIonVelocity()

bool Garfield::Medium::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 )
virtual

Negative ion drift velocity [cm / ns].

Definition at line 684 of file Medium.cc.

687 {
688 vx = vy = vz = 0.;
689 if (m_iMob.empty() && m_nMob.empty()) return false;
690
691 // Compute the magnitude of the electric field.
692 const double e = sqrt(ex * ex + ey * ey + ez * ez);
693 const double e0 = ScaleElectricField(e);
694 if (e < Small || e0 < Small) return true;
695 // Compute the magnitude of the electric field.
696 const double b = sqrt(bx * bx + by * by + bz * bz);
697
698 // Compute the angle between B field and E field.
699 const double ebang = m_tab2d ? GetAngle(ex, ey, ez, bx, by, bz, e, b) : 0.;
700 double mu = 0.;
701 if (m_nMob.empty()) {
702 if (!Interpolate(e0, b, ebang, m_iMob, mu, m_intpMob, m_extrMob)) mu = 0.;
703 } else {
704 if (!Interpolate(e0, b, ebang, m_nMob, mu, m_intpMob, m_extrMob)) mu = 0.;
705 }
706 constexpr double q = -1.;
707 mu *= q;
708 if (b < Small) {
709 vx = mu * ex;
710 vy = mu * ey;
711 vz = mu * ez;
712 } else {
713 Langevin(ex, ey, ez, bx, by, bz, mu, vx, vy, vz);
714 }
715 return true;
716}

◆ PlotAlphaEta()

void Garfield::Medium::PlotAlphaEta ( const std::string & carriers,
TPad * pad )

Plot Townsend and attachment coefficients.

Definition at line 189 of file Medium.cc.

189 {
190 ViewMedium view;
191 view.SetMedium(this);
192 if (pad) view.SetCanvas(pad);
193 view.PlotAlphaEta(opt, 'e');
194}

◆ PlotAttachment()

void Garfield::Medium::PlotAttachment ( const std::string & carriers,
TPad * pad )

Plot the attachment coefficient(s) as function of the electric field.

Definition at line 182 of file Medium.cc.

182 {
183 ViewMedium view;
184 view.SetMedium(this);
185 if (pad) view.SetCanvas(pad);
186 view.PlotAttachment(opt, 'e');
187}

◆ PlotDiffusion()

void Garfield::Medium::PlotDiffusion ( const std::string & carriers,
TPad * pad )

Plot the diffusion coefficients as function of the electric field.

Definition at line 168 of file Medium.cc.

168 {
169 ViewMedium view;
170 view.SetMedium(this);
171 if (pad) view.SetCanvas(pad);
172 view.PlotDiffusion(opt, 'e');
173}

◆ PlotTownsend()

void Garfield::Medium::PlotTownsend ( const std::string & carriers,
TPad * pad )

Plot the Townsend coefficient(s) as function of the electric field.

Definition at line 175 of file Medium.cc.

175 {
176 ViewMedium view;
177 view.SetMedium(this);
178 if (pad) view.SetCanvas(pad);
179 view.PlotTownsend(opt, 'e');
180}

◆ PlotVelocity()

void Garfield::Medium::PlotVelocity ( const std::string & carriers,
TPad * pad )

Plot the drift velocity as function of the electric field.

Definition at line 161 of file Medium.cc.

161 {
162 ViewMedium view;
163 view.SetMedium(this);
164 if (pad) view.SetCanvas(pad);
165 view.PlotVelocity(opt, 'e');
166}

◆ ResetElectronAttachment()

void Garfield::Medium::ResetElectronAttachment ( )
inline

Definition at line 451 of file Medium.hh.

451{ m_eAtt.clear(); }

Referenced by ResetTables().

◆ ResetElectronDiffusion()

void Garfield::Medium::ResetElectronDiffusion ( )
inline

Definition at line 445 of file Medium.hh.

445 {
446 m_eDifL.clear();
447 m_eDifT.clear();
448 m_eDifM.clear();
449 }

Referenced by ResetTables().

◆ ResetElectronLorentzAngle()

void Garfield::Medium::ResetElectronLorentzAngle ( )
inline

Definition at line 452 of file Medium.hh.

452{ m_eLor.clear(); }

Referenced by ResetTables().

◆ ResetElectronTownsend()

void Garfield::Medium::ResetElectronTownsend ( )
inline

Definition at line 450 of file Medium.hh.

450{ m_eAlp.clear(); }

Referenced by ResetTables().

◆ ResetElectronVelocity()

void Garfield::Medium::ResetElectronVelocity ( )
inline

Definition at line 440 of file Medium.hh.

440 {
441 m_eVelE.clear();
442 m_eVelB.clear();
443 m_eVelX.clear();
444 }

Referenced by ResetTables().

◆ ResetHoleAttachment()

void Garfield::Medium::ResetHoleAttachment ( )
inline

Definition at line 465 of file Medium.hh.

465{ m_hAtt.clear(); }

Referenced by ResetTables().

◆ ResetHoleDiffusion()

void Garfield::Medium::ResetHoleDiffusion ( )
inline

Definition at line 459 of file Medium.hh.

459 {
460 m_hDifL.clear();
461 m_hDifT.clear();
462 m_hDifM.clear();
463 }

Referenced by ResetTables().

◆ ResetHoleTownsend()

void Garfield::Medium::ResetHoleTownsend ( )
inline

Definition at line 464 of file Medium.hh.

464{ m_hAlp.clear(); }

Referenced by ResetTables().

◆ ResetHoleVelocity()

void Garfield::Medium::ResetHoleVelocity ( )
inline

Definition at line 454 of file Medium.hh.

454 {
455 m_hVelE.clear();
456 m_hVelB.clear();
457 m_hVelX.clear();
458 }

Referenced by ResetTables().

◆ ResetIonDiffusion()

void Garfield::Medium::ResetIonDiffusion ( )
inline

Definition at line 468 of file Medium.hh.

468 {
469 m_iDifL.clear();
470 m_iDifT.clear();
471 }

Referenced by ResetTables().

◆ ResetIonDissociation()

void Garfield::Medium::ResetIonDissociation ( )
inline

Definition at line 472 of file Medium.hh.

472{ m_iDis.clear(); }

Referenced by ResetTables().

◆ ResetIonMobility()

void Garfield::Medium::ResetIonMobility ( )
inline

Definition at line 467 of file Medium.hh.

467{ m_iMob.clear(); }

Referenced by ResetTables(), and SetIonMobility().

◆ ResetNegativeIonMobility()

void Garfield::Medium::ResetNegativeIonMobility ( )
inline

Definition at line 473 of file Medium.hh.

473{ m_nMob.clear(); }

Referenced by ResetTables(), and SetIonMobility().

◆ ResetTables()

void Garfield::Medium::ResetTables ( )
virtual

Reset all tables of transport parameters.

Reimplemented in Garfield::MediumGas.

Definition at line 999 of file Medium.cc.

999 {
1005
1010
1014
1016}
void ResetHoleAttachment()
Definition Medium.hh:465
void ResetIonMobility()
Definition Medium.hh:467
void ResetElectronAttachment()
Definition Medium.hh:451
void ResetIonDiffusion()
Definition Medium.hh:468
void ResetElectronLorentzAngle()
Definition Medium.hh:452
void ResetHoleVelocity()
Definition Medium.hh:454
void ResetIonDissociation()
Definition Medium.hh:472
void ResetElectronDiffusion()
Definition Medium.hh:445
void ResetNegativeIonMobility()
Definition Medium.hh:473
void ResetHoleDiffusion()
Definition Medium.hh:459
void ResetElectronTownsend()
Definition Medium.hh:450
void ResetHoleTownsend()
Definition Medium.hh:464
void ResetElectronVelocity()
Definition Medium.hh:440

Referenced by Garfield::MediumGas::ResetTables().

◆ ScaleAttachment()

virtual double Garfield::Medium::ScaleAttachment ( const double eta) const
inlinevirtual

Reimplemented in Garfield::MediumGas.

Definition at line 505 of file Medium.hh.

505{ return eta; }

Referenced by ElectronAttachment(), and HoleAttachment().

◆ ScaleDiffusion()

virtual double Garfield::Medium::ScaleDiffusion ( const double d) const
inlinevirtual

Reimplemented in Garfield::MediumGas.

Definition at line 502 of file Medium.hh.

502{ return d; }

Referenced by Diffusion().

◆ ScaleDiffusionTensor()

virtual double Garfield::Medium::ScaleDiffusionTensor ( const double d) const
inlinevirtual

Reimplemented in Garfield::MediumGas.

Definition at line 503 of file Medium.hh.

503{ return d; }

Referenced by Diffusion().

◆ ScaleDissociation()

virtual double Garfield::Medium::ScaleDissociation ( const double diss) const
inlinevirtual

Definition at line 507 of file Medium.hh.

507{ return diss; }

Referenced by IonDissociation().

◆ ScaleElectricField()

virtual double Garfield::Medium::ScaleElectricField ( const double e) const
inlinevirtual

Reimplemented in Garfield::MediumGas.

Definition at line 499 of file Medium.hh.

499{ return e; }

Referenced by Alpha(), Diffusion(), Diffusion(), ElectronLorentzAngle(), IonVelocity(), NegativeIonVelocity(), and Velocity().

◆ ScaleLorentzAngle()

virtual double Garfield::Medium::ScaleLorentzAngle ( const double lor) const
inlinevirtual

Reimplemented in Garfield::MediumGas.

Definition at line 506 of file Medium.hh.

506{ return lor; }

Referenced by ElectronLorentzAngle().

◆ ScaleTownsend()

virtual double Garfield::Medium::ScaleTownsend ( const double alpha) const
inlinevirtual

Reimplemented in Garfield::MediumGas.

Definition at line 504 of file Medium.hh.

504{ return alpha; }

Referenced by ElectronTownsend(), and HoleTownsend().

◆ ScaleVelocity()

virtual double Garfield::Medium::ScaleVelocity ( const double v) const
inlinevirtual

Definition at line 501 of file Medium.hh.

501{ return v; }

◆ SetAtomicNumber()

void Garfield::Medium::SetAtomicNumber ( const double z)
virtual

Set the effective atomic number.

Reimplemented in Garfield::MediumGas.

Definition at line 115 of file Medium.cc.

115 {
116 if (z < 1.) {
117 std::cerr << m_className << "::SetAtomicNumber:\n"
118 << " Atomic number must be >= 1.\n";
119 return;
120 }
121 m_z = z;
122 m_isChanged = true;
123}

Referenced by Garfield::MediumCdTe::MediumCdTe(), Garfield::MediumDiamond::MediumDiamond(), Garfield::MediumGaAs::MediumGaAs(), Garfield::MediumGaN::MediumGaN(), and Garfield::MediumSilicon::MediumSilicon().

◆ SetAtomicWeight()

void Garfield::Medium::SetAtomicWeight ( const double a)
virtual

Set the effective atomic weight.

Reimplemented in Garfield::MediumGas.

Definition at line 125 of file Medium.cc.

125 {
126 if (a <= 0.) {
127 std::cerr << m_className << "::SetAtomicWeight:\n"
128 << " Atomic weight must be greater than zero.\n";
129 return;
130 }
131 m_a = a;
132 m_isChanged = true;
133}

Referenced by Garfield::MediumCdTe::MediumCdTe(), Garfield::MediumDiamond::MediumDiamond(), Garfield::MediumGaAs::MediumGaAs(), Garfield::MediumGaN::MediumGaN(), and Garfield::MediumSilicon::MediumSilicon().

◆ SetDielectricConstant()

void Garfield::Medium::SetDielectricConstant ( const double eps)

Set the relative static dielectric constant.

Definition at line 92 of file Medium.cc.

92 {
93 if (eps < 1.) {
94 std::cerr << m_className << "::SetDielectricConstant:\n"
95 << " Dielectric constant must be >= 1.\n";
96 return;
97 }
98 m_epsilon = eps;
99 m_isChanged = true;
100}

Referenced by Garfield::MediumCdTe::MediumCdTe(), Garfield::MediumDiamond::MediumDiamond(), Garfield::MediumGaAs::MediumGaAs(), Garfield::MediumGaN::MediumGaN(), and Garfield::MediumSilicon::MediumSilicon().

◆ SetElectronAttachment()

bool Garfield::Medium::SetElectronAttachment ( const size_t ie,
const size_t ib,
const size_t ia,
const double eta )
inline

Set an entry in the table of attachment coefficients.

Definition at line 286 of file Medium.hh.

287 {
288 return SetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
289 }
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)
Definition Medium.cc:966

◆ SetElectronLongitudinalDiffusion()

bool Garfield::Medium::SetElectronLongitudinalDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
const double dl )
inline

Set an entry in the table of longitudinal diffusion coefficients.

Definition at line 256 of file Medium.hh.

257 {
258 return SetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
259 }

◆ SetElectronLorentzAngle()

bool Garfield::Medium::SetElectronLorentzAngle ( const size_t ie,
const size_t ib,
const size_t ia,
const double lor )
inline

Set an entry in the table of Lorentz angles.

Definition at line 297 of file Medium.hh.

298 {
299 return SetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
300 }

◆ SetElectronTownsend()

bool Garfield::Medium::SetElectronTownsend ( const size_t ie,
const size_t ib,
const size_t ia,
const double alpha )
inline

Set an entry in the table of Townsend coefficients.

Definition at line 276 of file Medium.hh.

277 {
278 return SetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
279 }

◆ SetElectronTransverseDiffusion()

bool Garfield::Medium::SetElectronTransverseDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
const double dt )
inline

Set an entry in the table of transverse diffusion coefficients.

Definition at line 266 of file Medium.hh.

267 {
268 return SetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
269 }

◆ SetElectronVelocityB()

bool Garfield::Medium::SetElectronVelocityB ( const size_t ie,
const size_t ib,
const size_t ia,
const double v )
inline

Set an entry in the table of drift speeds along Btrans.

Definition at line 246 of file Medium.hh.

247 {
248 return SetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
249 }

◆ SetElectronVelocityE()

bool Garfield::Medium::SetElectronVelocityE ( const size_t ie,
const size_t ib,
const size_t ia,
const double v )
inline

Set an entry in the table of drift speeds along E.

Definition at line 226 of file Medium.hh.

227 {
228 return SetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
229 }

◆ SetElectronVelocityExB()

bool Garfield::Medium::SetElectronVelocityExB ( const size_t ie,
const size_t ib,
const size_t ia,
const double v )
inline

Set an entry in the table of drift speeds along ExB.

Definition at line 236 of file Medium.hh.

237 {
238 return SetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
239 }

◆ SetEntry()

bool Garfield::Medium::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 )
protected

◆ SetExtrapolationMethod()

void Garfield::Medium::SetExtrapolationMethod ( const std::string & low,
const std::string & high,
std::pair< unsigned int, unsigned int > & extr,
const std::string & fcn )
protected

Definition at line 1215 of file Medium.cc.

1218 {
1219 unsigned int i = 0;
1220 if (GetExtrapolationIndex(low, i)) {
1221 extr.first = i;
1222 } else {
1223 std::cerr << m_className << "::SetExtrapolationMethod" << fcn << ":\n"
1224 << " Unknown extrapolation method (" << low << ")\n";
1225 }
1226 unsigned int j = 0;
1227 if (GetExtrapolationIndex(high, j)) {
1228 extr.second = j;
1229 } else {
1230 std::cerr << m_className << "::SetExtrapolationMethod" << fcn << ":\n"
1231 << " Unknown extrapolation method (" << high << ")\n";
1232 }
1233}
bool GetExtrapolationIndex(std::string str, unsigned int &nb) const
Definition Medium.cc:1235

Referenced by SetExtrapolationMethodAttachment(), SetExtrapolationMethodDiffusion(), Garfield::MediumGas::SetExtrapolationMethodExcitationRates(), SetExtrapolationMethodIonDissociation(), Garfield::MediumGas::SetExtrapolationMethodIonisationRates(), SetExtrapolationMethodIonMobility(), SetExtrapolationMethodTownsend(), and SetExtrapolationMethodVelocity().

◆ SetExtrapolationMethodAttachment()

void Garfield::Medium::SetExtrapolationMethodAttachment ( const std::string & extrLow,
const std::string & extrHigh )

Definition at line 1200 of file Medium.cc.

1201 {
1202 SetExtrapolationMethod(low, high, m_extrAtt, "Attachment");
1203}
void SetExtrapolationMethod(const std::string &low, const std::string &high, std::pair< unsigned int, unsigned int > &extr, const std::string &fcn)
Definition Medium.cc:1215

◆ SetExtrapolationMethodDiffusion()

void Garfield::Medium::SetExtrapolationMethodDiffusion ( const std::string & extrLow,
const std::string & extrHigh )

Definition at line 1190 of file Medium.cc.

1191 {
1192 SetExtrapolationMethod(low, high, m_extrDif, "Diffusion");
1193}

◆ SetExtrapolationMethodIonDissociation()

void Garfield::Medium::SetExtrapolationMethodIonDissociation ( const std::string & extrLow,
const std::string & extrHigh )

Definition at line 1210 of file Medium.cc.

1211 {
1212 SetExtrapolationMethod(low, high, m_extrDis, "IonDissociation");
1213}

◆ SetExtrapolationMethodIonMobility()

void Garfield::Medium::SetExtrapolationMethodIonMobility ( const std::string & extrLow,
const std::string & extrHigh )

Definition at line 1205 of file Medium.cc.

1206 {
1207 SetExtrapolationMethod(low, high, m_extrMob, "IonMobility");
1208}

◆ SetExtrapolationMethodTownsend()

void Garfield::Medium::SetExtrapolationMethodTownsend ( const std::string & extrLow,
const std::string & extrHigh )

Definition at line 1195 of file Medium.cc.

1196 {
1197 SetExtrapolationMethod(low, high, m_extrAlp, "Townsend");
1198}

◆ SetExtrapolationMethodVelocity()

void Garfield::Medium::SetExtrapolationMethodVelocity ( const std::string & extrLow,
const std::string & extrHigh )

Select the extrapolation method for fields below/above the table range. Possible options are "constant", "linear", and "exponential".

Definition at line 1185 of file Medium.cc.

1186 {
1187 SetExtrapolationMethod(low, high, m_extrVel, "Velocity");
1188}
std::pair< unsigned int, unsigned int > m_extrVel
Definition Medium.hh:615

◆ SetFanoFactor()

void Garfield::Medium::SetFanoFactor ( const double f)
inline

Set the Fano factor.

Definition at line 88 of file Medium.hh.

88{ m_fano = f; }

◆ SetFieldGrid() [1/2]

void Garfield::Medium::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.

Definition at line 882 of file Medium.cc.

884 {
885 const std::string hdr = m_className + "::SetFieldGrid";
886 if (!CheckFields(efields, hdr, "E-fields")) return;
887 if (!CheckFields(bfields, hdr, "B-fields")) return;
888 if (!CheckFields(angles, hdr, "angles")) return;
889
890 if (m_debug) {
891 std::cout << m_className << "::SetFieldGrid:\n E-fields:\n";
892 for (const auto efield : efields) std::cout << " " << efield << "\n";
893 std::cout << " B-fields:\n";
894 for (const auto bfield : bfields) std::cout << " " << bfield << "\n";
895 std::cout << " Angles:\n";
896 for (const auto angle : angles) std::cout << " " << angle << "\n";
897 }
898
899 // Clone the existing tables.
900 // Electrons
901 Clone(m_eVelE, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
902 "electron velocity along E");
903 Clone(m_eVelB, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
904 "electron velocity along Bt");
905 Clone(m_eVelX, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
906 "electron velocity along ExB");
907 Clone(m_eDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
908 "electron longitudinal diffusion");
909 Clone(m_eDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
910 "electron transverse diffusion");
911 Clone(m_eAlp, efields, bfields, angles, m_intpAlp, m_extrAlp, -30.,
912 "electron Townsend coefficient");
913 Clone(m_eAtt, efields, bfields, angles, m_intpAtt, m_extrAtt, -30.,
914 "electron attachment coefficient");
915 Clone(m_eLor, efields, bfields, angles, m_intpLor, m_extrLor, 0.,
916 "electron Lorentz angle");
917 if (!m_eDifM.empty()) {
918 Clone(m_eDifM, 6, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
919 "electron diffusion tensor");
920 }
921
922 // Holes
923 Clone(m_hVelE, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
924 "hole velocity along E");
925 Clone(m_hVelB, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
926 "hole velocity along Bt");
927 Clone(m_hVelX, efields, bfields, angles, m_intpVel, m_extrVel, 0.,
928 "hole velocity along ExB");
929 Clone(m_hDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
930 "hole longitudinal diffusion");
931 Clone(m_hDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
932 "hole transverse diffusion");
933 Clone(m_hAlp, efields, bfields, angles, m_intpAlp, m_extrAlp, -30.,
934 "hole Townsend coefficient");
935 Clone(m_hAtt, efields, bfields, angles, m_intpAtt, m_extrAtt, -30.,
936 "hole attachment coefficient");
937 if (!m_hDifM.empty()) {
938 Clone(m_hDifM, 6, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
939 "hole diffusion tensor");
940 }
941
942 // Ions
943 Clone(m_iMob, efields, bfields, angles, m_intpMob, m_extrMob, 0.,
944 "ion mobility");
945 Clone(m_iDifL, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
946 "ion longitudinal diffusion");
947 Clone(m_iDifT, efields, bfields, angles, m_intpDif, m_extrDif, 0.,
948 "ion transverse diffusion");
949 Clone(m_iDis, efields, bfields, angles, m_intpDis, m_extrDis, -30.,
950 "ion dissociation");
951
952 if (bfields.size() > 1 || angles.size() > 1) m_tab2d = true;
953 m_eFields = efields;
954 m_bFields = bfields;
955 m_bAngles = angles;
956}
unsigned int m_intpVel
Definition Medium.hh:624
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)
Definition Medium.cc:1018

◆ SetFieldGrid() [2/2]

void Garfield::Medium::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.

Definition at line 791 of file Medium.cc.

793 {
794 // Check if the requested E-field range makes sense.
795 if (ne <= 0) {
796 std::cerr << m_className << "::SetFieldGrid:\n"
797 << " Number of E-fields must be > 0.\n";
798 return;
799 }
800
801 if (emin < 0. || emax < 0.) {
802 std::cerr << m_className << "::SetFieldGrid:\n"
803 << " Electric fields must be positive.\n";
804 return;
805 }
806
807 if (emax < emin) {
808 std::cerr << m_className << "::SetFieldGrid: Swapping min./max. E-field.\n";
809 std::swap(emin, emax);
810 }
811
812 double estep = 0.;
813 if (logE) {
814 // Logarithmic scale
815 if (emin < Small) {
816 std::cerr << m_className << "::SetFieldGrid:\n"
817 << " Min. E-field must be non-zero for log. scale.\n";
818 return;
819 }
820 if (ne == 1) {
821 std::cerr << m_className << "::SetFieldGrid:\n"
822 << " Number of E-fields must be > 1 for log. scale.\n";
823 return;
824 }
825 estep = pow(emax / emin, 1. / (ne - 1.));
826 } else {
827 // Linear scale
828 if (ne > 1) estep = (emax - emin) / (ne - 1.);
829 }
830
831 // Check if the requested B-field range makes sense.
832 if (nb <= 0) {
833 std::cerr << m_className << "::SetFieldGrid:\n"
834 << " Number of B-fields must be > 0.\n";
835 return;
836 }
837 if (bmax < 0. || bmin < 0.) {
838 std::cerr << m_className << "::SetFieldGrid:\n"
839 << " Magnetic fields must be positive.\n";
840 return;
841 }
842 if (bmax < bmin) {
843 std::cerr << m_className << "::SetFieldGrid: Swapping min./max. B-field.\n";
844 std::swap(bmin, bmax);
845 }
846
847 const double bstep = nb > 1 ? (bmax - bmin) / (nb - 1.) : 0.;
848
849 // Check if the requested angular range makes sense.
850 if (na <= 0) {
851 std::cerr << m_className << "::SetFieldGrid:\n"
852 << " Number of angles must be > 0.\n";
853 return;
854 }
855 if (amax < 0. || amin < 0.) {
856 std::cerr << m_className << "::SetFieldGrid:"
857 << " Angles must be positive.\n";
858 return;
859 }
860 if (amax < amin) {
861 std::cerr << m_className << "::SetFieldGrid: Swapping min./max. angle.\n";
862 std::swap(amin, amax);
863 }
864 const double astep = na > 1 ? (amax - amin) / (na - 1.) : 0;
865
866 // Setup the field grids.
867 std::vector<double> eFields(ne);
868 std::vector<double> bFields(nb);
869 std::vector<double> bAngles(na);
870 for (size_t i = 0; i < ne; ++i) {
871 eFields[i] = logE ? emin * pow(estep, i) : emin + i * estep;
872 }
873 for (size_t i = 0; i < nb; ++i) {
874 bFields[i] = bmin + i * bstep;
875 }
876 for (size_t i = 0; i < na; ++i) {
877 bAngles[i] = amin + i * astep;
878 }
879 SetFieldGrid(eFields, bFields, bAngles);
880}
DoubleAc pow(const DoubleAc &f, double p)
Definition DoubleAc.cpp:337

Referenced by Medium(), and SetFieldGrid().

◆ SetHoleAttachment()

bool Garfield::Medium::SetHoleAttachment ( const size_t ie,
const size_t ib,
const size_t ia,
const double eta )
inline

Set an entry in the table of attachment coefficients.

Definition at line 369 of file Medium.hh.

370 {
371 return SetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
372 }

◆ SetHoleLongitudinalDiffusion()

bool Garfield::Medium::SetHoleLongitudinalDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
const double dl )
inline

Set an entry in the table of longitudinal diffusion coefficients.

Definition at line 339 of file Medium.hh.

340 {
341 return SetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
342 }

◆ SetHoleTownsend()

bool Garfield::Medium::SetHoleTownsend ( const size_t ie,
const size_t ib,
const size_t ia,
const double alpha )
inline

Set an entry in the table of Townsend coefficients.

Definition at line 359 of file Medium.hh.

360 {
361 return SetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
362 }

◆ SetHoleTransverseDiffusion()

bool Garfield::Medium::SetHoleTransverseDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
const double dt )
inline

Set an entry in the table of transverse diffusion coefficients.

Definition at line 349 of file Medium.hh.

350 {
351 return SetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
352 }

◆ SetHoleVelocityB()

bool Garfield::Medium::SetHoleVelocityB ( const size_t ie,
const size_t ib,
const size_t ia,
const double v )
inline

Set an entry in the table of drift speeds along Btrans.

Definition at line 328 of file Medium.hh.

329 {
330 return SetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
331 }

◆ SetHoleVelocityE()

bool Garfield::Medium::SetHoleVelocityE ( const size_t ie,
const size_t ib,
const size_t ia,
const double v )
inline

Set an entry in the table of drift speeds along E.

Definition at line 308 of file Medium.hh.

309 {
310 return SetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
311 }

◆ SetHoleVelocityExB()

bool Garfield::Medium::SetHoleVelocityExB ( const size_t ie,
const size_t ib,
const size_t ia,
const double v )
inline

Set an entry in the table of drift speeds along ExB.

Definition at line 318 of file Medium.hh.

319 {
320 return SetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
321 }

◆ SetInterpolationMethodAttachment()

void Garfield::Medium::SetInterpolationMethodAttachment ( const unsigned int intrp)

Definition at line 1288 of file Medium.cc.

1288 {
1289 if (intrp > 0) m_intpAtt = intrp;
1290}

◆ SetInterpolationMethodDiffusion()

void Garfield::Medium::SetInterpolationMethodDiffusion ( const unsigned int intrp)

Definition at line 1280 of file Medium.cc.

1280 {
1281 if (intrp > 0) m_intpDif = intrp;
1282}

◆ SetInterpolationMethodIonDissociation()

void Garfield::Medium::SetInterpolationMethodIonDissociation ( const unsigned int intrp)

Definition at line 1296 of file Medium.cc.

1296 {
1297 if (intrp > 0) m_intpDis = intrp;
1298}

◆ SetInterpolationMethodIonMobility()

void Garfield::Medium::SetInterpolationMethodIonMobility ( const unsigned int intrp)

Definition at line 1292 of file Medium.cc.

1292 {
1293 if (intrp > 0) m_intpMob = intrp;
1294}

◆ SetInterpolationMethodTownsend()

void Garfield::Medium::SetInterpolationMethodTownsend ( const unsigned int intrp)

Definition at line 1284 of file Medium.cc.

1284 {
1285 if (intrp > 0) m_intpAlp = intrp;
1286}

◆ SetInterpolationMethodVelocity()

void Garfield::Medium::SetInterpolationMethodVelocity ( const unsigned int intrp)

Set the degree of polynomial interpolation (usually 2).

Definition at line 1276 of file Medium.cc.

1276 {
1277 if (intrp > 0) m_intpVel = intrp;
1278}

◆ SetIonDissociation()

bool Garfield::Medium::SetIonDissociation ( const size_t ie,
const size_t ib,
const size_t ia,
const double diss )
inline

Set an entry in the table of dissociation coefficients.

Definition at line 416 of file Medium.hh.

417 {
418 return SetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
419 }

◆ SetIonLongitudinalDiffusion()

bool Garfield::Medium::SetIonLongitudinalDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
const double dl )
inline

Set an entry in the table of longitudinal diffusion coefficients.

Definition at line 396 of file Medium.hh.

397 {
398 return SetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
399 }

◆ SetIonMobility() [1/2]

bool Garfield::Medium::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.

Definition at line 1109 of file Medium.cc.

1110 {
1111 // Check the index.
1112 if (ie >= m_eFields.size() || ib >= m_bFields.size() ||
1113 ia >= m_bAngles.size()) {
1114 PrintOutOfRange(m_className, "SetIonMobility", ie, ib, ia);
1115 return false;
1116 }
1117
1118 if (m_iMob.empty()) {
1119 std::cerr << m_className << "::SetIonMobility:\n"
1120 << " Ion mobility table not initialised.\n";
1121 return false;
1122 }
1123
1124 if (mu == 0.) {
1125 std::cerr << m_className << "::SetIonMobility: Zero value not allowed.\n";
1126 return false;
1127 }
1128
1129 m_iMob[ia][ib][ie] = mu;
1130 if (m_debug) {
1131 std::cout << m_className << "::SetIonMobility:\n Ion mobility at E = "
1132 << m_eFields[ie] << " V/cm, B = "
1133 << m_bFields[ib] << " T, angle "
1134 << m_bAngles[ia] << " set to " << mu << " cm2/(V ns).\n";
1135 }
1136 return true;
1137}

◆ SetIonMobility() [2/2]

bool Garfield::Medium::SetIonMobility ( const std::vector< double > & fields,
const std::vector< double > & mobilities,
const bool negativeIons = false )

Initialise the table of ion mobilities from a list of electric fields and corresponding mobilities. The mobilities will be interpolated at the electric fields of the currently set grid.

Definition at line 1139 of file Medium.cc.

1141 {
1142 if (efields.size() != mobs.size()) {
1143 std::cerr << m_className << "::SetIonMobility:\n"
1144 << " E-field and mobility arrays have different sizes.\n";
1145 return false;
1146 }
1147
1148 if (negativeIons) {
1150 } else {
1152 }
1153 const auto nE = m_eFields.size();
1154 const auto nB = m_bFields.size();
1155 const auto nA = m_bAngles.size();
1156 if (negativeIons) {
1157 Init(nE, nB, nA, m_nMob, 0.);
1158 } else {
1159 Init(nE, nB, nA, m_iMob, 0.);
1160 }
1161 for (size_t i = 0; i < nE; ++i) {
1162 const double e = m_eFields[i];
1163 const double mu = Interpolate1D(e, mobs, efields, m_intpMob, m_extrMob);
1164 if (negativeIons) {
1165 m_nMob[0][0][i] = mu;
1166 } else {
1167 m_iMob[0][0][i] = mu;
1168 }
1169 }
1170 if (!m_tab2d) return true;
1171 for (size_t i = 0; i < nA; ++i) {
1172 for (size_t j = 0; j < nB; ++j) {
1173 for (size_t k = 0; k < nE; ++k) {
1174 if (negativeIons) {
1175 m_nMob[i][j][k] = m_nMob[0][0][k];
1176 } else {
1177 m_iMob[i][j][k] = m_iMob[0][0][k];
1178 }
1179 }
1180 }
1181 }
1182 return true;
1183}

Referenced by Garfield::MediumGas::LoadMobility().

◆ SetIonTransverseDiffusion()

bool Garfield::Medium::SetIonTransverseDiffusion ( const size_t ie,
const size_t ib,
const size_t ia,
const double dt )
inline

Set an entry in the table of transverse diffusion coefficients.

Definition at line 406 of file Medium.hh.

407 {
408 return SetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
409 }

◆ SetMassDensity()

void Garfield::Medium::SetMassDensity ( const double rho)
virtual

Set the mass density [g/cm3].

Reimplemented in Garfield::MediumGas.

Definition at line 145 of file Medium.cc.

145 {
146 if (rho <= 0.) {
147 std::cerr << m_className << "::SetMassDensity:\n"
148 << " Density [g/cm3] must be greater than zero.\n";
149 return;
150 }
151
152 if (m_a <= 0.) {
153 std::cerr << m_className << "::SetMassDensity:\n"
154 << " Atomic weight is not defined.\n";
155 return;
156 }
157 m_density = rho / (AtomicMassUnit * m_a);
158 m_isChanged = true;
159}

Referenced by Garfield::MediumCdTe::MediumCdTe(), Garfield::MediumDiamond::MediumDiamond(), Garfield::MediumGaAs::MediumGaAs(), Garfield::MediumGaN::MediumGaN(), and Garfield::MediumSilicon::MediumSilicon().

◆ SetNegativeIonMobility()

bool Garfield::Medium::SetNegativeIonMobility ( const size_t ie,
const size_t ib,
const size_t ia,
const double mu )
inline

Set an entry in the table of negative ion mobilities.

Definition at line 427 of file Medium.hh.

428 {
429 return SetEntry(ie, ib, ia, "NegativeIonMobility", m_nMob, mu);
430 }

◆ SetNumberDensity()

void Garfield::Medium::SetNumberDensity ( const double n)
virtual

Set the number density [cm-3].

Reimplemented in Garfield::MediumGas.

Definition at line 135 of file Medium.cc.

135 {
136 if (n <= 0.) {
137 std::cerr << m_className << "::SetNumberDensity:\n"
138 << " Density [cm-3] must be greater than zero.\n";
139 return;
140 }
141 m_density = n;
142 m_isChanged = true;
143}

◆ SetPressure()

void Garfield::Medium::SetPressure ( const double p)

Definition at line 82 of file Medium.cc.

82 {
83 if (p <= 0.) {
84 std::cerr << m_className << "::SetPressure:\n"
85 << " Pressure [Torr] must be greater than zero.\n";
86 return;
87 }
88 m_pressure = p;
89 m_isChanged = true;
90}

Referenced by main().

◆ SetTemperature()

void Garfield::Medium::SetTemperature ( const double t)

Set the temperature [K].

Definition at line 72 of file Medium.cc.

72 {
73 if (t <= 0.) {
74 std::cerr << m_className << "::SetTemperature:\n"
75 << " Temperature [K] must be greater than zero.\n";
76 return;
77 }
78 m_temperature = t;
79 m_isChanged = true;
80}

Referenced by main(), Garfield::MediumCdTe::MediumCdTe(), Garfield::MediumDiamond::MediumDiamond(), Garfield::MediumGaAs::MediumGaAs(), Garfield::MediumGaN::MediumGaN(), and Garfield::MediumSilicon::MediumSilicon().

◆ SetThreshold()

size_t Garfield::Medium::SetThreshold ( const std::vector< std::vector< std::vector< double > > > & tab) const
protected

Definition at line 1252 of file Medium.cc.

1253 {
1254
1255 if (tab.empty()) return 0;
1256 const auto nE = m_eFields.size();
1257 const auto nB = m_bFields.size();
1258 const auto nA = m_bAngles.size();
1259 for (size_t i = 0; i < nE; ++i) {
1260 bool below = false;
1261 for (size_t k = 0; k < nA; ++k) {
1262 for (size_t j = 0; j < nB; ++j) {
1263 if (tab[k][j][i] < -20.) {
1264 below = true;
1265 break;
1266 }
1267 }
1268 if (below) break;
1269 }
1270 if (below) continue;
1271 return i;
1272 }
1273 return nE - 1;
1274}

Referenced by Garfield::MediumGas::AdjustTownsendCoefficient(), Garfield::MediumMagboltz::GenerateGasTable(), and Garfield::MediumGas::MergeGasFile().

◆ SetW()

void Garfield::Medium::SetW ( const double w)
inline

Set the W value (average energy to produce an electron/ion or e/h pair).

Definition at line 84 of file Medium.hh.

84{ m_w = w; }

◆ UnScaleElectricField()

virtual double Garfield::Medium::UnScaleElectricField ( const double e) const
inlinevirtual

Reimplemented in Garfield::MediumGas.

Definition at line 500 of file Medium.hh.

500{ return e; }

Referenced by ElectronMobility(), and HoleMobility().

◆ Velocity()

bool Garfield::Medium::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
protected

Definition at line 196 of file Medium.cc.

201 {
202
203 vx = vy = vz = 0.;
204 // Make sure there is at least a table of velocities along E.
205 if (velE.empty()) return false;
206
207 // Compute the magnitude of the electric field.
208 const double e = sqrt(ex * ex + ey * ey + ez * ez);
209 const double e0 = ScaleElectricField(e);
210 if (e < Small || e0 < Small) return false;
211
212 // Compute the magnitude of the magnetic field.
213 const double b = sqrt(bx * bx + by * by + bz * bz);
214 // Compute the angle between B field and E field.
215 const double ebang = GetAngle(ex, ey, ez, bx, by, bz, e, b);
216
217 // Calculate the velocity along E.
218 double ve = 0.;
219 if (!Interpolate(e0, b, ebang, velE, ve, m_intpVel, m_extrVel)) {
220 std::cerr << m_className << "::Velocity: Interpolation along E failed.\n";
221 return false;
222 }
223 if (b < Small) {
224 // No magnetic field.
225 const double mu = q * ve / e;
226 vx = mu * ex;
227 vy = mu * ey;
228 vz = mu * ez;
229 return true;
230 } else if (velX.empty() || velB.empty() ||
231 (m_bFields.size() == 1 && fabs(m_bFields[0]) < Small)) {
232 // Magnetic field, velocities along ExB, Bt not available.
233 Langevin(ex, ey, ez, bx, by, bz, q * ve / e, vx, vy, vz);
234 return true;
235 }
236
237 // Magnetic field, velocities along ExB and Bt available.
238 // Compute unit vectors along E, E x B and Bt.
239 double ue[3] = {ex / e, ey / e, ez / e};
240 double uexb[3] = {ey * bz - ez * by, ez * bx - ex * bz, ex * by - ey * bx};
241 const double exb =
242 sqrt(uexb[0] * uexb[0] + uexb[1] * uexb[1] + uexb[2] * uexb[2]);
243 if (exb > 0.) {
244 uexb[0] /= exb;
245 uexb[1] /= exb;
246 uexb[2] /= exb;
247 } else {
248 uexb[0] = ue[0];
249 uexb[1] = ue[1];
250 uexb[2] = ue[2];
251 }
252
253 double ubt[3] = {uexb[1] * ez - uexb[2] * ey,
254 uexb[2] * ex - uexb[0] * ez,
255 uexb[0] * ey - uexb[1] * ex};
256 const double bt = sqrt(ubt[0] * ubt[0] + ubt[1] * ubt[1] + ubt[2] * ubt[2]);
257 if (bt > 0.) {
258 ubt[0] /= bt;
259 ubt[1] /= bt;
260 ubt[2] /= bt;
261 } else {
262 ubt[0] = ue[0];
263 ubt[1] = ue[1];
264 ubt[2] = ue[2];
265 }
266
267 if (m_debug) {
268 std::cout << m_className << "::Velocity:\n";
269 std::printf(" unit vector along E: (%15.5f, %15.5f, %15.5f)\n",
270 ue[0], ue[1], ue[2]);
271 std::printf(" unit vector along E x B: (%15.5f, %15.5f, %15.5f)\n",
272 uexb[0], uexb[1], uexb[2]);
273 std::printf(" unit vector along Bt: (%15.5f, %15.5f, %15.5f)\n",
274 ubt[0], ubt[1], ubt[2]);
275 }
276
277 // Calculate the velocities in all directions.
278 double vexb = 0.;
279 if (!Interpolate(e0, b, ebang, velX, vexb, m_intpVel, m_extrVel)) {
280 std::cerr << m_className << "::Velocity: Interpolation along ExB failed.\n";
281 return false;
282 }
283 vexb *= q;
284 double vbt = 0.;
285 if (!Interpolate(e0, b, ebang, velB, vbt, m_intpVel, m_extrVel)) {
286 std::cerr << m_className << "::Velocity: Interpolation along Bt failed.\n";
287 return false;
288 }
289 if (ex * bx + ey * by + ez * bz > 0.) {
290 vbt = fabs(vbt);
291 } else {
292 vbt = -fabs(vbt);
293 }
294 vbt *= q * q;
295 vx = q * (ve * ue[0] + vbt * ubt[0] + vexb * uexb[0]);
296 vy = q * (ve * ue[1] + vbt * ubt[1] + vexb * uexb[1]);
297 vz = q * (ve * ue[2] + vbt * ubt[2] + vexb * uexb[2]);
298 return true;
299}

Referenced by ElectronVelocity(), and HoleVelocity().

Member Data Documentation

◆ m_a

double Garfield::Medium::m_a = 0.
protected

Definition at line 548 of file Medium.hh.

Referenced by GetAtomicWeight(), GetMassDensity(), SetAtomicWeight(), and SetMassDensity().

◆ m_bAngles

◆ m_bFields

◆ m_className

std::string Garfield::Medium::m_className = "Medium"
protected

Definition at line 529 of file Medium.hh.

Referenced by Garfield::MediumGas::AdjustTownsendCoefficient(), Clone(), Clone(), Garfield::MediumMagboltz::ComputeDeexcitation(), Garfield::MediumGas::DisablePenningTransfer(), Garfield::MediumMagboltz::DisablePenningTransfer(), Garfield::MediumSilicon::ElectronAttachment(), ElectronCollision(), Garfield::MediumMagboltz::ElectronCollision(), Garfield::MediumSilicon::ElectronCollision(), Garfield::MediumSilicon::ElectronTownsend(), Garfield::MediumSilicon::ElectronVelocity(), Garfield::MediumMagboltz::EnableDeexcitation(), Garfield::MediumGas::EnablePenningTransfer(), Garfield::MediumGas::EnablePenningTransfer(), Garfield::MediumGas::EnablePenningTransfer(), Garfield::MediumMagboltz::EnablePenningTransfer(), Garfield::MediumMagboltz::EnablePenningTransfer(), Garfield::MediumMagboltz::EnableRadiationTrapping(), Garfield::MediumMagboltz::GenerateGasTable(), GetComponent(), Garfield::MediumCdTe::GetComponent(), Garfield::MediumDiamond::GetComponent(), Garfield::MediumGas::GetComponent(), Garfield::MediumSilicon::GetConductionBandDensityOfStates(), GetDeexcitationProduct(), GetDielectricFunction(), Garfield::MediumSilicon::GetDielectricFunction(), Garfield::MediumSilicon::GetElectronBandPopulation(), GetElectronCollisionRate(), Garfield::MediumMagboltz::GetElectronCollisionRate(), Garfield::MediumMagboltz::GetElectronCollisionRate(), Garfield::MediumSilicon::GetElectronCollisionRate(), GetElectronEnergy(), Garfield::MediumSilicon::GetElectronEnergy(), Garfield::MediumGas::GetElectronExcitationRate(), Garfield::MediumGas::GetElectronIonisationRate(), Garfield::MediumSilicon::GetElectronMomentum(), GetElectronNullCollisionRate(), Garfield::MediumSilicon::GetElectronNullCollisionRate(), GetEntry(), Garfield::MediumGas::GetExcitationLevel(), Garfield::MediumGas::GetIonisationLevel(), Garfield::MediumMagboltz::GetLevel(), Garfield::MediumGas::GetMixture(), Garfield::MediumMagboltz::GetNumberOfElectronCollisions(), Garfield::MediumSilicon::GetNumberOfElectronCollisions(), GetOpticalDataRange(), Garfield::MediumSilicon::GetOpticalDataRange(), Garfield::MediumGas::GetPenningTransfer(), GetPhotoAbsorptionCrossSection(), Garfield::MediumGas::GetPhotoAbsorptionCrossSection(), Garfield::MediumMagboltz::GetPhotonCollision(), Garfield::MediumMagboltz::GetPhotonCollisionRate(), Garfield::MediumSilicon::GetValenceBandDensityOfStates(), Garfield::MediumSilicon::HoleAttachment(), Garfield::MediumSilicon::HoleTownsend(), Garfield::MediumSilicon::HoleVelocity(), Init(), Init(), Garfield::MediumMagboltz::Initialise(), Garfield::MediumSilicon::Initialise(), Interpolate1D(), Garfield::MediumGas::LoadGasFile(), Garfield::MediumGas::LoadMobility(), Garfield::MediumCdTe::MediumCdTe(), Garfield::MediumConductor::MediumConductor(), Garfield::MediumDiamond::MediumDiamond(), Garfield::MediumGaAs::MediumGaAs(), Garfield::MediumGaN::MediumGaN(), Garfield::MediumGas::MediumGas(), Garfield::MediumMagboltz::MediumMagboltz(), Garfield::MediumPlastic::MediumPlastic(), Garfield::MediumSilicon::MediumSilicon(), Garfield::MediumGas::MergeGasFile(), Garfield::MediumGas::PrintGas(), Garfield::MediumGas::ReadHeader(), Garfield::MediumMagboltz::RunMagboltz(), SetAtomicNumber(), Garfield::MediumGas::SetAtomicNumber(), SetAtomicWeight(), Garfield::MediumGas::SetAtomicWeight(), Garfield::MediumGas::SetComposition(), SetDielectricConstant(), Garfield::MediumSilicon::SetDoping(), Garfield::MediumGaN::SetElectronConcentration(), SetEntry(), Garfield::MediumMagboltz::SetExcitationScaling(), SetExtrapolationMethod(), SetFieldGrid(), SetFieldGrid(), SetIonMobility(), SetIonMobility(), Garfield::MediumCdTe::SetLowFieldMobility(), Garfield::MediumDiamond::SetLowFieldMobility(), Garfield::MediumGaAs::SetLowFieldMobility(), Garfield::MediumGaN::SetLowFieldMobility(), Garfield::MediumSilicon::SetLowFieldMobility(), SetMassDensity(), Garfield::MediumGas::SetMassDensity(), Garfield::MediumMagboltz::SetMaxElectronEnergy(), Garfield::MediumSilicon::SetMaxElectronEnergy(), Garfield::MediumMagboltz::SetMaxPhotonEnergy(), SetNumberDensity(), Garfield::MediumGas::SetNumberDensity(), SetPressure(), Garfield::MediumDiamond::SetSaturationVelocity(), Garfield::MediumSilicon::SetSaturationVelocity(), Garfield::MediumMagboltz::SetSplittingFunctionGreenSawada(), SetTemperature(), Garfield::MediumSilicon::SetTrapCrossSection(), Garfield::MediumSilicon::SetTrapDensity(), Garfield::MediumSilicon::SetTrappingTime(), Velocity(), and Garfield::MediumGas::WriteGasFile().

◆ m_debug

◆ m_density

double Garfield::Medium::m_density = 0.
protected

◆ m_driftable

◆ m_eAlp

◆ m_eAtt

◆ m_eDifL

◆ m_eDifM

◆ m_eDifT

◆ m_eFields

◆ m_eLor

◆ m_epsilon

double Garfield::Medium::m_epsilon = 1.
protected

◆ m_eThrAlp

unsigned int Garfield::Medium::m_eThrAlp = 0
protected

◆ m_eThrAtt

unsigned int Garfield::Medium::m_eThrAtt = 0
protected

◆ m_eVelB

◆ m_eVelE

◆ m_eVelX

◆ m_extrAlp

std::pair<unsigned int, unsigned int> Garfield::Medium::m_extrAlp = {0, 1}
protected

◆ m_extrAtt

std::pair<unsigned int, unsigned int> Garfield::Medium::m_extrAtt = {0, 1}
protected

◆ m_extrDif

std::pair<unsigned int, unsigned int> Garfield::Medium::m_extrDif = {0, 1}
protected

◆ m_extrDis

std::pair<unsigned int, unsigned int> Garfield::Medium::m_extrDis = {0, 1}
protected

◆ m_extrLor

std::pair<unsigned int, unsigned int> Garfield::Medium::m_extrLor = {0, 1}
protected

◆ m_extrMob

◆ m_extrVel

std::pair<unsigned int, unsigned int> Garfield::Medium::m_extrVel = {0, 1}
protected

◆ m_fano

◆ m_hAlp

◆ m_hAtt

◆ m_hDifL

std::vector<std::vector<std::vector<double> > > Garfield::Medium::m_hDifL
protected

◆ m_hDifM

std::vector<std::vector<std::vector<std::vector<double> > > > Garfield::Medium::m_hDifM
protected

Definition at line 597 of file Medium.hh.

Referenced by HoleDiffusion(), ResetHoleDiffusion(), and SetFieldGrid().

◆ m_hDifT

std::vector<std::vector<std::vector<double> > > Garfield::Medium::m_hDifT
protected

◆ m_hThrAlp

unsigned int Garfield::Medium::m_hThrAlp = 0
protected

Definition at line 610 of file Medium.hh.

Referenced by HoleTownsend().

◆ m_hThrAtt

unsigned int Garfield::Medium::m_hThrAtt = 0
protected

Definition at line 611 of file Medium.hh.

Referenced by HoleAttachment().

◆ m_hVelB

std::vector<std::vector<std::vector<double> > > Garfield::Medium::m_hVelB
protected

◆ m_hVelE

◆ m_hVelX

std::vector<std::vector<std::vector<double> > > Garfield::Medium::m_hVelX
protected

◆ m_id

int Garfield::Medium::m_id
protected

Definition at line 534 of file Medium.hh.

Referenced by GetId(), and Medium().

◆ m_idCounter

int Garfield::Medium::m_idCounter = -1
staticprotected

Definition at line 531 of file Medium.hh.

Referenced by Medium().

◆ m_iDifL

◆ m_iDifT

◆ m_iDis

◆ m_iMob

◆ m_intpAlp

◆ m_intpAtt

◆ m_intpDif

◆ m_intpDis

◆ m_intpLor

◆ m_intpMob

◆ m_intpVel

◆ m_ionisable

◆ m_isChanged

bool Garfield::Medium::m_isChanged = true
protected

Definition at line 563 of file Medium.hh.

Referenced by Garfield::MediumSilicon::ElectronAttachment(), Garfield::MediumSilicon::ElectronCollision(), Garfield::MediumGaAs::ElectronTownsend(), Garfield::MediumGaN::ElectronTownsend(), Garfield::MediumSilicon::ElectronTownsend(), Garfield::MediumCdTe::ElectronVelocity(), Garfield::MediumDiamond::ElectronVelocity(), Garfield::MediumGaAs::ElectronVelocity(), Garfield::MediumGaN::ElectronVelocity(), Garfield::MediumSilicon::ElectronVelocity(), Garfield::MediumMagboltz::EnableAnisotropicScattering(), Garfield::MediumMagboltz::EnableDeexcitation(), Garfield::MediumMagboltz::EnableRadiationTrapping(), Garfield::MediumSilicon::GetElectronCollisionRate(), Garfield::MediumSilicon::GetElectronNullCollisionRate(), Garfield::MediumSilicon::HoleAttachment(), Garfield::MediumGaAs::HoleTownsend(), Garfield::MediumGaN::HoleTownsend(), Garfield::MediumSilicon::HoleTownsend(), Garfield::MediumCdTe::HoleVelocity(), Garfield::MediumDiamond::HoleVelocity(), Garfield::MediumGaAs::HoleVelocity(), Garfield::MediumGaN::HoleVelocity(), Garfield::MediumSilicon::HoleVelocity(), Garfield::MediumMagboltz::Initialise(), Garfield::MediumSilicon::Initialise(), Garfield::MediumGas::LoadGasFile(), Garfield::MediumGas::MediumGas(), Garfield::MediumMagboltz::MediumMagboltz(), Garfield::MediumMagboltz::PrintGas(), SetAtomicNumber(), SetAtomicWeight(), Garfield::MediumGas::SetComposition(), SetDielectricConstant(), Garfield::MediumSilicon::SetDoping(), Garfield::MediumSilicon::SetDopingMobilityModelMasetti(), Garfield::MediumSilicon::SetDopingMobilityModelMinimos(), Garfield::MediumMagboltz::SetExcitationScaling(), Garfield::MediumSilicon::SetHighFieldMobilityModelCanali(), Garfield::MediumSilicon::SetHighFieldMobilityModelMinimos(), Garfield::MediumSilicon::SetHighFieldMobilityModelReggiani(), Garfield::MediumSilicon::SetImpactIonisationModelGrant(), Garfield::MediumSilicon::SetImpactIonisationModelMassey(), Garfield::MediumSilicon::SetImpactIonisationModelOkutoCrowell(), Garfield::MediumSilicon::SetImpactIonisationModelVanOverstraetenDeMan(), Garfield::MediumSilicon::SetLatticeMobilityModelMinimos(), Garfield::MediumSilicon::SetLatticeMobilityModelReggiani(), Garfield::MediumSilicon::SetLatticeMobilityModelSentaurus(), Garfield::MediumCdTe::SetLowFieldMobility(), Garfield::MediumDiamond::SetLowFieldMobility(), Garfield::MediumGaAs::SetLowFieldMobility(), Garfield::MediumGaN::SetLowFieldMobility(), Garfield::MediumSilicon::SetLowFieldMobility(), SetMassDensity(), Garfield::MediumMagboltz::SetMaxElectronEnergy(), Garfield::MediumSilicon::SetMaxElectronEnergy(), Garfield::MediumMagboltz::SetMaxPhotonEnergy(), SetNumberDensity(), SetPressure(), Garfield::MediumSilicon::SetSaturationVelocity(), Garfield::MediumSilicon::SetSaturationVelocityModelCanali(), Garfield::MediumSilicon::SetSaturationVelocityModelMinimos(), Garfield::MediumSilicon::SetSaturationVelocityModelReggiani(), Garfield::MediumMagboltz::SetSplittingFunctionGreenSawada(), SetTemperature(), Garfield::MediumSilicon::SetTrapCrossSection(), Garfield::MediumSilicon::SetTrapDensity(), Garfield::MediumSilicon::SetTrappingTime(), Garfield::MediumCdTe::UnsetLowFieldMobility(), Garfield::MediumDiamond::UnsetLowFieldMobility(), Garfield::MediumGaAs::UnsetLowFieldMobility(), and Garfield::MediumGaN::UnsetLowFieldMobility().

◆ m_iThrDis

unsigned int Garfield::Medium::m_iThrDis = 0
protected

◆ m_microscopic

◆ m_name

◆ m_nComponents

◆ m_nMob

std::vector<std::vector<std::vector<double> > > Garfield::Medium::m_nMob
protected

◆ m_pressure

◆ m_tab2d

◆ m_temperature

◆ m_w

◆ m_z

double Garfield::Medium::m_z = 1.
protected

Definition at line 546 of file Medium.hh.

Referenced by GetAtomicNumber(), and SetAtomicNumber().


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