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

#include <MediumGaAs.hh>

+ Inheritance diagram for Garfield::MediumGaAs:

Public Member Functions

 MediumGaAs ()
 
 ~MediumGaAs ()
 
bool IsSemiconductor () const
 
void GetComponent (const unsigned int &i, std::string &label, double &f)
 
void SetTrapCrossSection (const double ecs, const double hcs)
 
void SetTrapDensity (const double n)
 
void SetTrappingTime (const double etau, const double htau)
 
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)
 
bool ElectronTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
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)
 
bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
void SetLowFieldMobility (const double mue, const double muh)
 
bool GetOpticalDataRange (double &emin, double &emax, const unsigned int &i=0)
 
bool GetDielectricFunction (const double &e, double &eps1, double &eps2, const unsigned int &i=0)
 
- Public Member Functions inherited from Garfield::Medium
 Medium ()
 
virtual ~Medium ()
 
int GetId () const
 
std::string GetName () const
 
virtual bool IsGas () const
 
virtual bool IsSemiconductor () const
 
void SetTemperature (const double &t)
 
double GetTemperature () const
 
void SetPressure (const double &p)
 
double GetPressure () const
 
void SetDielectricConstant (const double &eps)
 
double GetDielectricConstant () const
 
unsigned int GetNumberOfComponents () const
 
virtual void GetComponent (const unsigned int &i, std::string &label, double &f)
 
virtual void SetAtomicNumber (const double &z)
 
virtual double GetAtomicNumber () const
 
virtual void SetAtomicWeight (const double &a)
 
virtual double GetAtomicWeight () const
 
virtual void SetNumberDensity (const double &n)
 
virtual double GetNumberDensity () const
 
virtual void SetMassDensity (const double &rho)
 
virtual double GetMassDensity () const
 
virtual void EnableDrift ()
 
void DisableDrift ()
 
virtual void EnablePrimaryIonisation ()
 
void DisablePrimaryIonisation ()
 
bool IsDriftable () const
 
bool IsMicroscopic () const
 
bool IsIonisable () const
 
void SetW (const double &w)
 
double GetW ()
 
void SetFanoFactor (const double &f)
 
double GetFanoFactor ()
 
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)
 
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)
 
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)
 
virtual bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
virtual double GetElectronEnergy (const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
 
virtual void GetElectronMomentum (const double e, double &px, double &py, double &pz, int &band)
 
virtual double GetElectronNullCollisionRate (const int band=0)
 
virtual double GetElectronCollisionRate (const double e, const int band=0)
 
virtual bool GetElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, int &nion, int &ndxc, int &band)
 
virtual int GetNumberOfIonisationProducts ()
 
virtual bool GetIonisationProduct (const int i, int &type, double &energy)
 
virtual int GetNumberOfDeexcitationProducts ()
 
virtual bool GetDeexcitationProduct (const int i, double &t, double &s, int &type, double &energy)
 
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)
 
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)
 
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])
 
virtual bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
virtual bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
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)
 
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)
 
virtual bool IonDissociation (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &diss)
 
void SetFieldGrid (double emin, double emax, int ne, bool logE, double bmin=0., double bmax=0., int nb=1, double amin=0., double amax=0., int na=1)
 
void SetFieldGrid (const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles)
 
void GetFieldGrid (std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles)
 
bool GetElectronVelocityE (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &v)
 
bool GetElectronVelocityExB (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &v)
 
bool GetElectronVelocityB (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &v)
 
bool GetElectronLongitudinalDiffusion (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &dl)
 
bool GetElectronTransverseDiffusion (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &dt)
 
bool GetElectronTownsend (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &alpha)
 
bool GetElectronAttachment (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &eta)
 
bool GetHoleVelocityE (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &v)
 
bool GetHoleVelocityExB (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &v)
 
bool GetHoleVelocityB (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &v)
 
bool GetHoleLongitudinalDiffusion (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &dl)
 
bool GetHoleTransverseDiffusion (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &dt)
 
bool GetHoleTownsend (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &alpha)
 
bool GetHoleAttachment (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &eta)
 
bool GetIonMobility (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &mu)
 
bool GetIonLongitudinalDiffusion (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &dl)
 
bool GetIonTransverseDiffusion (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &dt)
 
bool GetIonDissociation (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, double &diss)
 
void ResetElectronVelocity ()
 
void ResetElectronDiffusion ()
 
void ResetElectronTownsend ()
 
void ResetElectronAttachment ()
 
void ResetHoleVelocity ()
 
void ResetHoleDiffusion ()
 
void ResetHoleTownsend ()
 
void ResetHoleAttachment ()
 
void ResetIonMobility ()
 
void ResetIonDiffusion ()
 
void ResetIonDissociation ()
 
bool SetIonMobility (const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, const double &mu)
 
bool SetIonMobility (const std::vector< double > &fields, const std::vector< double > &mobilities)
 
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)
 
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 ScaleDissociation (const double &diss) const
 
virtual bool GetOpticalDataRange (double &emin, double &emax, const unsigned int &i=0)
 
virtual bool GetDielectricFunction (const double &e, double &eps1, double &eps2, const unsigned int &i=0)
 
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 ()
 
void DisableDebugging ()
 

Additional Inherited Members

- Protected Member Functions inherited from Garfield::Medium
double Interpolate1D (const double &e, const std::vector< double > &table, const std::vector< double > &fields, const unsigned int &intpMeth, const int &jExtr, const int &iExtr)
 
bool GetExtrapolationIndex (std::string extrStr, unsigned int &extrNb)
 
void CloneTable (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 unsigned int &extrLow, const unsigned int &extrHigh, const double init, const std::string label)
 
void CloneTensor (std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const unsigned int &n, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int &intp, const unsigned int &extrLow, const unsigned int &extrHigh, const double &init, const std::string &label)
 
void InitParamArrays (const unsigned int &eRes, const unsigned int &bRes, const unsigned int &aRes, std::vector< std::vector< std::vector< double > > > &tab, const double &val)
 
void InitParamTensor (const unsigned int &eRes, const unsigned int &bRes, const unsigned int &aRes, const unsigned int &tRes, std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const double &val)
 
- Protected Attributes inherited from Garfield::Medium
std::string m_className
 
int m_id
 
std::string m_name
 
double m_temperature
 
double m_pressure
 
double m_epsilon
 
unsigned int m_nComponents
 
double m_z
 
double m_a
 
double m_density
 
bool m_driftable
 
bool m_microscopic
 
bool m_ionisable
 
double m_w
 
double m_fano
 
bool m_isChanged
 
bool m_debug
 
unsigned int m_nEfields
 
unsigned int m_nBfields
 
unsigned int m_nAngles
 
std::vector< double > eFields
 
std::vector< double > bFields
 
std::vector< double > bAngles
 
bool m_map2d
 
bool m_hasElectronVelocityE
 
bool m_hasElectronVelocityB
 
bool m_hasElectronVelocityExB
 
bool m_hasElectronDiffLong
 
bool m_hasElectronDiffTrans
 
bool m_hasElectronDiffTens
 
bool m_hasElectronTownsend
 
bool m_hasElectronAttachment
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityE
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityExB
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityB
 
std::vector< std::vector< std::vector< double > > > tabElectronDiffLong
 
std::vector< std::vector< std::vector< double > > > tabElectronDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabElectronTownsend
 
std::vector< std::vector< std::vector< double > > > tabElectronAttachment
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabElectronDiffTens
 
bool m_hasHoleVelocityE
 
bool m_hasHoleVelocityB
 
bool m_hasHoleVelocityExB
 
bool m_hasHoleDiffLong
 
bool m_hasHoleDiffTrans
 
bool m_hasHoleDiffTens
 
bool m_hasHoleTownsend
 
bool m_hasHoleAttachment
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityE
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityExB
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityB
 
std::vector< std::vector< std::vector< double > > > tabHoleDiffLong
 
std::vector< std::vector< std::vector< double > > > tabHoleDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabHoleTownsend
 
std::vector< std::vector< std::vector< double > > > tabHoleAttachment
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabHoleDiffTens
 
bool m_hasIonMobility
 
bool m_hasIonDiffLong
 
bool m_hasIonDiffTrans
 
bool m_hasIonDissociation
 
std::vector< std::vector< std::vector< double > > > tabIonMobility
 
std::vector< std::vector< std::vector< double > > > tabIonDiffLong
 
std::vector< std::vector< std::vector< double > > > tabIonDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabIonDissociation
 
int thrElectronTownsend
 
int thrElectronAttachment
 
int thrHoleTownsend
 
int thrHoleAttachment
 
int thrIonDissociation
 
unsigned int m_extrLowVelocity
 
unsigned int m_extrHighVelocity
 
unsigned int m_extrLowDiffusion
 
unsigned int m_extrHighDiffusion
 
unsigned int m_extrLowTownsend
 
unsigned int m_extrHighTownsend
 
unsigned int m_extrLowAttachment
 
unsigned int m_extrHighAttachment
 
unsigned int m_extrLowMobility
 
unsigned int m_extrHighMobility
 
unsigned int m_extrLowDissociation
 
unsigned int m_extrHighDissociation
 
unsigned int m_intpVelocity
 
unsigned int m_intpDiffusion
 
unsigned int m_intpTownsend
 
unsigned int m_intpAttachment
 
unsigned int m_intpMobility
 
unsigned int m_intpDissociation
 
- Static Protected Attributes inherited from Garfield::Medium
static int m_idCounter = -1
 

Detailed Description

Definition at line 9 of file MediumGaAs.hh.

Constructor & Destructor Documentation

◆ MediumGaAs()

Garfield::MediumGaAs::MediumGaAs ( )

Definition at line 14 of file MediumGaAs.cc.

15 : Medium(),
16 // bandGap(1.42),
17 eMobility(8.8e-6),
18 hMobility(3.2e-6),
19 eHallFactor(1.05),
20 hHallFactor(1.25),
21 eTrapCs(1.e-15),
22 hTrapCs(1.e-15),
23 eTrapDensity(1.e13),
24 hTrapDensity(1.e13),
25 eTrapTime(0.),
26 hTrapTime(0.),
27 trappingModel(0),
28 m_hasUserMobility(false),
29 m_hasOpticalData(false),
30 opticalDataFile("OpticalData_GaAs.txt") {
31
32 m_className = "MediumGaAs";
33 m_name = "GaAs";
34
35 SetTemperature(300.);
38 SetAtomicWeight(72.32);
39 SetMassDensity(5.317);
40
43 m_microscopic = false;
44
45 m_w = 4.35;
46 m_fano = 0.1;
47}
bool m_microscopic
Definition: Medium.hh:309
virtual void SetAtomicWeight(const double &a)
Definition: Medium.cc:178
virtual void SetMassDensity(const double &rho)
Definition: Medium.cc:200
double m_fano
Definition: Medium.hh:313
std::string m_name
Definition: Medium.hh:291
virtual void SetAtomicNumber(const double &z)
Definition: Medium.cc:167
virtual void EnableDrift()
Definition: Medium.hh:52
void SetTemperature(const double &t)
Definition: Medium.cc:117
void SetDielectricConstant(const double &eps)
Definition: Medium.cc:139
virtual void EnablePrimaryIonisation()
Definition: Medium.hh:54
std::string m_className
Definition: Medium.hh:284

◆ ~MediumGaAs()

Garfield::MediumGaAs::~MediumGaAs ( )
inline

Definition at line 15 of file MediumGaAs.hh.

15{}

Member Function Documentation

◆ ElectronAttachment()

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

Reimplemented from Garfield::Medium.

Definition at line 158 of file MediumGaAs.cc.

161 {
162
163 eta = 0.;
165 // Interpolation in user table.
166 return Medium::ElectronAttachment(ex, ey, ez, bx, by, bz, eta);
167 }
168
169 switch (trappingModel) {
170 case 0:
171 eta = eTrapCs * eTrapDensity;
172 break;
173 case 1:
174 double vx, vy, vz;
175 ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
176 eta = eTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
177 if (eta > 0.) eta = 1. / eta;
178 break;
179 default:
180 std::cerr << m_className << "::ElectronAttachment:\n";
181 std::cerr << " Unknown model activated. Program bug!\n";
182 return false;
183 break;
184 }
185
186 return true;
187}
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:313
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)
Definition: MediumGaAs.cc:114
virtual bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Definition: Medium.cc:609
bool m_hasElectronAttachment
Definition: Medium.hh:335

◆ ElectronTownsend()

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

Reimplemented from Garfield::Medium.

Definition at line 145 of file MediumGaAs.cc.

148 {
149
150 alpha = 0.;
152 // Interpolation in user table.
153 return Medium::ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
154 }
155 return false;
156}
bool m_hasElectronTownsend
Definition: Medium.hh:335
virtual bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Definition: Medium.cc:542

◆ ElectronVelocity()

bool Garfield::MediumGaAs::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

Reimplemented from Garfield::Medium.

Definition at line 114 of file MediumGaAs.cc.

117 {
118
119 vx = vy = vz = 0.;
121 // Interpolation in user table.
122 return Medium::ElectronVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
123 }
124 // Calculate the mobility
125 double mu = eMobility;
126 mu = -mu;
127 const double b = sqrt(bx * bx + by * by + bz * bz);
128 if (b < Small) {
129 vx = mu * ex;
130 vy = mu * ey;
131 vz = mu * ez;
132 } else {
133 // Hall mobility
134 const double muH = eHallFactor * mu;
135 const double eb = bx * ex + by * ey + bz * ez;
136 const double nom = 1. + pow(muH * b, 2);
137 // Compute the drift velocity using the Langevin equation.
138 vx = mu * (ex + muH * (ey * bz - ez * by) + muH * muH * bx * eb) / nom;
139 vy = mu * (ey + muH * (ez * bx - ex * bz) + muH * muH * by * eb) / nom;
140 vz = mu * (ez + muH * (ex * by - ey * bx) + muH * muH * bz * eb) / nom;
141 }
142 return true;
143}
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:336
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)
Definition: Medium.cc:217
bool m_hasElectronVelocityE
Definition: Medium.hh:333

Referenced by ElectronAttachment().

◆ GetComponent()

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

Reimplemented from Garfield::Medium.

Definition at line 49 of file MediumGaAs.cc.

49 {
50
51 if (i == 0) {
52 label = "Ga";
53 f = 0.5;
54 } else if (i == 1) {
55 label = "As";
56 f = 0.5;
57 }
58}

◆ GetDielectricFunction()

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

Reimplemented from Garfield::Medium.

Definition at line 299 of file MediumGaAs.cc.

300 {
301
302 if (i != 0) {
303 std::cerr << m_className << "::GetDielectricFunction:\n";
304 std::cerr << " Medium has only one component.\n";
305 return false;
306 }
307
308 // Make sure the optical data table has been loaded.
309 if (!m_hasOpticalData) {
310 if (!LoadOpticalData(opticalDataFile)) {
311 std::cerr << m_className << "::GetDielectricFunction:\n";
312 std::cerr << " Optical data table could not be loaded.\n";
313 return false;
314 }
315 m_hasOpticalData = true;
316 }
317
318 // Make sure the requested energy is within the range of the table.
319 const double emin = opticalDataTable[0].energy;
320 const double emax = opticalDataTable.back().energy;
321 if (e < emin || e > emax) {
322 std::cerr << m_className << "::GetDielectricFunction:\n";
323 std::cerr << " Requested energy (" << e << " eV) "
324 << " is outside the range of the optical data table.\n";
325 std::cerr << " " << emin << " < E [eV] < " << emax << "\n";
326 eps1 = eps2 = 0.;
327 return false;
328 }
329
330 // Locate the requested energy in the table.
331 int iLow = 0;
332 int iUp = opticalDataTable.size() - 1;
333 int iM;
334 while (iUp - iLow > 1) {
335 iM = (iUp + iLow) >> 1;
336 if (e >= opticalDataTable[iM].energy) {
337 iLow = iM;
338 } else {
339 iUp = iM;
340 }
341 }
342
343 // Interpolate the real part of dielectric function.
344 // Use linear interpolation if one of the values is negative,
345 // Otherwise use log-log interpolation.
346 const double logX0 = log(opticalDataTable[iLow].energy);
347 const double logX1 = log(opticalDataTable[iUp].energy);
348 const double logX = log(e);
349 if (opticalDataTable[iLow].eps1 <= 0. || opticalDataTable[iUp].eps1 <= 0.) {
350 eps1 = opticalDataTable[iLow].eps1 +
351 (e - opticalDataTable[iLow].energy) *
352 (opticalDataTable[iUp].eps1 - opticalDataTable[iLow].eps1) /
353 (opticalDataTable[iUp].energy - opticalDataTable[iLow].energy);
354 } else {
355 const double logY0 = log(opticalDataTable[iLow].eps1);
356 const double logY1 = log(opticalDataTable[iUp].eps1);
357 eps1 = logY0 + (logX - logX0) * (logY1 - logY0) / (logX1 - logX0);
358 eps1 = exp(eps1);
359 }
360
361 // Interpolate the imaginary part of dielectric function,
362 // using log-log interpolation.
363 const double logY0 = log(opticalDataTable[iLow].eps2);
364 const double logY1 = log(opticalDataTable[iUp].eps2);
365 eps2 = logY0 + (log(e) - logX0) * (logY1 - logY0) / (logX1 - logX0);
366 eps2 = exp(eps2);
367 return true;
368}
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:376

◆ GetOpticalDataRange()

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

Reimplemented from Garfield::Medium.

Definition at line 272 of file MediumGaAs.cc.

273 {
274
275 if (i != 0) {
276 std::cerr << m_className << "::GetOpticalDataRange:\n";
277 std::cerr << " Medium has only one component.\n";
278 }
279
280 // Make sure the optical data table has been loaded.
281 if (!m_hasOpticalData) {
282 if (!LoadOpticalData(opticalDataFile)) {
283 std::cerr << m_className << "::GetOpticalDataRange:\n";
284 std::cerr << " Optical data table could not be loaded.\n";
285 return false;
286 }
287 m_hasOpticalData = true;
288 }
289
290 emin = opticalDataTable[0].energy;
291 emax = opticalDataTable.back().energy;
292 if (m_debug) {
293 std::cout << m_className << "::GetOpticalDataRange:\n";
294 std::cout << " " << emin << " < E [eV] < " << emax << "\n";
295 }
296 return true;
297}

◆ HoleAttachment()

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

Reimplemented from Garfield::Medium.

Definition at line 230 of file MediumGaAs.cc.

232 {
233
234 eta = 0.;
236 // Interpolation in user table.
237 return Medium::HoleAttachment(ex, ey, ez, bx, by, bz, eta);
238 }
239 switch (trappingModel) {
240 case 0:
241 eta = hTrapCs * hTrapDensity;
242 break;
243 case 1:
244 double vx, vy, vz;
245 HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
246 eta = hTrapTime * sqrt(vx * vx + vy * vy + vz * vz);
247 if (eta > 0.) eta = 1. / eta;
248 break;
249 default:
250 std::cerr << m_className << "::HoleAttachment:\n";
251 std::cerr << " Unknown model activated. Program bug!\n";
252 return false;
253 break;
254 }
255 return true;
256}
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)
Definition: MediumGaAs.cc:189
bool m_hasHoleAttachment
Definition: Medium.hh:350
virtual bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Definition: Medium.cc:1144

◆ HoleTownsend()

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

Reimplemented from Garfield::Medium.

Definition at line 218 of file MediumGaAs.cc.

220 {
221
222 alpha = 0.;
223 if (m_hasHoleTownsend) {
224 // Interpolation in user table.
225 return Medium::HoleTownsend(ex, ey, ez, bx, by, bz, alpha);
226 }
227 return false;
228}
virtual bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Definition: Medium.cc:1079
bool m_hasHoleTownsend
Definition: Medium.hh:350

◆ HoleVelocity()

bool Garfield::MediumGaAs::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

Reimplemented from Garfield::Medium.

Definition at line 189 of file MediumGaAs.cc.

191 {
192
193 vx = vy = vz = 0.;
194 if (m_hasHoleVelocityE) {
195 // Interpolation in user table.
196 return Medium::HoleVelocity(ex, ey, ez, bx, by, bz, vx, vy, vz);
197 }
198 // Calculate the mobility
199 double mu = hMobility;
200 const double b = sqrt(bx * bx + by * by + bz * bz);
201 if (b < Small) {
202 vx = mu * ex;
203 vy = mu * ey;
204 vz = mu * ez;
205 } else {
206 // Hall mobility
207 const double muH = hHallFactor * mu;
208 const double eb = bx * ex + by * ey + bz * ez;
209 const double nom = 1. + pow(muH * b, 2);
210 // Compute the drift velocity using the Langevin equation.
211 vx = mu * (ex + muH * (ey * bz - ez * by) + muH * muH * bx * eb) / nom;
212 vy = mu * (ey + muH * (ez * bx - ex * bz) + muH * muH * by * eb) / nom;
213 vz = mu * (ez + muH * (ex * by - ey * bx) + muH * muH * bz * eb) / nom;
214 }
215 return true;
216}
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)
Definition: Medium.cc:773
bool m_hasHoleVelocityE
Definition: Medium.hh:348

Referenced by HoleAttachment().

◆ IsSemiconductor()

bool Garfield::MediumGaAs::IsSemiconductor ( ) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 17 of file MediumGaAs.hh.

17{ return true; }

◆ SetLowFieldMobility()

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

Definition at line 258 of file MediumGaAs.cc.

258 {
259
260 if (mue <= 0. || muh <= 0.) {
261 std::cerr << m_className << "::SetLowFieldMobility:\n";
262 std::cerr << " Mobility must be greater than zero.\n";
263 return;
264 }
265
266 eMobility = mue;
267 hMobility = muh;
268 m_hasUserMobility = true;
269 m_isChanged = true;
270}
bool m_isChanged
Definition: Medium.hh:316

◆ SetTrapCrossSection()

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

Definition at line 60 of file MediumGaAs.cc.

60 {
61
62 if (ecs < 0.) {
63 std::cerr << m_className << "::SetTrapCrossSection:\n";
64 std::cerr << " Capture cross-section [cm2] must positive.\n";
65 } else {
66 eTrapCs = ecs;
67 }
68
69 if (hcs < 0.) {
70 std::cerr << m_className << "::SetTrapCrossSection:\n";
71 std::cerr << " Capture cross-section [cm2] must be positive.n";
72 } else {
73 hTrapCs = hcs;
74 }
75
76 trappingModel = 0;
77 m_isChanged = true;
78}

◆ SetTrapDensity()

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

Definition at line 80 of file MediumGaAs.cc.

80 {
81
82 if (n < 0.) {
83 std::cerr << m_className << "::SetTrapDensity:\n";
84 std::cerr << " Trap density [cm-3] must be greater than zero.\n";
85 } else {
86 eTrapDensity = n;
87 hTrapDensity = n;
88 }
89
90 trappingModel = 0;
91 m_isChanged = true;
92}

◆ SetTrappingTime()

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

Definition at line 94 of file MediumGaAs.cc.

94 {
95
96 if (etau <= 0.) {
97 std::cerr << m_className << "::SetTrappingTime:\n";
98 std::cerr << " Trapping time [ns-1] must be positive.\n";
99 } else {
100 eTrapTime = etau;
101 }
102
103 if (htau <= 0.) {
104 std::cerr << m_className << "::SetTrappingTime:\n";
105 std::cerr << " Trapping time [ns-1] must be positive.\n";
106 } else {
107 hTrapTime = htau;
108 }
109
110 trappingModel = 1;
111 m_isChanged = true;
112}

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