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

#include <MediumMagboltz.hh>

+ Inheritance diagram for Garfield::MediumMagboltz:

Public Member Functions

 MediumMagboltz ()
 
virtual ~MediumMagboltz ()
 
bool SetMaxElectronEnergy (const double e)
 
double GetMaxElectronEnergy () const
 
bool SetMaxPhotonEnergy (const double e)
 
double GetMaxPhotonEnergy () const
 
void EnableEnergyRangeAdjustment ()
 
void DisableEnergyRangeAdjustment ()
 
void EnableAnisotropicScattering ()
 
void DisableAnisotropicScattering ()
 
void SetSplittingFunctionOpalBeaty ()
 
void SetSplittingFunctionGreenSawada ()
 
void SetSplittingFunctionFlat ()
 
void EnableDeexcitation ()
 
void DisableDeexcitation ()
 
void EnableRadiationTrapping ()
 
void DisableRadiationTrapping ()
 
void EnablePenningTransfer (const double r, const double lambda)
 
void EnablePenningTransfer (const double r, const double lambda, std::string gasname)
 
void DisablePenningTransfer ()
 
void DisablePenningTransfer (std::string gasname)
 
void EnableCrossSectionOutput ()
 
void DisableCrossSectionOutput ()
 
void SetExcitationScalingFactor (const double r, std::string gasname)
 
bool Initialise (const bool verbose=false)
 
void PrintGas ()
 
double GetElectronNullCollisionRate (const int band)
 
double GetElectronCollisionRate (const double e, const int band)
 
double GetElectronCollisionRate (const double e, const unsigned int level, const int band)
 
bool GetElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, int &nion, int &ndxc, int &band)
 
unsigned int GetNumberOfIonisationProducts () const
 
bool GetIonisationProduct (const unsigned int i, int &type, double &energy) const
 
void ComputeDeexcitation (int iLevel, int &fLevel)
 
unsigned int GetNumberOfDeexcitationProducts () const
 
bool GetDeexcitationProduct (const unsigned int i, double &t, double &s, int &type, double &energy) const
 
double GetPhotonCollisionRate (const double e)
 
bool GetPhotonCollision (const double e, int &type, int &level, double &e1, double &ctheta, int &nsec, double &esec)
 
void ResetCollisionCounters ()
 
unsigned int GetNumberOfElectronCollisions () const
 
unsigned int GetNumberOfElectronCollisions (int &nElastic, int &nIonising, int &nAttachment, int &nInelastic, int &nExcitation, int &nSuperelastic) const
 
int GetNumberOfLevels ()
 
bool GetLevel (const unsigned int i, int &ngas, int &type, std::string &descr, double &e)
 
unsigned int GetNumberOfElectronCollisions (const unsigned int level) const
 
int GetNumberOfPenningTransfers () const
 
int GetNumberOfPhotonCollisions () const
 
int GetNumberOfPhotonCollisions (int &nElastic, int &nIonising, int &nInelastic) const
 
void RunMagboltz (const double e, const double b, const double btheta, const int ncoll, bool verbose, double &vx, double &vy, double &vz, double &dl, double &dt, double &alpha, double &eta, double &lor, double &vxerr, double &vyerr, double &vzerr, double &dlerr, double &dterr, double &alphaerr, double &etaerr, double &lorerr, double &alphatof)
 
void GenerateGasTable (const int numCollisions=10, const bool verbose=true)
 
- Public Member Functions inherited from Garfield::MediumGas
 MediumGas ()
 
virtual ~MediumGas ()
 
bool IsGas () const
 
bool SetComposition (const std::string &gas1, const double f1=1., const std::string &gas2="", const double f2=0., const std::string &gas3="", const double f3=0., const std::string &gas4="", const double f4=0., const std::string &gas5="", const double f5=0., const std::string &gas6="", const double f6=0.)
 
void GetComposition (std::string &gas1, double &f1, std::string &gas2, double &f2, std::string &gas3, double &f3, std::string &gas4, double &f4, std::string &gas5, double &f5, std::string &gas6, double &f6)
 
void GetComponent (const unsigned int i, std::string &label, double &f)
 
void SetAtomicNumber (const double z)
 
double GetAtomicNumber () const
 
void SetAtomicWeight (const double a)
 
double GetAtomicWeight () const
 
void SetNumberDensity (const double n)
 
double GetNumberDensity () const
 
void SetMassDensity (const double rho)
 
double GetMassDensity () const
 
bool LoadGasFile (const std::string &filename)
 
bool WriteGasFile (const std::string &filename)
 
void PrintGas ()
 
bool LoadIonMobility (const std::string &filename)
 
void SetExtrapolationMethodExcitationRates (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodIonisationRates (const std::string &extrLow, const std::string &extrHigh)
 
void SetInterpolationMethodExcitationRates (const int intrp)
 
void SetInterpolationMethodIonisationRates (const int intrp)
 
double ScaleElectricField (const double e) const
 
double UnScaleElectricField (const double e) const
 
double ScaleDiffusion (const double d) const
 
double ScaleDiffusionTensor (const double d) const
 
double ScaleTownsend (const double alpha) const
 
double ScaleAttachment (const double eta) const
 
double ScaleLorentzAngle (const double lor) const
 
bool GetPhotoabsorptionCrossSection (const double &e, double &sigma, const unsigned int &i)
 
- Public Member Functions inherited from Garfield::Medium
 Medium ()
 
virtual ~Medium ()
 
int GetId () const
 
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 bool ElectronLorentzAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
 
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 unsigned int GetNumberOfIonisationProducts () const
 
virtual bool GetIonisationProduct (const unsigned int i, int &type, double &energy) const
 
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)
 
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 GetElectronLorentzAngle (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &lor)
 
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 ResetElectronLorentzAngle ()
 
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 ScaleLorentzAngle (const double lor) 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 ()
 

Public Attributes

double fit3d4p
 
double fitHigh4p
 
double fit3dQCO2
 
double fit3dEtaCO2
 
double fit3dQCH4
 
double fit3dEtaCH4
 
double fit3dQC2H6
 
double fit3dEtaC2H6
 
double fit4pEtaCH4
 
double fit4pEtaC2H6
 
double fit4sEtaC2H6
 
double fitLineCut
 

Additional Inherited Members

- Protected Member Functions inherited from Garfield::MediumGas
bool GetGasInfo (const std::string &gasname, double &a, double &z) const
 
bool GetGasName (const int gasnumber, const int version, std::string &gasname)
 
bool GetGasName (std::string input, std::string &gasname) const
 
bool GetGasNumberGasFile (const std::string &input, int &number) const
 
- Protected Member Functions inherited from Garfield::Medium
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
 
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::MediumGas
std::string m_gas [m_nMaxGases]
 
double m_fraction [m_nMaxGases]
 
double m_atWeight [m_nMaxGases]
 
double m_atNum [m_nMaxGases]
 
bool m_usePenning
 
double m_rPenningGlobal
 
double m_rPenningGas [m_nMaxGases]
 
double m_lambdaPenningGlobal
 
double m_lambdaPenningGas [m_nMaxGases]
 
double m_pressureTable
 
double m_temperatureTable
 
std::vector< std::vector< std::vector< double > > > m_tabTownsendNoPenning
 
bool m_hasExcRates
 
bool m_hasIonRates
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_tabExcRates
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_tabIonRates
 
std::vector< excListElementm_excitationList
 
std::vector< ionListElementm_ionisationList
 
unsigned int m_extrLowExcRates
 
unsigned int m_extrHighExcRates
 
unsigned int m_extrLowIonRates
 
unsigned int m_extrHighIonRates
 
unsigned int m_intpExcRates
 
unsigned int m_intpIonRates
 
- 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
 
std::vector< double > m_eFields
 
std::vector< double > m_bFields
 
std::vector< double > m_bAngles
 
bool m_map2d
 
bool m_hasElectronVelocityE
 
bool m_hasElectronVelocityB
 
bool m_hasElectronVelocityExB
 
bool m_hasElectronDiffLong
 
bool m_hasElectronDiffTrans
 
bool m_hasElectronDiffTens
 
bool m_hasElectronAttachment
 
bool m_hasElectronLorentzAngle
 
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< double > > > tabElectronLorentzAngle
 
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_extrLowLorentzAngle
 
unsigned int m_extrHighLorentzAngle
 
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_intpLorentzAngle
 
unsigned int m_intpMobility
 
unsigned int m_intpDissociation
 
- Static Protected Attributes inherited from Garfield::MediumGas
static const unsigned int m_nMaxGases = 6
 
- Static Protected Attributes inherited from Garfield::Medium
static int m_idCounter = -1
 

Detailed Description

Interface to Magboltz (version 9).

Definition at line 11 of file MediumMagboltz.hh.

Constructor & Destructor Documentation

◆ MediumMagboltz()

Garfield::MediumMagboltz::MediumMagboltz ( )

Definition at line 23 of file MediumMagboltz.cc.

24 : MediumGas(),
25 m_eFinal(40.),
26 m_eStep(m_eFinal / nEnergySteps),
27 m_eHigh(1.e4),
28 m_eHighLog(log(m_eHigh)),
29 m_lnStep(1.),
30 m_useAutoAdjust(true),
31 m_useCsOutput(false),
32 m_nTerms(0),
33 m_useAnisotropic(true),
34 m_nPenning(0),
35 m_useDeexcitation(false),
36 m_useRadTrap(true),
37 m_useOpalBeaty(true),
38 m_useGreenSawada(false),
39 m_eFinalGamma(20.),
40 m_eStepGamma(m_eFinalGamma / nEnergyStepsGamma) {
41
42 fit3d4p = fitHigh4p = 1.;
46 fit4sEtaC2H6 = 0.5;
47 fitLineCut = 1000;
48
49 m_className = "MediumMagboltz";
50
51 // Set physical constants in Magboltz common blocks.
52 Magboltz::cnsts_.echarg = ElementaryCharge * 1.e-15;
53 Magboltz::cnsts_.emass = ElectronMassGramme;
54 Magboltz::cnsts_.amu = AtomicMassUnit;
55 Magboltz::cnsts_.pir2 = BohrRadius * BohrRadius * Pi;
56 Magboltz::inpt_.ary = RydbergEnergy;
57
58 // Set parameters in Magboltz common blocks.
60 Magboltz::inpt_.nStep = nEnergySteps;
61 // Select the scattering model.
62 Magboltz::inpt_.nAniso = 2;
63 // Max. energy [eV]
64 Magboltz::inpt_.efinal = m_eFinal;
65 // Energy step size [eV]
66 Magboltz::inpt_.estep = m_eStep;
67 // Temperature and pressure
68 Magboltz::inpt_.akt = BoltzmannConstant * m_temperature;
69 Magboltz::inpt_.tempc = m_temperature - ZeroCelsius;
71 // Disable Penning transfer.
72 Magboltz::inpt_.ipen = 0;
73
74 // Initialise Penning parameters
75 for (int i = nMaxLevels; i--;) {
76 m_rPenning[i] = 0.;
77 m_lambdaPenning[i] = 0.;
78 }
79
80 m_isChanged = true;
81
84 m_microscopic = true;
85
86 // Initialize the collision counters.
87 m_nCollisionsDetailed.clear();
88 for (int i = nCsTypes; i--;) m_nCollisions[i] = 0;
89 for (int i = nCsTypesGamma; i--;) m_nPhotonCollisions[i] = 0;
90
91 m_ionProducts.clear();
92 m_dxcProducts.clear();
93
94 for (unsigned int i = 0; i < m_nMaxGases; ++i) m_scaleExc[i] = 1.;
95}
static const unsigned int m_nMaxGases
Definition: MediumGas.hh:87
bool m_microscopic
Definition: Medium.hh:319
double m_pressure
Definition: Medium.hh:305
virtual void EnableDrift()
Definition: Medium.hh:52
virtual void EnablePrimaryIonisation()
Definition: Medium.hh:54
unsigned int m_nComponents
Definition: Medium.hh:309
std::string m_className
Definition: Medium.hh:294
bool m_isChanged
Definition: Medium.hh:326
double m_temperature
Definition: Medium.hh:303
struct Garfield::Magboltz::@1 inpt_
struct Garfield::Magboltz::@3 cnsts_

◆ ~MediumMagboltz()

virtual Garfield::MediumMagboltz::~MediumMagboltz ( )
inlinevirtual

Definition at line 17 of file MediumMagboltz.hh.

17{}

Member Function Documentation

◆ ComputeDeexcitation()

void Garfield::MediumMagboltz::ComputeDeexcitation ( int  iLevel,
int &  fLevel 
)

Definition at line 4981 of file MediumMagboltz.cc.

4981 {
4982
4983 if (!m_useDeexcitation) {
4984 std::cerr << m_className << "::ComputeDeexcitation:\n";
4985 std::cerr << " Deexcitation is disabled.\n";
4986 return;
4987 }
4988
4989 // Make sure that the tables are updated.
4990 if (m_isChanged) {
4991 if (!Mixer()) {
4992 std::cerr << m_className << "::ComputeDeexcitation:\n";
4993 std::cerr << " Error calculating the collision rates table.\n";
4994 return;
4995 }
4996 m_isChanged = false;
4997 }
4998
4999 if (iLevel < 0 || iLevel >= (int)m_nTerms) {
5000 std::cerr << m_className << "::ComputeDeexcitation:\n";
5001 std::cerr << " Level index is out of range.\n";
5002 return;
5003 }
5004
5005 iLevel = m_iDeexcitation[iLevel];
5006 if (iLevel < 0 || iLevel >= (int)m_deexcitations.size()) {
5007 std::cerr << m_className << "::ComputeDeexcitation:\n";
5008 std::cerr << " Level is not deexcitable.\n";
5009 return;
5010 }
5011
5012 ComputeDeexcitationInternal(iLevel, fLevel);
5013 if (fLevel >= 0 && fLevel < (int)m_deexcitations.size()) {
5014 fLevel = m_deexcitations[fLevel].level;
5015 }
5016}

◆ DisableAnisotropicScattering()

void Garfield::MediumMagboltz::DisableAnisotropicScattering ( )
inline

Definition at line 39 of file MediumMagboltz.hh.

39 {
40 m_useAnisotropic = false;
41 m_isChanged = true;
42 }

◆ DisableCrossSectionOutput()

void Garfield::MediumMagboltz::DisableCrossSectionOutput ( )
inline

Definition at line 67 of file MediumMagboltz.hh.

67{ m_useCsOutput = false; }

◆ DisableDeexcitation()

void Garfield::MediumMagboltz::DisableDeexcitation ( )
inline

Definition at line 51 of file MediumMagboltz.hh.

51{ m_useDeexcitation = false; }

◆ DisableEnergyRangeAdjustment()

void Garfield::MediumMagboltz::DisableEnergyRangeAdjustment ( )
inline

Definition at line 32 of file MediumMagboltz.hh.

32{ m_useAutoAdjust = false; }

◆ DisablePenningTransfer() [1/2]

void Garfield::MediumMagboltz::DisablePenningTransfer ( )

Definition at line 313 of file MediumMagboltz.cc.

313 {
314
315 for (unsigned int i = 0; i < m_nTerms; ++i) {
316 m_rPenning[i] = 0.;
317 m_lambdaPenning[i] = 0.;
318 }
319 m_rPenningGlobal = 0.;
321
322 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
323 m_rPenningGas[i] = 0.;
324 m_lambdaPenningGas[i] = 0.;
325 }
326
327 m_usePenning = false;
328}
double m_rPenningGlobal
Definition: MediumGas.hh:99
double m_lambdaPenningGlobal
Definition: MediumGas.hh:102
double m_rPenningGas[m_nMaxGases]
Definition: MediumGas.hh:100
double m_lambdaPenningGas[m_nMaxGases]
Definition: MediumGas.hh:103

◆ DisablePenningTransfer() [2/2]

void Garfield::MediumMagboltz::DisablePenningTransfer ( std::string  gasname)

Definition at line 330 of file MediumMagboltz.cc.

330 {
331
332 // Get the "standard" name of this gas.
333 if (!GetGasName(gasname, gasname)) {
334 std::cerr << m_className << "::DisablePenningTransfer:\n";
335 std::cerr << " Gas " << gasname << " is not defined.\n";
336 return;
337 }
338
339 // Look for this gas in the present gas mixture.
340 bool found = false;
341 int iGas = -1;
342 for (unsigned int i = 0; i < m_nComponents; ++i) {
343 if (m_gas[i] == gasname) {
344 m_rPenningGas[i] = 0.;
345 m_lambdaPenningGas[i] = 0.;
346 found = true;
347 iGas = i;
348 break;
349 }
350 }
351
352 if (!found) {
353 std::cerr << m_className << "::DisablePenningTransfer:\n";
354 std::cerr << " Specified gas (" << gasname
355 << ") is not part of the present gas mixture.\n";
356 return;
357 }
358
359 unsigned int nLevelsFound = 0;
360 for (unsigned int i = 0; i < m_nTerms; ++i) {
361 if (int(m_csType[i] / nCsTypes) == iGas) {
362 m_rPenning[i] = 0.;
363 m_lambdaPenning[i] = 0.;
364 } else {
365 if (m_csType[i] % nCsTypes == ElectronCollisionTypeExcitation &&
366 m_rPenning[i] > Small) {
367 ++nLevelsFound;
368 }
369 }
370 }
371
372 if (nLevelsFound == 0) {
373 // There are no more excitation levels with r > 0.
374 std::cout << m_className << "::DisablePenningTransfer:\n"
375 << " Penning transfer globally switched off.\n";
376 m_usePenning = false;
377 }
378}
std::string m_gas[m_nMaxGases]
Definition: MediumGas.hh:90
bool GetGasName(const int gasnumber, const int version, std::string &gasname)
Definition: MediumGas.cc:2031

◆ DisableRadiationTrapping()

void Garfield::MediumMagboltz::DisableRadiationTrapping ( )
inline

Definition at line 54 of file MediumMagboltz.hh.

54{ m_useRadTrap = false; }

◆ EnableAnisotropicScattering()

void Garfield::MediumMagboltz::EnableAnisotropicScattering ( )
inline

Definition at line 35 of file MediumMagboltz.hh.

35 {
36 m_useAnisotropic = true;
37 m_isChanged = true;
38 }

◆ EnableCrossSectionOutput()

void Garfield::MediumMagboltz::EnableCrossSectionOutput ( )
inline

Definition at line 66 of file MediumMagboltz.hh.

66{ m_useCsOutput = true; }

◆ EnableDeexcitation()

void Garfield::MediumMagboltz::EnableDeexcitation ( )

Definition at line 174 of file MediumMagboltz.cc.

174 {
175
176 if (m_usePenning) {
177 std::cout << m_className << "::EnableDeexcitation:\n";
178 std::cout << " Penning transfer will be switched off.\n";
179 }
180 // if (m_useRadTrap) {
181 // std::cout << " Radiation trapping is switched on.\n";
182 // } else {
183 // std::cout << " Radiation trapping is switched off.\n";
184 // }
185 m_usePenning = false;
186 m_useDeexcitation = true;
187 m_isChanged = true;
188 m_dxcProducts.clear();
189}

◆ EnableEnergyRangeAdjustment()

void Garfield::MediumMagboltz::EnableEnergyRangeAdjustment ( )
inline

Definition at line 31 of file MediumMagboltz.hh.

31{ m_useAutoAdjust = true; }

◆ EnablePenningTransfer() [1/2]

void Garfield::MediumMagboltz::EnablePenningTransfer ( const double  r,
const double  lambda 
)

Definition at line 203 of file MediumMagboltz.cc.

204 {
205
206 if (r < 0. || r > 1.) {
207 std::cerr << m_className << "::EnablePenningTransfer:\n";
208 std::cerr << " Penning transfer probability must be "
209 << " in the range [0, 1].\n";
210 return;
211 }
212
214 if (lambda < Small) {
216 } else {
217 m_lambdaPenningGlobal = lambda;
218 }
219
220 std::cout << m_className << "::EnablePenningTransfer:\n";
221 std::cout << " Global Penning transfer parameters set to: \n";
222 std::cout << " r = " << m_rPenningGlobal << "\n";
223 std::cout << " lambda = " << m_lambdaPenningGlobal << " cm\n";
224
225 for (unsigned int i = 0; i < m_nTerms; ++i) {
226 m_rPenning[i] = m_rPenningGlobal;
227 m_lambdaPenning[i] = m_lambdaPenningGlobal;
228 }
229
230 if (m_useDeexcitation) {
231 std::cout << m_className << "::EnablePenningTransfer:\n";
232 std::cout << " Deexcitation handling will be switched off.\n";
233 }
234 m_usePenning = true;
235}

Referenced by GarfieldPhysics::InitializePhysics().

◆ EnablePenningTransfer() [2/2]

void Garfield::MediumMagboltz::EnablePenningTransfer ( const double  r,
const double  lambda,
std::string  gasname 
)

Definition at line 237 of file MediumMagboltz.cc.

238 {
239
240 if (r < 0. || r > 1.) {
241 std::cerr << m_className << "::EnablePenningTransfer:\n";
242 std::cerr << " Penning transfer probability must be "
243 << " in the range [0, 1].\n";
244 return;
245 }
246
247 // Get the "standard" name of this gas.
248 if (!GetGasName(gasname, gasname)) {
249 std::cerr << m_className << "::EnablePenningTransfer:\n";
250 std::cerr << " Unknown gas name.\n";
251 return;
252 }
253
254 // Look for this gas in the present gas mixture.
255 bool found = false;
256 int iGas = -1;
257 for (unsigned int i = 0; i < m_nComponents; ++i) {
258 if (m_gas[i] == gasname) {
259 m_rPenningGas[i] = r;
260 if (lambda < Small) {
261 m_lambdaPenningGas[i] = 0.;
262 } else {
263 m_lambdaPenningGas[i] = lambda;
264 }
265 found = true;
266 iGas = i;
267 break;
268 }
269 }
270
271 if (!found) {
272 std::cerr << m_className << "::EnablePenningTransfer:\n";
273 std::cerr << " Specified gas (" << gasname
274 << ") is not part of the present gas mixture.\n";
275 return;
276 }
277
278 // Make sure that the collision rate table is updated.
279 if (m_isChanged) {
280 if (!Mixer()) {
281 std::cerr << m_className << "::EnablePenningTransfer:\n";
282 std::cerr << " Error calculating the collision rates table.\n";
283 return;
284 }
285 m_isChanged = false;
286 }
287
288 unsigned int nLevelsFound = 0;
289 for (unsigned int i = 0; i < m_nTerms; ++i) {
290 if (int(m_csType[i] / nCsTypes) != iGas) continue;
291 if (m_csType[i] % nCsTypes == ElectronCollisionTypeExcitation) {
292 ++nLevelsFound;
293 }
294 m_rPenning[i] = m_rPenningGas[iGas];
295 m_lambdaPenning[i] = m_lambdaPenningGas[iGas];
296 }
297
298 if (nLevelsFound > 0) {
299 std::cout << m_className << "::EnablePenningTransfer:\n";
300 std::cout << " Penning transfer parameters for " << nLevelsFound
301 << " excitation levels set to:\n";
302 std::cout << " r = " << m_rPenningGas[iGas] << "\n";
303 std::cout << " lambda = " << m_lambdaPenningGas[iGas] << " cm\n";
304 } else {
305 std::cerr << m_className << "::EnablePenningTransfer:\n";
306 std::cerr << " Specified gas (" << gasname
307 << ") has no excitation levels in the present energy range.\n";
308 }
309
310 m_usePenning = true;
311}

◆ EnableRadiationTrapping()

void Garfield::MediumMagboltz::EnableRadiationTrapping ( )

Definition at line 191 of file MediumMagboltz.cc.

191 {
192
193 m_useRadTrap = true;
194 if (!m_useDeexcitation) {
195 std::cout << m_className << "::EnableRadiationTrapping:\n";
196 std::cout << " Radiation trapping is enabled"
197 << " but de-excitation is not.\n";
198 } else {
199 m_isChanged = true;
200 }
201}

◆ GenerateGasTable()

void Garfield::MediumMagboltz::GenerateGasTable ( const int  numCollisions = 10,
const bool  verbose = true 
)

Definition at line 5479 of file MediumMagboltz.cc.

5479 {
5480
5481 // Set the reference pressure and temperature.
5484
5485 // Initialize the parameter arrays.
5486 const unsigned int nEfields = m_eFields.size();
5487 const unsigned int nBfields = m_bFields.size();
5488 const unsigned int nAngles = m_bAngles.size();
5489 InitParamArrays(nEfields, nBfields, nAngles, tabElectronVelocityE, 0.);
5490 InitParamArrays(nEfields, nBfields, nAngles, tabElectronVelocityB, 0.);
5491 InitParamArrays(nEfields, nBfields, nAngles, tabElectronVelocityExB, 0.);
5492 InitParamArrays(nEfields, nBfields, nAngles, tabElectronDiffLong, 0.);
5493 InitParamArrays(nEfields, nBfields, nAngles, tabElectronDiffTrans, 0.);
5494 InitParamArrays(nEfields, nBfields, nAngles, tabElectronLorentzAngle, 0.);
5495 InitParamArrays(nEfields, nBfields, nAngles, tabElectronTownsend, -30.);
5496 InitParamArrays(nEfields, nBfields, nAngles, m_tabTownsendNoPenning, -30.);
5497 InitParamArrays(nEfields, nBfields, nAngles, tabElectronAttachment, -30.);
5498
5502 m_hasElectronDiffLong = true;
5506
5507 m_hasExcRates = false;
5508 m_tabExcRates.clear();
5509 m_excitationList.clear();
5510 m_hasIonRates = false;
5511 m_tabIonRates.clear();
5512 m_ionisationList.clear();
5513
5514 m_hasIonMobility = false;
5515 m_hasIonDissociation = false;
5516 m_hasIonDiffLong = false;
5517 m_hasIonDiffTrans = false;
5518
5519 // gasBits = "TFTTFTFTTTFFFFFF";
5520 // The version number is 11 because there are slight
5521 // differences between the way these gas files are written
5522 // and the ones from Garfield. This is mainly in the way
5523 // the gas tables are stored.
5524 // versionNumber = 11;
5525
5526 double vx = 0., vy = 0., vz = 0.;
5527 double difl = 0., dift = 0.;
5528 double alpha = 0., eta = 0.;
5529 double lor = 0.;
5530 double vxerr = 0., vyerr = 0., vzerr = 0.;
5531 double diflerr = 0., difterr = 0.;
5532 double alphaerr = 0., etaerr = 0.;
5533 double alphatof = 0.;
5534 double lorerr = 0.;
5535
5536 // Run through the grid of E- and B-fields and angles.
5537 for (unsigned int i = 0; i < nEfields; ++i) {
5538 for (unsigned int j = 0; j < nAngles; ++j) {
5539 for (unsigned int k = 0; k < nBfields; ++k) {
5540 if (m_debug) {
5541 std::cout << m_className << "::GenerateGasTable:\n"
5542 << " E = " << m_eFields[i] << " V/cm, B = "
5543 << m_bFields[k] << " T, angle: " << m_bAngles[j] << " rad\n";
5544 }
5545 RunMagboltz(m_eFields[i], m_bFields[k], m_bAngles[j], numColl, verbose, vx,
5546 vy, vz, difl, dift, alpha, eta, lor, vxerr, vyerr, vzerr,
5547 diflerr, difterr, alphaerr, etaerr, lorerr, alphatof);
5548 tabElectronVelocityE[j][k][i] = vz;
5549 tabElectronVelocityExB[j][k][i] = vy;
5550 tabElectronVelocityB[j][k][i] = vx;
5551 tabElectronDiffLong[j][k][i] = difl;
5552 tabElectronDiffTrans[j][k][i] = dift;
5553 tabElectronLorentzAngle[j][k][i] = lor;
5554 if (alpha > 0.) {
5555 tabElectronTownsend[j][k][i] = log(alpha);
5556 m_tabTownsendNoPenning[j][k][i] = log(alpha);
5557 } else {
5558 tabElectronTownsend[j][k][i] = -30.;
5559 m_tabTownsendNoPenning[j][k][i] = -30.;
5560 }
5561 if (eta > 0.) {
5562 tabElectronAttachment[j][k][i] = log(eta);
5563 } else {
5564 tabElectronAttachment[j][k][i] = -30.;
5565 }
5566 }
5567 }
5568 }
5569}
std::vector< excListElement > m_excitationList
Definition: MediumGas.hh:126
std::vector< std::vector< std::vector< double > > > m_tabTownsendNoPenning
Definition: MediumGas.hh:111
double m_temperatureTable
Definition: MediumGas.hh:108
std::vector< std::vector< std::vector< std::vector< double > > > > m_tabIonRates
Definition: MediumGas.hh:116
std::vector< ionListElement > m_ionisationList
Definition: MediumGas.hh:132
std::vector< std::vector< std::vector< std::vector< double > > > > m_tabExcRates
Definition: MediumGas.hh:115
void RunMagboltz(const double e, const double b, const double btheta, const int ncoll, bool verbose, double &vx, double &vy, double &vz, double &dl, double &dt, double &alpha, double &eta, double &lor, double &vxerr, double &vyerr, double &vzerr, double &dlerr, double &dterr, double &alphaerr, double &etaerr, double &lorerr, double &alphatof)
std::vector< double > m_bFields
Definition: Medium.hh:333
bool m_hasElectronDiffTrans
Definition: Medium.hh:340
std::vector< std::vector< std::vector< double > > > tabElectronVelocityB
Definition: Medium.hh:345
bool m_hasElectronVelocityB
Definition: Medium.hh:339
std::vector< std::vector< std::vector< double > > > tabElectronVelocityE
Definition: Medium.hh:343
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)
Definition: Medium.cc:2640
std::vector< std::vector< std::vector< double > > > tabElectronDiffLong
Definition: Medium.hh:346
std::vector< std::vector< std::vector< double > > > tabElectronLorentzAngle
Definition: Medium.hh:350
std::vector< std::vector< std::vector< double > > > tabElectronAttachment
Definition: Medium.hh:349
std::vector< double > m_eFields
Definition: Medium.hh:332
std::vector< double > m_bAngles
Definition: Medium.hh:334
bool m_hasElectronVelocityExB
Definition: Medium.hh:339
bool m_hasElectronLorentzAngle
Definition: Medium.hh:342
bool m_hasIonDiffTrans
Definition: Medium.hh:371
bool m_hasIonMobility
Definition: Medium.hh:370
bool m_hasElectronDiffLong
Definition: Medium.hh:340
bool m_hasIonDissociation
Definition: Medium.hh:372
bool m_hasIonDiffLong
Definition: Medium.hh:371
bool m_hasElectronVelocityE
Definition: Medium.hh:339
std::vector< std::vector< std::vector< double > > > tabElectronTownsend
Definition: Medium.hh:348
bool m_hasElectronAttachment
Definition: Medium.hh:341
std::vector< std::vector< std::vector< double > > > tabElectronDiffTrans
Definition: Medium.hh:347
std::vector< std::vector< std::vector< double > > > tabElectronVelocityExB
Definition: Medium.hh:344

◆ GetDeexcitationProduct()

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

Reimplemented from Garfield::Medium.

Definition at line 895 of file MediumMagboltz.cc.

896 {
897
898 if (i >= m_dxcProducts.size() || !(m_useDeexcitation || m_usePenning)) {
899 return false;
900 }
901 t = m_dxcProducts[i].t;
902 s = m_dxcProducts[i].s;
903 type = m_dxcProducts[i].type;
904 energy = m_dxcProducts[i].energy;
905 return true;
906}

◆ GetElectronCollision()

bool Garfield::MediumMagboltz::GetElectronCollision ( const double  e,
int &  type,
int &  level,
double &  e1,
double &  dx,
double &  dy,
double &  dz,
int &  nion,
int &  ndxc,
int &  band 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 644 of file MediumMagboltz.cc.

647 {
648
649 // Check if the electron energy is within the currently set range.
650 if (e > m_eFinal && m_useAutoAdjust) {
651 std::cerr << m_className << "::GetElectronCollision:\n";
652 std::cerr << " Provided electron energy (" << e
653 << " eV) exceeds current energy range (" << m_eFinal << " eV).\n";
654 std::cerr << " Increasing energy range to " << 1.05 * e << " eV.\n";
655 SetMaxElectronEnergy(1.05 * e);
656 } else if (e <= 0.) {
657 std::cerr << m_className << "::GetElectronCollision:\n";
658 std::cerr << " Electron energy must be greater than zero.\n";
659 return false;
660 }
661
662 // If necessary, update the collision rates table.
663 if (m_isChanged) {
664 if (!Mixer()) {
665 std::cerr << m_className << "::GetElectronCollision:\n";
666 std::cerr << " Error calculating the collision rates table.\n";
667 return false;
668 }
669 m_isChanged = false;
670 }
671
672 if (m_debug && band > 0) {
673 std::cerr << m_className << "::GetElectronCollision:\n";
674 std::cerr << " Warning: unexpected band index.\n";
675 }
676
677 double angCut = 1.;
678 double angPar = 0.5;
679
680 if (e <= m_eHigh) {
681 // Linear binning
682 // Get the energy interval.
683 int iE = int(e / m_eStep);
684 if (iE >= nEnergySteps) iE = nEnergySteps - 1;
685 if (iE < 0) iE = 0;
686
687 // Sample the scattering process.
688 const double r = RndmUniform();
689 int iLow = 0;
690 int iUp = m_nTerms - 1;
691 if (r <= m_cf[iE][iLow]) {
692 level = iLow;
693 } else if (r >= m_cf[iE][iUp]) {
694 level = iUp;
695 } else {
696 int iMid;
697 while (iUp - iLow > 1) {
698 iMid = (iLow + iUp) >> 1;
699 if (r < m_cf[iE][iMid]) {
700 iUp = iMid;
701 } else {
702 iLow = iMid;
703 }
704 }
705 level = iUp;
706 }
707 // Get the angular distribution parameters.
708 angCut = m_scatCut[iE][level];
709 angPar = m_scatParameter[iE][level];
710 } else {
711 // Logarithmic binning
712 // Get the energy interval.
713 int iE = int(log(e / m_eHigh) / m_lnStep);
714 if (iE < 0) iE = 0;
715 if (iE >= nEnergyStepsLog) iE = nEnergyStepsLog - 1;
716 // Sample the scattering process.
717 const double r = RndmUniform();
718 int iLow = 0;
719 int iUp = m_nTerms - 1;
720 if (r <= m_cfLog[iE][iLow]) {
721 level = iLow;
722 } else if (r >= m_cfLog[iE][iUp]) {
723 level = iUp;
724 } else {
725 int iMid;
726 while (iUp - iLow > 1) {
727 iMid = (iLow + iUp) >> 1;
728 if (r < m_cfLog[iE][iMid]) {
729 iUp = iMid;
730 } else {
731 iLow = iMid;
732 }
733 }
734 level = iUp;
735 }
736 // Get the angular distribution parameters.
737 angCut = m_scatCutLog[iE][level];
738 angPar = m_scatParameterLog[iE][level];
739 }
740
741 // Extract the collision type.
742 type = m_csType[level] % nCsTypes;
743 const int igas = int(m_csType[level] / nCsTypes);
744 // Increase the collision counters.
745 ++m_nCollisions[type];
746 ++m_nCollisionsDetailed[level];
747
748 // Get the energy loss for this process.
749 double loss = m_energyLoss[level];
750 nion = ndxc = 0;
751
752 if (type == ElectronCollisionTypeIonisation) {
753 // Sample the secondary electron energy according to
754 // the Opal-Beaty-Peterson parameterisation.
755 double esec = 0.;
756 if (m_useOpalBeaty) {
757 // Get the splitting parameter.
758 const double w = m_wOpalBeaty[level];
759 esec = w * tan(RndmUniform() * atan(0.5 * (e - loss) / w));
760 // Rescaling (SST)
761 // esec = w * pow(esec / w, 0.9524);
762 } else if (m_useGreenSawada) {
763 const double w = m_gsGreenSawada[igas] * e / (e + m_gbGreenSawada[igas]);
764 const double esec0 =
765 m_tsGreenSawada[igas] - m_taGreenSawada[igas] / (e + m_tbGreenSawada[igas]);
766 const double r = RndmUniform();
767 esec = esec0 + w * tan((r - 1.) * atan(esec0 / w) +
768 r * atan((0.5 * (e - loss) - esec0) / w));
769 } else {
770 esec = RndmUniform() * (e - loss);
771 }
772 if (esec <= 0) esec = Small;
773 loss += esec;
774 m_ionProducts.clear();
775 // Add the secondary electron.
776 ionProd newIonProd;
777 newIonProd.type = IonProdTypeElectron;
778 newIonProd.energy = esec;
779 m_ionProducts.push_back(newIonProd);
780 // Add the ion.
781 newIonProd.type = IonProdTypeIon;
782 newIonProd.energy = 0.;
783 m_ionProducts.push_back(newIonProd);
784 nion = 2;
785 } else if (type == ElectronCollisionTypeExcitation) {
786 // if (m_gas[igas] == "CH4" && loss * m_rgas[igas] < 13.35 && e > 12.65) {
787 // if (RndmUniform() < 0.5) {
788 // loss = 8.55 + RndmUniform() * (13.3 - 8.55);
789 // loss /= m_rgas[igas];
790 // } else {
791 // loss = std::max(Small, RndmGaussian(loss * m_rgas[igas], 1.));
792 // loss /= m_rgas[igas];
793 // }
794 // }
795 // Follow the de-excitation cascade (if switched on).
796 if (m_useDeexcitation && m_iDeexcitation[level] >= 0) {
797 int fLevel = 0;
798 ComputeDeexcitationInternal(m_iDeexcitation[level], fLevel);
799 ndxc = m_dxcProducts.size();
800 } else if (m_usePenning) {
801 m_dxcProducts.clear();
802 // Simplified treatment of Penning ionisation.
803 // If the energy threshold of this level exceeds the
804 // ionisation potential of one of the gases,
805 // create a new electron (with probability m_rPenning).
806 if (m_energyLoss[level] * m_rgas[igas] > m_minIonPot &&
807 RndmUniform() < m_rPenning[level]) {
808 // The energy of the secondary electron is assumed to be given by
809 // the difference of excitation and ionisation threshold.
810 double esec = m_energyLoss[level] * m_rgas[igas] - m_minIonPot;
811 if (esec <= 0) esec = Small;
812 // Add the secondary electron to the list.
813 dxcProd newDxcProd;
814 newDxcProd.t = 0.;
815 newDxcProd.s = 0.;
816 if (m_lambdaPenning[level] > Small) {
817 // Uniform distribution within a sphere of radius lambda
818 newDxcProd.s = m_lambdaPenning[level] * pow(RndmUniformPos(), 1. / 3.);
819 }
820 newDxcProd.energy = esec;
821 newDxcProd.type = DxcProdTypeElectron;
822 m_dxcProducts.push_back(newDxcProd);
823 ndxc = 1;
824 ++m_nPenning;
825 }
826 }
827 }
828
829 // Make sure the energy loss is smaller than the energy.
830 if (e < loss) loss = e - 0.0001;
831
832 // Determine the scattering angle.
833 double ctheta0 = 1. - 2. * RndmUniform();
834 if (m_useAnisotropic) {
835 switch (m_scatModel[level]) {
836 case 0:
837 break;
838 case 1:
839 ctheta0 = 1. - RndmUniform() * angCut;
840 if (RndmUniform() > angPar) ctheta0 = -ctheta0;
841 break;
842 case 2:
843 ctheta0 = (ctheta0 + angPar) / (1. + angPar * ctheta0);
844 break;
845 default:
846 std::cerr << m_className << "::GetElectronCollision:\n";
847 std::cerr << " Unknown scattering model. \n";
848 std::cerr << " Using isotropic distribution.\n";
849 break;
850 }
851 }
852
853 const double s1 = m_rgas[igas];
854 const double s2 = (s1 * s1) / (s1 - 1.);
855 const double theta0 = acos(ctheta0);
856 const double arg = std::max(1. - s1 * loss / e, Small);
857 const double d = 1. - ctheta0 * sqrt(arg);
858
859 // Update the energy.
860 e1 = std::max(e * (1. - loss / (s1 * e) - 2. * d / s2), Small);
861 double q = std::min(sqrt((e / e1) * arg) / s1, 1.);
862 const double theta = asin(q * sin(theta0));
863 double ctheta = cos(theta);
864 if (ctheta0 < 0.) {
865 const double u = (s1 - 1.) * (s1 - 1.) / arg;
866 if (ctheta0 * ctheta0 > u) ctheta = -ctheta;
867 }
868 const double stheta = sin(theta);
869 // Calculate the direction after the collision.
870 dz = std::min(dz, 1.);
871 const double argZ = sqrt(dx * dx + dy * dy);
872
873 // Azimuth is chosen at random.
874 const double phi = TwoPi * RndmUniform();
875 const double cphi = cos(phi);
876 const double sphi = sin(phi);
877
878 if (argZ == 0.) {
879 dz = ctheta;
880 dx = cphi * stheta;
881 dy = sphi * stheta;
882 } else {
883 const double a = stheta / argZ;
884 const double dz1 = dz * ctheta + argZ * stheta * sphi;
885 const double dy1 = dy * ctheta + a * (dx * cphi - dy * dz * sphi);
886 const double dx1 = dx * ctheta - a * (dy * cphi + dx * dz * sphi);
887 dz = dz1;
888 dy = dy1;
889 dx = dx1;
890 }
891
892 return true;
893}
bool SetMaxElectronEnergy(const double e)
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition: Random.hh:17
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337
DoubleAc acos(const DoubleAc &f)
Definition: DoubleAc.cpp:490
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384
DoubleAc asin(const DoubleAc &f)
Definition: DoubleAc.cpp:470
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

◆ GetElectronCollisionRate() [1/2]

double Garfield::MediumMagboltz::GetElectronCollisionRate ( const double  e,
const int  band 
)
virtual

Reimplemented from Garfield::Medium.

Definition at line 546 of file MediumMagboltz.cc.

547 {
548
549 // Check if the electron energy is within the currently set range.
550 if (e <= 0.) {
551 std::cerr << m_className << "::GetElectronCollisionRate:\n";
552 std::cerr << " Electron energy must be greater than zero.\n";
553 return m_cfTot[0];
554 }
555 if (e > m_eFinal && m_useAutoAdjust) {
556 std::cerr << m_className << "::GetElectronCollisionRate:\n";
557 std::cerr << " Collision rate at " << e
558 << " eV is not included in the current table.\n";
559 std::cerr << " Increasing energy range to " << 1.05 * e << " eV.\n";
560 SetMaxElectronEnergy(1.05 * e);
561 }
562
563 // If necessary, update the collision rates table.
564 if (m_isChanged) {
565 if (!Mixer()) {
566 std::cerr << m_className << "::GetElectronCollisionRate:\n";
567 std::cerr << " Error calculating the collision rates table.\n";
568 return 0.;
569 }
570 m_isChanged = false;
571 }
572
573 if (m_debug && band > 0) {
574 std::cerr << m_className << "::GetElectronCollisionRate:\n";
575 std::cerr << " Warning: unexpected band index.\n";
576 }
577
578 // Get the energy interval.
579 int iE = 0;
580 if (e <= m_eHigh) {
581 // Linear binning
582 iE = int(e / m_eStep);
583 if (iE >= nEnergySteps) return m_cfTot[nEnergySteps - 1];
584 if (iE < 0) return m_cfTot[0];
585 return m_cfTot[iE];
586 }
587
588 // Logarithmic binning
589 const double eLog = log(e);
590 iE = int((eLog - m_eHighLog) / m_lnStep);
591 // Calculate the collision rate by log-log interpolation.
592 const double fmax = m_cfTotLog[iE];
593 const double fmin = iE == 0 ? log(m_cfTot[nEnergySteps - 1]) : m_cfTotLog[iE - 1];
594 const double emin = m_eHighLog + iE * m_lnStep;
595 const double f = fmin + (eLog - emin) * (fmax - fmin) / m_lnStep;
596 return exp(f);
597}
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377

Referenced by GetElectronCollisionRate().

◆ GetElectronCollisionRate() [2/2]

double Garfield::MediumMagboltz::GetElectronCollisionRate ( const double  e,
const unsigned int  level,
const int  band 
)

Definition at line 599 of file MediumMagboltz.cc.

601 {
602
603 // Check if the electron energy is within the currently set range.
604 if (e <= 0.) {
605 std::cerr << m_className << "::GetElectronCollisionRate:\n";
606 std::cerr << " Electron energy must be greater than zero.\n";
607 return 0.;
608 }
609
610 // Check if the level exists.
611 if (level >= m_nTerms) {
612 std::cerr << m_className << "::GetElectronCollisionRate:\n";
613 std::cerr << " Level " << level << " does not exist.\n";
614 std::cerr << " The present gas mixture has " << m_nTerms
615 << " cross-section terms.\n";
616 return 0.;
617 }
618
619 // Get the total scattering rate.
620 double rate = GetElectronCollisionRate(e, band);
621 // Get the energy interval.
622 int iE = 0;
623 if (e <= m_eHigh) {
624 // Linear binning
625 iE = int(e / m_eStep);
626 if (iE >= nEnergySteps) return m_cfTot[nEnergySteps - 1];
627 if (level == 0) {
628 rate *= m_cf[iE][0];
629 } else {
630 rate *= m_cf[iE][level] - m_cf[iE][level - 1];
631 }
632 } else {
633 // Logarithmic binning
634 iE = int((log(e) - m_eHighLog) / m_lnStep);
635 if (level == 0) {
636 rate *= m_cfLog[iE][0];
637 } else {
638 rate *= m_cfLog[iE][level] - m_cfLog[iE][level - 1];
639 }
640 }
641 return rate;
642}
double GetElectronCollisionRate(const double e, const int band)

◆ GetElectronNullCollisionRate()

double Garfield::MediumMagboltz::GetElectronNullCollisionRate ( const int  band)
virtual

Reimplemented from Garfield::Medium.

Definition at line 526 of file MediumMagboltz.cc.

526 {
527
528 // If necessary, update the collision rates table.
529 if (m_isChanged) {
530 if (!Mixer()) {
531 std::cerr << m_className << "::GetElectronNullCollisionRate:\n";
532 std::cerr << " Error calculating the collision rates table.\n";
533 return 0.;
534 }
535 m_isChanged = false;
536 }
537
538 if (m_debug && band > 0) {
539 std::cerr << m_className << "::GetElectronNullCollisionRate:\n";
540 std::cerr << " Warning: unexpected band index.\n";
541 }
542
543 return m_cfNull;
544}

◆ GetIonisationProduct()

bool Garfield::MediumMagboltz::GetIonisationProduct ( const unsigned int  i,
int &  type,
double &  energy 
) const
virtual

Reimplemented from Garfield::Medium.

Definition at line 908 of file MediumMagboltz.cc.

909 {
910
911 if (i >= m_ionProducts.size()) {
912 std::cerr << m_className << "::GetIonisationProduct:\n"
913 << " Index (" << i << ") out of range.\n";
914 return false;
915 }
916
917 type = m_ionProducts[i].type;
918 energy = m_ionProducts[i].energy;
919 return true;
920}

◆ GetLevel()

bool Garfield::MediumMagboltz::GetLevel ( const unsigned int  i,
int &  ngas,
int &  type,
std::string &  descr,
double &  e 
)

Definition at line 1122 of file MediumMagboltz.cc.

1123 {
1124
1125 if (m_isChanged) {
1126 if (!Mixer()) {
1127 std::cerr << m_className << "::GetLevel:\n";
1128 std::cerr << " Error calculating the collision rates table.\n";
1129 return false;
1130 }
1131 m_isChanged = false;
1132 }
1133
1134 if (i >= m_nTerms) {
1135 std::cerr << m_className << "::GetLevel:\n";
1136 std::cerr << " Requested level (" << i << ") does not exist.\n";
1137 return false;
1138 }
1139
1140 // Collision type
1141 type = m_csType[i] % nCsTypes;
1142 ngas = int(m_csType[i] / nCsTypes);
1143 // Description (from Magboltz)
1144 descr = std::string(50, ' ');
1145 for (int j = 50; j--;) descr[j] = m_description[i][j];
1146 // Threshold energy
1147 e = m_rgas[ngas] * m_energyLoss[i];
1148 if (m_debug) {
1149 std::cout << m_className << "::GetLevel:\n";
1150 std::cout << " Level " << i << ": " << descr << "\n";
1151 std::cout << " Type " << type << "\n",
1152 std::cout << " Threshold energy: " << e << " eV\n";
1153 if (type == ElectronCollisionTypeExcitation && m_usePenning &&
1154 e > m_minIonPot) {
1155 std::cout << " Penning transfer coefficient: " << m_rPenning[i] << "\n";
1156 } else if (type == ElectronCollisionTypeExcitation && m_useDeexcitation) {
1157 const int idxc = m_iDeexcitation[i];
1158 if (idxc < 0 || idxc >= (int)m_deexcitations.size()) {
1159 std::cout << " Deexcitation cascade not implemented.\n";
1160 return true;
1161 }
1162 if (m_deexcitations[idxc].osc > 0.) {
1163 std::cout << " Oscillator strength: " << m_deexcitations[idxc].osc
1164 << "\n";
1165 }
1166 std::cout << " Decay channels:\n";
1167 for (int j = 0; j < m_deexcitations[idxc].nChannels; ++j) {
1168 if (m_deexcitations[idxc].type[j] == DxcTypeRad) {
1169 std::cout << " Radiative decay to ";
1170 if (m_deexcitations[idxc].final[j] < 0) {
1171 std::cout << "ground state: ";
1172 } else {
1173 std::cout << m_deexcitations[m_deexcitations[idxc].final[j]].label
1174 << ": ";
1175 }
1176 } else if (m_deexcitations[idxc].type[j] == DxcTypeCollIon) {
1177 if (m_deexcitations[idxc].final[j] < 0) {
1178 std::cout << " Penning ionisation: ";
1179 } else {
1180 std::cout << " Associative ionisation: ";
1181 }
1182 } else if (m_deexcitations[idxc].type[j] == DxcTypeCollNonIon) {
1183 if (m_deexcitations[idxc].final[j] >= 0) {
1184 std::cout << " Collision-induced transition to "
1185 << m_deexcitations[m_deexcitations[idxc].final[j]].label
1186 << ": ";
1187 } else {
1188 std::cout << " Loss: ";
1189 }
1190 }
1191 if (j == 0) {
1192 std::cout << std::setprecision(5) << m_deexcitations[idxc].p[j] * 100.
1193 << "%\n";
1194 } else {
1195 std::cout << std::setprecision(5) << (m_deexcitations[idxc].p[j] -
1196 m_deexcitations[idxc].p[j - 1]) *
1197 100. << "%\n";
1198 }
1199 }
1200 }
1201 }
1202
1203 return true;
1204}

◆ GetMaxElectronEnergy()

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

Definition at line 22 of file MediumMagboltz.hh.

22{ return m_eFinal; }

◆ GetMaxPhotonEnergy()

double Garfield::MediumMagboltz::GetMaxPhotonEnergy ( ) const
inline

Definition at line 27 of file MediumMagboltz.hh.

27{ return m_eFinalGamma; }

◆ GetNumberOfDeexcitationProducts()

unsigned int Garfield::MediumMagboltz::GetNumberOfDeexcitationProducts ( ) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 92 of file MediumMagboltz.hh.

92 {
93 return m_dxcProducts.size();
94 }

◆ GetNumberOfElectronCollisions() [1/3]

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

Definition at line 1087 of file MediumMagboltz.cc.

1087 {
1088
1089 unsigned int ncoll = 0;
1090 for (int j = nCsTypes; j--;) ncoll += m_nCollisions[j];
1091 return ncoll;
1092}

◆ GetNumberOfElectronCollisions() [2/3]

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

Definition at line 1206 of file MediumMagboltz.cc.

1206 {
1207
1208 if (level >= m_nTerms) {
1209 std::cerr << m_className << "::GetNumberOfElectronCollisions:\n"
1210 << " Cross-section term (" << level << ") does not exist.\n";
1211 return 0;
1212 }
1213 return m_nCollisionsDetailed[level];
1214}

◆ GetNumberOfElectronCollisions() [3/3]

unsigned int Garfield::MediumMagboltz::GetNumberOfElectronCollisions ( int &  nElastic,
int &  nIonising,
int &  nAttachment,
int &  nInelastic,
int &  nExcitation,
int &  nSuperelastic 
) const

Definition at line 1094 of file MediumMagboltz.cc.

1096 {
1097
1098 nElastic = m_nCollisions[ElectronCollisionTypeElastic];
1099 nIonisation = m_nCollisions[ElectronCollisionTypeIonisation];
1100 nAttachment = m_nCollisions[ElectronCollisionTypeAttachment];
1101 nInelastic = m_nCollisions[ElectronCollisionTypeInelastic];
1102 nExcitation = m_nCollisions[ElectronCollisionTypeExcitation];
1103 nSuperelastic = m_nCollisions[ElectronCollisionTypeSuperelastic];
1104 return nElastic + nIonisation + nAttachment + nInelastic + nExcitation +
1105 nSuperelastic;
1106}

◆ GetNumberOfIonisationProducts()

unsigned int Garfield::MediumMagboltz::GetNumberOfIonisationProducts ( ) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 86 of file MediumMagboltz.hh.

86 {
87 return m_ionProducts.size();
88 }

◆ GetNumberOfLevels()

int Garfield::MediumMagboltz::GetNumberOfLevels ( )

Definition at line 1108 of file MediumMagboltz.cc.

1108 {
1109
1110 if (m_isChanged) {
1111 if (!Mixer()) {
1112 std::cerr << m_className << "::GetNumberOfLevels:\n";
1113 std::cerr << " Error calculating the collision rates table.\n";
1114 return 0;
1115 }
1116 m_isChanged = false;
1117 }
1118
1119 return m_nTerms;
1120}

◆ GetNumberOfPenningTransfers()

int Garfield::MediumMagboltz::GetNumberOfPenningTransfers ( ) const
inline

Definition at line 119 of file MediumMagboltz.hh.

119{ return m_nPenning; }

◆ GetNumberOfPhotonCollisions() [1/2]

int Garfield::MediumMagboltz::GetNumberOfPhotonCollisions ( ) const

Definition at line 1216 of file MediumMagboltz.cc.

1216 {
1217
1218 int ncoll = 0;
1219 for (int j = nCsTypesGamma; j--;) ncoll += m_nPhotonCollisions[j];
1220 return ncoll;
1221}

◆ GetNumberOfPhotonCollisions() [2/2]

int Garfield::MediumMagboltz::GetNumberOfPhotonCollisions ( int &  nElastic,
int &  nIonising,
int &  nInelastic 
) const

Definition at line 1223 of file MediumMagboltz.cc.

1224 {
1225
1226 nElastic = m_nPhotonCollisions[0];
1227 nIonising = m_nPhotonCollisions[1];
1228 nInelastic = m_nPhotonCollisions[2];
1229 return nElastic + nIonising + nInelastic;
1230}

◆ GetPhotonCollision()

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

Reimplemented from Garfield::Medium.

Definition at line 968 of file MediumMagboltz.cc.

970 {
971
972 if (e > m_eFinalGamma && m_useAutoAdjust) {
973 std::cerr << m_className << "::GetPhotonCollision:\n";
974 std::cerr << " Provided electron energy (" << e
975 << " eV) exceeds current energy range (" << m_eFinalGamma
976 << " eV).\n";
977 std::cerr << " Increasing energy range to " << 1.05 * e << " eV.\n";
978 SetMaxPhotonEnergy(1.05 * e);
979 } else if (e <= 0.) {
980 std::cerr << m_className << "::GetPhotonCollision:\n";
981 std::cerr << " Photon energy must be greater than zero.\n";
982 return false;
983 }
984
985 if (m_isChanged) {
986 if (!Mixer()) {
987 std::cerr << m_className << "::GetPhotonCollision:\n";
988 std::cerr << " Error calculating the collision rates table.\n";
989 return false;
990 }
991 m_isChanged = false;
992 }
993
994 // Energy interval
995 int iE = int(e / m_eStepGamma);
996 if (iE >= nEnergyStepsGamma) iE = nEnergyStepsGamma - 1;
997 if (iE < 0) iE = 0;
998
999 double r = m_cfTotGamma[iE];
1000 if (m_useDeexcitation && m_useRadTrap && !m_deexcitations.empty()) {
1001 int nLines = 0;
1002 std::vector<double> pLine(0);
1003 std::vector<int> iLine(0);
1004 // Loop over the excitations.
1005 const unsigned int nDeexcitations = m_deexcitations.size();
1006 for (unsigned int i = 0; i < nDeexcitations; ++i) {
1007 if (m_deexcitations[i].cf > 0. &&
1008 fabs(e - m_deexcitations[i].energy) <= m_deexcitations[i].width) {
1009 r += m_deexcitations[i].cf * TMath::Voigt(e - m_deexcitations[i].energy,
1010 m_deexcitations[i].sDoppler,
1011 2 * m_deexcitations[i].gPressure);
1012 pLine.push_back(r);
1013 iLine.push_back(i);
1014 ++nLines;
1015 }
1016 }
1017 r *= RndmUniform();
1018 if (nLines > 0 && r >= m_cfTotGamma[iE]) {
1019 // Photon is absorbed by a discrete line.
1020 for (int i = 0; i < nLines; ++i) {
1021 if (r <= pLine[i]) {
1022 ++m_nPhotonCollisions[PhotonCollisionTypeExcitation];
1023 int fLevel = 0;
1024 ComputeDeexcitationInternal(iLine[i], fLevel);
1025 type = PhotonCollisionTypeExcitation;
1026 nsec = nDeexcitationProducts;
1027 return true;
1028 }
1029 }
1030 std::cerr << m_className << "::GetPhotonCollision:\n";
1031 std::cerr << " Random sampling of deexcitation line failed.\n";
1032 std::cerr << " Program bug!\n";
1033 return false;
1034 }
1035 } else {
1036 r *= RndmUniform();
1037 }
1038
1039 int iLow = 0;
1040 int iUp = nPhotonTerms - 1;
1041 if (r <= m_cfGamma[iE][iLow]) {
1042 level = iLow;
1043 } else if (r >= m_cfGamma[iE][iUp]) {
1044 level = iUp;
1045 } else {
1046 int iMid;
1047 while (iUp - iLow > 1) {
1048 iMid = (iLow + iUp) >> 1;
1049 if (r < m_cfGamma[iE][iMid]) {
1050 iUp = iMid;
1051 } else {
1052 iLow = iMid;
1053 }
1054 }
1055 level = iUp;
1056 }
1057
1058 nsec = 0;
1059 esec = e1 = 0.;
1060 type = csTypeGamma[level];
1061 // Collision type
1062 type = type % nCsTypesGamma;
1063 int ngas = int(csTypeGamma[level] / nCsTypesGamma);
1064 ++m_nPhotonCollisions[type];
1065 // Ionising collision
1066 if (type == 1) {
1067 esec = e - m_ionPot[ngas];
1068 if (esec < Small) esec = Small;
1069 nsec = 1;
1070 }
1071
1072 // Determine the scattering angle
1073 ctheta = 2 * RndmUniform() - 1.;
1074
1075 return true;
1076}
bool SetMaxPhotonEnergy(const double e)
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

◆ GetPhotonCollisionRate()

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

Reimplemented from Garfield::Medium.

Definition at line 922 of file MediumMagboltz.cc.

922 {
923
924 if (e <= 0.) {
925 std::cerr << m_className << "::GetPhotonCollisionRate:\n";
926 std::cerr << " Photon energy must be greater than zero.\n";
927 return m_cfTotGamma[0];
928 }
929 if (e > m_eFinalGamma && m_useAutoAdjust) {
930 std::cerr << m_className << "::GetPhotonCollisionRate:\n";
931 std::cerr << " Collision rate at " << e
932 << " eV is not included in the current table.\n";
933 std::cerr << " Increasing energy range to " << 1.05 * e << " eV.\n";
934 SetMaxPhotonEnergy(1.05 * e);
935 }
936
937 if (m_isChanged) {
938 if (!Mixer()) {
939 std::cerr << m_className << "::GetPhotonCollisionRate:\n";
940 std::cerr << " Error calculating the collision rates table.\n";
941 return 0.;
942 }
943 m_isChanged = false;
944 }
945
946 int iE = int(e / m_eStepGamma);
947 if (iE >= nEnergyStepsGamma) iE = nEnergyStepsGamma - 1;
948 if (iE < 0) iE = 0;
949
950 double cfSum = m_cfTotGamma[iE];
951 if (m_useDeexcitation && m_useRadTrap && !m_deexcitations.empty()) {
952 // Loop over the excitations.
953 const unsigned int nDeexcitations = m_deexcitations.size();
954 for (unsigned int i = 0; i < nDeexcitations; ++i) {
955 if (m_deexcitations[i].cf > 0. &&
956 fabs(e - m_deexcitations[i].energy) <= m_deexcitations[i].width) {
957 cfSum +=
958 m_deexcitations[i].cf * TMath::Voigt(e - m_deexcitations[i].energy,
959 m_deexcitations[i].sDoppler,
960 2 * m_deexcitations[i].gPressure);
961 }
962 }
963 }
964
965 return cfSum;
966}

◆ Initialise()

bool Garfield::MediumMagboltz::Initialise ( const bool  verbose = false)

Definition at line 417 of file MediumMagboltz.cc.

417 {
418
419 if (!m_isChanged) {
420 if (m_debug) {
421 std::cerr << m_className << "::Initialise:\n";
422 std::cerr << " Nothing changed.\n";
423 }
424 return true;
425 }
426 if (!Mixer(verbose)) {
427 std::cerr << m_className << "::Initialise:\n";
428 std::cerr << " Error calculating the collision rates table.\n";
429 return false;
430 }
431 m_isChanged = false;
432 return true;
433}

Referenced by GarfieldPhysics::InitializePhysics(), and PrintGas().

◆ PrintGas()

void Garfield::MediumMagboltz::PrintGas ( )

Definition at line 435 of file MediumMagboltz.cc.

435 {
436
438
439 if (m_isChanged) {
440 if (!Initialise()) return;
441 }
442
443 std::cout << m_className << "::PrintGas:\n";
444 for (unsigned int i = 0; i < m_nTerms; ++i) {
445 // Collision type
446 int type = m_csType[i] % nCsTypes;
447 int ngas = int(m_csType[i] / nCsTypes);
448 // Description (from Magboltz)
449 std::string descr = std::string(50, ' ');
450 for (int j = 50; j--;) descr[j] = m_description[i][j];
451 // Threshold energy
452 double e = m_rgas[ngas] * m_energyLoss[i];
453 std::cout << " Level " << i << ": " << descr << "\n";
454 std::cout << " Type " << type;
455 if (type == ElectronCollisionTypeElastic) {
456 std::cout << " (elastic)\n";
457 } else if (type == ElectronCollisionTypeIonisation) {
458 std::cout << " (ionisation)\n";
459 std::cout << " Ionisation threshold: " << e << " eV\n";
460 } else if (type == ElectronCollisionTypeAttachment) {
461 std::cout << " (attachment)\n";
462 } else if (type == ElectronCollisionTypeInelastic) {
463 std::cout << " (inelastic)\n";
464 std::cout << " Energy loss: " << e << " eV\n";
465 } else if (type == ElectronCollisionTypeExcitation) {
466 std::cout << " (excitation)\n";
467 std::cout << " Excitation energy: " << e << " eV\n";
468 } else if (type == ElectronCollisionTypeSuperelastic) {
469 std::cout << " (super-elastic)\n";
470 std::cout << " Energy gain: " << -e << " eV\n";
471 } else {
472 std::cout << " (unknown)\n";
473 }
474 if (type == ElectronCollisionTypeExcitation && m_usePenning &&
475 e > m_minIonPot) {
476 std::cout << " Penning transfer coefficient: " << m_rPenning[i]
477 << "\n";
478 } else if (type == ElectronCollisionTypeExcitation && m_useDeexcitation) {
479 const int idxc = m_iDeexcitation[i];
480 if (idxc < 0 || idxc >= (int)m_deexcitations.size()) {
481 std::cout << " Deexcitation cascade not implemented.\n";
482 continue;
483 }
484 if (m_deexcitations[idxc].osc > 0.) {
485 std::cout << " Oscillator strength: " << m_deexcitations[idxc].osc
486 << "\n";
487 }
488 std::cout << " Decay channels:\n";
489 for (int j = 0; j < m_deexcitations[idxc].nChannels; ++j) {
490 if (m_deexcitations[idxc].type[j] == DxcTypeRad) {
491 std::cout << " Radiative decay to ";
492 if (m_deexcitations[idxc].final[j] < 0) {
493 std::cout << "ground state: ";
494 } else {
495 std::cout << m_deexcitations[m_deexcitations[idxc].final[j]].label
496 << ": ";
497 }
498 } else if (m_deexcitations[idxc].type[j] == DxcTypeCollIon) {
499 if (m_deexcitations[idxc].final[j] < 0) {
500 std::cout << " Penning ionisation: ";
501 } else {
502 std::cout << " Associative ionisation: ";
503 }
504 } else if (m_deexcitations[idxc].type[j] == DxcTypeCollNonIon) {
505 if (m_deexcitations[idxc].final[j] >= 0) {
506 std::cout << " Collision-induced transition to "
507 << m_deexcitations[m_deexcitations[idxc].final[j]].label
508 << ": ";
509 } else {
510 std::cout << " Loss: ";
511 }
512 }
513 if (j == 0) {
514 std::cout << std::setprecision(5) << m_deexcitations[idxc].p[j] * 100.
515 << "%\n";
516 } else {
517 std::cout << std::setprecision(5) << (m_deexcitations[idxc].p[j] -
518 m_deexcitations[idxc].p[j - 1]) *
519 100. << "%\n";
520 }
521 }
522 }
523 }
524}
bool Initialise(const bool verbose=false)

◆ ResetCollisionCounters()

void Garfield::MediumMagboltz::ResetCollisionCounters ( )

Definition at line 1078 of file MediumMagboltz.cc.

1078 {
1079
1080 for (int j = nCsTypes; j--;) m_nCollisions[j] = 0;
1081 m_nCollisionsDetailed.resize(m_nTerms);
1082 for (unsigned int j = 0; j < m_nTerms; ++j) m_nCollisionsDetailed[j] = 0;
1083 m_nPenning = 0;
1084 for (int j = nCsTypesGamma; j--;) m_nPhotonCollisions[j] = 0;
1085}

◆ RunMagboltz()

void Garfield::MediumMagboltz::RunMagboltz ( const double  e,
const double  b,
const double  btheta,
const int  ncoll,
bool  verbose,
double &  vx,
double &  vy,
double &  vz,
double &  dl,
double &  dt,
double &  alpha,
double &  eta,
double &  lor,
double &  vxerr,
double &  vyerr,
double &  vzerr,
double &  dlerr,
double &  dterr,
double &  alphaerr,
double &  etaerr,
double &  lorerr,
double &  alphatof 
)

Definition at line 5278 of file MediumMagboltz.cc.

5286 {
5287
5288 // Initialize the values.
5289 vx = vy = vz = 0.;
5290 dl = dt = 0.;
5291 alpha = eta = alphatof = 0.;
5292 lor = 0.;
5293 vxerr = vyerr = vzerr = 0.;
5294 dlerr = dterr = 0.;
5295 alphaerr = etaerr = 0.;
5296 lorerr = 0.;
5297
5298 // Set input parameters in Magboltz common blocks.
5300 Magboltz::inpt_.nStep = 4000;
5301 Magboltz::inpt_.nAniso = 2;
5302
5303 Magboltz::inpt_.tempc = m_temperature - ZeroCelsius;
5305 Magboltz::inpt_.ipen = 0;
5306 Magboltz::setp_.nmax = ncoll;
5307
5308 Magboltz::setp_.efield = e;
5309 // Convert from Tesla to kGauss.
5310 Magboltz::bfld_.bmag = bmag * 10.;
5311 // Convert from radians to degree.
5312 Magboltz::bfld_.btheta = btheta * 180. / Pi;
5313
5314 // Set the gas composition in Magboltz.
5315 for (unsigned int i = 0; i < m_nComponents; ++i) {
5316 int ng = 0;
5317 if (!GetGasNumberMagboltz(m_gas[i], ng)) {
5318 std::cerr << m_className << "::RunMagboltz:\n"
5319 << " Gas " << m_gas[i] << " has no corresponding"
5320 << " gas number in Magboltz.\n";
5321 return;
5322 }
5323 Magboltz::gasn_.ngasn[i] = ng;
5324 Magboltz::ratio_.frac[i] = 100 * m_fraction[i];
5325 }
5326
5327 // Call Magboltz internal setup routine.
5329
5330 // Calculate the max. energy in the table.
5331 if (e * m_temperature / (293.15 * m_pressure) > 15) {
5332 // If E/p > 15 start with 8 eV.
5333 Magboltz::inpt_.efinal = 8.;
5334 } else {
5335 Magboltz::inpt_.efinal = 0.5;
5336 }
5337 Magboltz::setp_.estart = Magboltz::inpt_.efinal / 50.;
5338
5339 long long ielow = 1;
5340 while (ielow == 1) {
5342 if (bmag == 0. || btheta == 0. || fabs(btheta) == Pi) {
5343 Magboltz::elimit_(&ielow);
5344 } else if (btheta == HalfPi) {
5345 Magboltz::elimitb_(&ielow);
5346 } else {
5347 Magboltz::elimitc_(&ielow);
5348 }
5349 if (ielow == 1) {
5350 // Increase the max. energy.
5351 Magboltz::inpt_.efinal *= sqrt(2.);
5352 Magboltz::setp_.estart = Magboltz::inpt_.efinal / 50.;
5353 }
5354 }
5355
5356 if (m_debug || verbose) Magboltz::prnter_();
5357
5358 // Run the Monte Carlo calculation.
5359 if (bmag == 0.) {
5361 } else if (btheta == 0. || btheta == Pi) {
5363 } else if (btheta == HalfPi) {
5365 } else {
5367 }
5368 if (m_debug || verbose) Magboltz::output_();
5369
5370 // If attachment or ionisation rate is greater than sstmin,
5371 // include spatial gradients in the solution.
5372 const double sstmin = 30.;
5373 const double epscale = 760. * m_temperature / (m_pressure * 293.15);
5374 double alpp = Magboltz::ctowns_.alpha * epscale;
5375 double attp = Magboltz::ctowns_.att * epscale;
5376 bool useSST = false;
5377 if (fabs(alpp - attp) > sstmin || alpp > sstmin || attp > sstmin) {
5378 useSST = true;
5379 if (bmag == 0.) {
5381 } else if (btheta == 0. || btheta == Pi) {
5383 } else if (btheta == HalfPi) {
5385 } else {
5387 }
5388 // Calculate the (effective) TOF Townsend coefficient.
5389 double alphapt = Magboltz::tofout_.ralpha;
5390 double etapt = Magboltz::tofout_.rattof;
5391 double fc1 =
5392 1.e5 * Magboltz::tofout_.tofwr / (2. * Magboltz::tofout_.tofdl);
5393 double fc2 = 1.e12 * (alphapt - etapt) / Magboltz::tofout_.tofdl;
5394 alphatof = fc1 - sqrt(fc1 * fc1 - fc2);
5395 }
5396 if (m_debug || verbose) Magboltz::output2_();
5397
5398 // Velocities. Convert to cm / ns.
5399 vx = Magboltz::vel_.wx * 1.e-9;
5400 vxerr = Magboltz::velerr_.dwx;
5401 vy = Magboltz::vel_.wy * 1.e-9;
5402 vyerr = Magboltz::velerr_.dwy;
5403 vz = Magboltz::vel_.wz * 1.e-9;
5404 vzerr = Magboltz::velerr_.dwz;
5405
5406 // Calculate the Lorentz angle.
5407 const double forcalc = vx * vx + vy * vy;
5408 double elvel = sqrt(forcalc + vz * vz);
5409 if (forcalc != 0 && elvel != 0) {
5410 lor = acos(vz / elvel);
5411 const double ainlorerr = sqrt(forcalc * forcalc * vzerr * vzerr +
5412 vx * vx * vx * vx * vxerr * vxerr +
5413 vy * vy * vy * vy * vyerr * vyerr);
5414 lorerr = vz * ainlorerr/ elvel / elvel / sqrt (forcalc) / lor;
5415 }
5416
5417 // Diffusion coefficients.
5418 // dt = sqrt(0.2 * Magboltz::difvel_.diftr / vz) * 1.e-4;
5419 dt = sqrt(0.2 * 0.5 * (Magboltz::diflab_.difxx + Magboltz::diflab_.difyy) /
5420 vz) * 1.e-4;
5421 dterr = Magboltz::diferl_.dfter;
5422 // dl = sqrt(0.2 * Magboltz::difvel_.difln / vz) * 1.e-4;
5423 dl = sqrt(0.2 * Magboltz::diflab_.difzz / vz) * 1.e-4;
5424 dlerr = Magboltz::diferl_.dfler;
5425 // Diffusion tensor.
5426 // SBOL(1)=2D-6*DIFZZ*TORR/VBOL
5427 // SBOL(2)=2D-6*DIFXX*TORR/VBOL
5428 // SBOL(3)=2D-6*DIFYY*TORR/VBOL
5429 // SBOL(4)=2D-6*DIFXZ*TORR/VBOL
5430 // SBOL(5)=2D-6*DIFYZ*TORR/VBOL
5431 // SBOL(6)=2D-6*DIFXY*TORR/VBOL
5432 alpha = Magboltz::ctowns_.alpha;
5433 alphaerr = Magboltz::ctwner_.alper;
5434 eta = Magboltz::ctowns_.att;
5435 etaerr = Magboltz::ctwner_.atter;
5436
5437 // Print the results.
5438 if (m_debug) {
5439 std::cout << m_className << "::RunMagboltz:\n Results:\n";
5440 std::cout << " Drift velocity along E: " << std::right
5441 << std::setw(10) << std::setprecision(6) << vz << " cm/ns +/- "
5442 << std::setprecision(2) << vzerr << "%\n";
5443 std::cout << " Drift velocity along Bt: " << std::right
5444 << std::setw(10) << std::setprecision(6) << vx << " cm/ns +/- "
5445 << std::setprecision(2) << vxerr << "%\n";
5446 std::cout << " Drift velocity along ExB: " << std::right
5447 << std::setw(10) << std::setprecision(6) << vy << " cm/ns +/- "
5448 << std::setprecision(2) << vyerr << "%\n";
5449 std::cout << " Longitudinal diffusion: " << std::right
5450 << std::setw(10) << std::setprecision(6) << dl << " cm1/2 +/- "
5451 << std::setprecision(2) << dlerr << "%\n";
5452 std::cout << " Transverse diffusion: " << std::right
5453 << std::setw(10) << std::setprecision(6) << dt << " cm1/2 +/- "
5454 << std::setprecision(2) << dterr << "%\n";
5455 std::cout << " Lorentz Angle: " << std::right
5456 << std::setw(10) << std::setprecision(6) << (lor / Pi * 180.)
5457 << " degree +/- " << std::setprecision(2) << lorerr << "%\n";
5458 if (useSST) {
5459 std::cout << " Townsend coefficient (SST): " << std::right
5460 << std::setw(10) << std::setprecision(6) << alpha
5461 << " cm-1 +/- " << std::setprecision(2) << alphaerr << "%\n";
5462 std::cout << " Attachment coefficient (SST): " << std::right
5463 << std::setw(10) << std::setprecision(6) << eta << " cm-1 +/- "
5464 << std::setprecision(2) << etaerr << "%\n";
5465 std::cout << " Eff. Townsend coefficient (TOF): " << std::right
5466 << std::setw(10) << std::setprecision(6) << alphatof
5467 << " cm-1\n";
5468 } else {
5469 std::cout << " Townsend coefficient: " << std::right
5470 << std::setw(10) << std::setprecision(6) << alpha
5471 << " cm-1 +/- " << std::setprecision(2) << alphaerr << "%\n";
5472 std::cout << " Attachment coefficient: " << std::right
5473 << std::setw(10) << std::setprecision(6) << eta << " cm-1 +/- "
5474 << std::setprecision(2) << etaerr << "%\n";
5475 }
5476 }
5477}
double m_fraction[m_nMaxGases]
Definition: MediumGas.hh:91
struct Garfield::Magboltz::@8 diflab_
struct Garfield::Magboltz::@5 ratio_
struct Garfield::Magboltz::@12 ctowns_
void elimitc_(long long *ielow)
struct Garfield::Magboltz::@13 ctwner_
struct Garfield::Magboltz::@4 gasn_
struct Garfield::Magboltz::@0 bfld_
struct Garfield::Magboltz::@7 velerr_
struct Garfield::Magboltz::@2 setp_
void elimit_(long long *ielow)
void elimitb_(long long *ielow)
struct Garfield::Magboltz::@14 tofout_
struct Garfield::Magboltz::@6 vel_
struct Garfield::Magboltz::@11 diferl_

Referenced by GenerateGasTable().

◆ SetExcitationScalingFactor()

void Garfield::MediumMagboltz::SetExcitationScalingFactor ( const double  r,
std::string  gasname 
)

Definition at line 380 of file MediumMagboltz.cc.

381 {
382
383 if (r <= 0.) {
384 std::cerr << m_className << "::SetScalingFactor:\n";
385 std::cerr << " Incorrect value for scaling factor: " << r << "\n";
386 return;
387 }
388
389 // Get the "standard" name of this gas.
390 if (!GetGasName(gasname, gasname)) {
391 std::cerr << m_className << "::SetExcitationScalingFactor:\n";
392 std::cerr << " Unknown gas name.\n";
393 return;
394 }
395
396 // Look for this gas in the present gas mixture.
397 bool found = false;
398 for (unsigned int i = 0; i < m_nComponents; ++i) {
399 if (m_gas[i] == gasname) {
400 m_scaleExc[i] = r;
401 found = true;
402 break;
403 }
404 }
405
406 if (!found) {
407 std::cerr << m_className << "::SetExcitationScalingFactor:\n";
408 std::cerr << " Specified gas (" << gasname
409 << ") is not part of the present gas mixture.\n";
410 return;
411 }
412
413 // Make sure that the collision rate table is updated.
414 m_isChanged = true;
415}

◆ SetMaxElectronEnergy()

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

Definition at line 97 of file MediumMagboltz.cc.

97 {
98
99 if (e <= Small) {
100 std::cerr << m_className << "::SetMaxElectronEnergy:\n";
101 std::cerr << " Provided upper electron energy limit (" << e
102 << " eV) is too small.\n";
103 return false;
104 }
105 m_eFinal = e;
106
107 // Determine the energy interval size.
108 if (m_eFinal <= m_eHigh) {
109 m_eStep = m_eFinal / nEnergySteps;
110 } else {
111 m_eStep = m_eHigh / nEnergySteps;
112 }
113
114 // Set max. energy and step size also in Magboltz common block.
115 Magboltz::inpt_.efinal = m_eFinal;
116 Magboltz::inpt_.estep = m_eStep;
117
118 // Force recalculation of the scattering rates table.
119 m_isChanged = true;
120
121 return true;
122}

Referenced by GetElectronCollision(), and GetElectronCollisionRate().

◆ SetMaxPhotonEnergy()

bool Garfield::MediumMagboltz::SetMaxPhotonEnergy ( const double  e)

Definition at line 124 of file MediumMagboltz.cc.

124 {
125
126 if (e <= Small) {
127 std::cerr << m_className << "::SetMaxPhotonEnergy:\n";
128 std::cerr << " Provided upper photon energy limit (" << e
129 << " eV) is too small.\n";
130 return false;
131 }
132 m_eFinalGamma = e;
133
134 // Determine the energy interval size.
135 m_eStepGamma = m_eFinalGamma / nEnergyStepsGamma;
136
137 // Force recalculation of the scattering rates table.
138 m_isChanged = true;
139
140 return true;
141}

Referenced by GetPhotonCollision(), and GetPhotonCollisionRate().

◆ SetSplittingFunctionFlat()

void Garfield::MediumMagboltz::SetSplittingFunctionFlat ( )

Definition at line 168 of file MediumMagboltz.cc.

168 {
169
170 m_useOpalBeaty = false;
171 m_useGreenSawada = false;
172}

◆ SetSplittingFunctionGreenSawada()

void Garfield::MediumMagboltz::SetSplittingFunctionGreenSawada ( )

Definition at line 149 of file MediumMagboltz.cc.

149 {
150
151 m_useOpalBeaty = false;
152 m_useGreenSawada = true;
153 if (m_isChanged) return;
154
155 bool allset = true;
156 for (unsigned int i = 0; i < m_nComponents; ++i) {
157 if (!m_hasGreenSawada[i]) {
158 if (allset) {
159 std::cout << m_className << "::SetSplittingFunctionGreenSawada:\n";
160 allset = false;
161 }
162 std::cout << " Fit parameters for " << m_gas[i] << " not available.\n";
163 std::cout << " Opal-Beaty formula is used instead.\n";
164 }
165 }
166}

◆ SetSplittingFunctionOpalBeaty()

void Garfield::MediumMagboltz::SetSplittingFunctionOpalBeaty ( )

Definition at line 143 of file MediumMagboltz.cc.

143 {
144
145 m_useOpalBeaty = true;
146 m_useGreenSawada = false;
147}

Member Data Documentation

◆ fit3d4p

double Garfield::MediumMagboltz::fit3d4p

Definition at line 140 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit3dEtaC2H6

double Garfield::MediumMagboltz::fit3dEtaC2H6

Definition at line 143 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit3dEtaCH4

double Garfield::MediumMagboltz::fit3dEtaCH4

Definition at line 142 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit3dEtaCO2

double Garfield::MediumMagboltz::fit3dEtaCO2

Definition at line 141 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit3dQC2H6

double Garfield::MediumMagboltz::fit3dQC2H6

Definition at line 143 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit3dQCH4

double Garfield::MediumMagboltz::fit3dQCH4

Definition at line 142 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit3dQCO2

double Garfield::MediumMagboltz::fit3dQCO2

Definition at line 141 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit4pEtaC2H6

double Garfield::MediumMagboltz::fit4pEtaC2H6

Definition at line 145 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit4pEtaCH4

double Garfield::MediumMagboltz::fit4pEtaCH4

Definition at line 144 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fit4sEtaC2H6

double Garfield::MediumMagboltz::fit4sEtaC2H6

Definition at line 146 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fitHigh4p

double Garfield::MediumMagboltz::fitHigh4p

Definition at line 140 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().

◆ fitLineCut

double Garfield::MediumMagboltz::fitLineCut

Definition at line 147 of file MediumMagboltz.hh.

Referenced by MediumMagboltz().


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