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

#include <MediumGas.hh>

+ Inheritance diagram for Garfield::MediumGas:

Classes

struct  excListElement
 
struct  ionListElement
 

Public Member Functions

 MediumGas ()
 
 ~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
 
bool GetPhotoabsorptionCrossSection (const double &e, double &sigma, const unsigned int &i)
 
- 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 ()
 

Protected Member Functions

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

std::string gas [m_nMaxGases]
 
double fraction [m_nMaxGases]
 
double atWeight [m_nMaxGases]
 
double atNum [m_nMaxGases]
 
bool usePenning
 
double rPenningGlobal
 
double rPenningGas [m_nMaxGases]
 
double lambdaPenningGlobal
 
double lambdaPenningGas [m_nMaxGases]
 
double pressureTable
 
double temperatureTable
 
std::vector< std::vector< std::vector< double > > > tabTownsendNoPenning
 
bool m_hasExcRates
 
bool m_hasIonRates
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabExcRates
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabIonRates
 
int nExcListElements
 
std::vector< excListElementexcitationList
 
int nIonListElements
 
std::vector< ionListElementionisationList
 
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
 
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

static const unsigned int m_nMaxGases = 6
 
- Static Protected Attributes inherited from Garfield::Medium
static int m_idCounter = -1
 

Detailed Description

Definition at line 11 of file MediumGas.hh.

Constructor & Destructor Documentation

◆ MediumGas()

Garfield::MediumGas::MediumGas ( )

Definition at line 16 of file MediumGas.cc.

17 : Medium(),
18 usePenning(false),
23 m_hasExcRates(false),
24 m_hasIonRates(false),
33
34 m_className = "MediumGas";
35
36 // Default gas mixture: pure argon
37 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
38 fraction[i] = 0.;
39 gas[i] = "";
40 atWeight[i] = 0.;
41 atNum[i] = 0.;
42 }
43 gas[0] = "Ar";
44 fraction[0] = 1.;
45 m_name = gas[0];
46 GetGasInfo(gas[0], atWeight[0], atNum[0]);
47
48 m_isChanged = true;
49
52
53 // Initialise Penning parameters
54 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
55 rPenningGas[i] = 0.;
56 lambdaPenningGas[i] = 0.;
57 }
58
59 tabElectronTownsend.clear();
60 excitationList.clear();
61 ionisationList.clear();
62}
std::vector< ionListElement > ionisationList
Definition: MediumGas.hh:128
bool GetGasInfo(const std::string gasname, double &a, double &z) const
Definition: MediumGas.cc:1875
double lambdaPenningGas[m_nMaxGases]
Definition: MediumGas.hh:98
double lambdaPenningGlobal
Definition: MediumGas.hh:97
double fraction[m_nMaxGases]
Definition: MediumGas.hh:86
unsigned int m_intpExcRates
Definition: MediumGas.hh:133
unsigned int m_extrHighExcRates
Definition: MediumGas.hh:131
std::string gas[m_nMaxGases]
Definition: MediumGas.hh:85
static const unsigned int m_nMaxGases
Definition: MediumGas.hh:82
double rPenningGas[m_nMaxGases]
Definition: MediumGas.hh:95
double atNum[m_nMaxGases]
Definition: MediumGas.hh:88
std::vector< excListElement > excitationList
Definition: MediumGas.hh:121
unsigned int m_intpIonRates
Definition: MediumGas.hh:134
unsigned int m_extrHighIonRates
Definition: MediumGas.hh:132
unsigned int m_extrLowIonRates
Definition: MediumGas.hh:132
unsigned int m_extrLowExcRates
Definition: MediumGas.hh:131
double atWeight[m_nMaxGases]
Definition: MediumGas.hh:87
double m_pressure
Definition: Medium.hh:295
std::string m_name
Definition: Medium.hh:291
virtual void EnableDrift()
Definition: Medium.hh:52
virtual void EnablePrimaryIonisation()
Definition: Medium.hh:54
std::string m_className
Definition: Medium.hh:284
bool m_isChanged
Definition: Medium.hh:316
std::vector< std::vector< std::vector< double > > > tabElectronTownsend
Definition: Medium.hh:341
double m_temperature
Definition: Medium.hh:293

◆ ~MediumGas()

Garfield::MediumGas::~MediumGas ( )
inline

Definition at line 17 of file MediumGas.hh.

17{}

Member Function Documentation

◆ GetAtomicNumber()

double Garfield::MediumGas::GetAtomicNumber ( ) const
virtual

Reimplemented from Garfield::Medium.

Definition at line 273 of file MediumGas.cc.

273 {
274
275 // Effective Z, weighted by the fractions of the components.
276 double z = 0.;
277 for (unsigned int i = 0; i < m_nComponents; ++i) {
278 z += atNum[i] * fraction[i];
279 }
280 return z;
281}
unsigned int m_nComponents
Definition: Medium.hh:299

◆ GetAtomicWeight()

double Garfield::MediumGas::GetAtomicWeight ( ) const
virtual

Reimplemented from Garfield::Medium.

Definition at line 251 of file MediumGas.cc.

251 {
252
253 // Effective A, weighted by the fractions of the components.
254 double a = 0.;
255 for (unsigned int i = 0; i < m_nComponents; ++i) {
256 a += atWeight[i] * fraction[i];
257 }
258 return a;
259}

Referenced by GetMassDensity().

◆ GetComponent()

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

Reimplemented from Garfield::Medium.

Definition at line 205 of file MediumGas.cc.

206 {
207
208 if (i >= m_nComponents) {
209 std::cerr << m_className << "::GetComponent:\n";
210 std::cerr << " Index out of range.\n";
211 label = "";
212 f = 0.;
213 return;
214 }
215
216 label = gas[i];
217 f = fraction[i];
218}

◆ GetComposition()

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

Definition at line 186 of file MediumGas.cc.

189 {
190
191 gas1 = gas[0];
192 f1 = fraction[0];
193 gas2 = gas[1];
194 f2 = fraction[1];
195 gas3 = gas[2];
196 f3 = fraction[2];
197 gas4 = gas[3];
198 f4 = fraction[3];
199 gas5 = gas[4];
200 f5 = fraction[4];
201 gas6 = gas[5];
202 f6 = fraction[5];
203}

◆ GetGasInfo()

bool Garfield::MediumGas::GetGasInfo ( const std::string  gasname,
double &  a,
double &  z 
) const
protected

Definition at line 1875 of file MediumGas.cc.

1876 {
1877
1878 if (gasname == "CF4") {
1879 a = 12.0107 + 4 * 18.9984032;
1880 z = 6 + 4 * 9;
1881 return true;
1882 } else if (gasname == "Ar") {
1883 a = 39.948;
1884 z = 18;
1885 } else if (gasname == "He") {
1886 a = 4.002602;
1887 z = 2;
1888 } else if (gasname == "He-3") {
1889 a = 3.01602931914;
1890 z = 2;
1891 } else if (gasname == "Ne") {
1892 a = 20.1797;
1893 z = 10;
1894 } else if (gasname == "Kr") {
1895 a = 37.798;
1896 z = 36;
1897 } else if (gasname == "Xe") {
1898 a = 131.293;
1899 z = 54;
1900 } else if (gasname == "CH4") {
1901 a = 12.0107 + 4 * 1.00794;
1902 z = 6 + 4;
1903 } else if (gasname == "C2H6") {
1904 a = 2 * 12.0107 + 6 * 1.00794;
1905 z = 2 * 6 + 6;
1906 } else if (gasname == "C3H8") {
1907 a = 3 * 12.0107 + 8 * 1.00794;
1908 z = 3 * 6 + 8;
1909 } else if (gasname == "iC4H10") {
1910 a = 4 * 12.0107 + 10 * 1.00794;
1911 z = 4 * 6 + 10;
1912 } else if (gasname == "CO2") {
1913 a = 12.0107 + 2 * 15.9994;
1914 z = 6 + 2 * 8;
1915 } else if (gasname == "neoC5H12") {
1916 a = 5 * 12.0107 + 12 * 1.00794;
1917 z = 5 * 6 + 12;
1918 } else if (gasname == "H2O") {
1919 a = 2 * 1.00794 + 15.9994;
1920 z = 2 + 8;
1921 } else if (gasname == "O2") {
1922 a = 2 * 15.9994;
1923 z = 2 * 8;
1924 } else if (gasname == "N2") {
1925 a = 2 * 14.0067;
1926 z = 2 * 7;
1927 } else if (gasname == "NO") {
1928 a = 14.0067 + 15.9994;
1929 z = 7 + 8;
1930 } else if (gasname == "N2O") {
1931 a = 2 * 14.0067 + 15.9994;
1932 z = 2 * 7 + 8;
1933 } else if (gasname == "C2H4") {
1934 a = 2 * 12.0107 + 4 * 1.00794;
1935 z = 2 * 6 + 4;
1936 } else if (gasname == "C2H2") {
1937 a = 2 * 12.0107 + 2 * 1.00794;
1938 z = 2 * 6 + 2;
1939 } else if (gasname == "H2") {
1940 a = 2 * 1.00794;
1941 z = 2;
1942 } else if (gasname == "D2") {
1943 a = 2 * 2.01410177785;
1944 z = 2;
1945 } else if (gasname == "CO") {
1946 a = 12.0107 + 15.9994;
1947 z = 6 + 8;
1948 } else if (gasname == "Methylal") {
1949 a = 3 * 12.0107 + 8 * 1.00794 + 2 * 15.9994;
1950 z = 3 * 6 + 8 + 2 * 8;
1951 } else if (gasname == "DME") {
1952 a = 4 * 12.0107 + 10 * 1.00794 + 2 * 15.9994;
1953 z = 4 * 6 + 10 + 2 * 8;
1954 } else if (gasname == "Reid-Step" || gasname == "Mawell-Model" ||
1955 gasname == "Reid-Ramp") {
1956 a = 1.;
1957 z = 1.;
1958 } else if (gasname == "C2F6") {
1959 a = 2 * 12.0107 + 6 * 18.9984032;
1960 z = 2 * 6 + 6 * 9;
1961 } else if (gasname == "SF6") {
1962 a = 32.065 + 6 * 18.9984032;
1963 z = 16 + 6 * 9;
1964 } else if (gasname == "NH3") {
1965 a = 14.0067 + 3 * 1.00794;
1966 z = 7 + 3;
1967 } else if (gasname == "C3H6") {
1968 a = 3 * 12.0107 + 6 * 1.00794;
1969 z = 3 * 6 + 6;
1970 } else if (gasname == "cC3H6") {
1971 a = 3 * 12.0107 + 6 * 1.00794;
1972 z = 3 * 6 + 6;
1973 } else if (gasname == "CH3OH") {
1974 a = 12.0107 + 4 * 1.00794 + 15.9994;
1975 z = 6 + 4 + 8;
1976 } else if (gasname == "C2H5OH") {
1977 a = 2 * 12.0107 + 6 * 1.00794 + 15.9994;
1978 z = 2 * 6 + 6 + 8;
1979 } else if (gasname == "C3H7OH") {
1980 a = 3 * 12.0107 + 8 * 1.00794 + 15.9994;
1981 z = 3 * 6 + 8 * 8;
1982 } else if (gasname == "Cs") {
1983 a = 132.9054519;
1984 z = 55;
1985 } else if (gasname == "F2") {
1986 a = 2 * 18.9984032;
1987 z = 2 * 9;
1988 } else if (gasname == "CS2") {
1989 a = 12.0107 + 2 * 32.065;
1990 z = 6 + 2 * 16;
1991 } else if (gasname == "COS") {
1992 a = 12.0107 + 15.9994 + 32.065;
1993 z = 6 + 8 + 16;
1994 } else if (gasname == "CD4") {
1995 a = 12.0107 + 4 * 2.01410177785;
1996 z = 6 + 4;
1997 } else if (gasname == "BF3") {
1998 a = 10.811 + 3 * 18.9984032;
1999 z = 5 + 3 * 9;
2000 } else if (gasname == "C2H2F4") {
2001 a = 2 * 12.0107 + 2 * 1.00794 + 4 * 18.9984032;
2002 z = 2 * 6 + 2 + 4 * 9;
2003 } else if (gasname == "CHF3") {
2004 a = 12.0107 + 1.00794 + 3 * 18.9984032;
2005 z = 6 + 1 + 3 * 9;
2006 } else if (gasname == "CF3Br") {
2007 a = 12.0107 + 3 * 18.9984032 + 79.904;
2008 z = 6 + 3 * 9 + 35;
2009 } else if (gasname == "C3F8") {
2010 a = 3 * 12.0107 + 8 * 18.9984032;
2011 z = 3 * 6 + 8 * 9;
2012 } else if (gasname == "O3") {
2013 a = 3 * 15.9994;
2014 z = 3 * 8;
2015 } else if (gasname == "Hg") {
2016 a = 2 * 200.59;
2017 z = 80;
2018 } else if (gasname == "H2S") {
2019 a = 2 * 1.00794 + 32.065;
2020 z = 2 + 16;
2021 } else if (gasname == "nC4H10") {
2022 a = 4 * 12.0107 + 10 * 1.00794;
2023 z = 4 * 6 + 10;
2024 } else if (gasname == "nC5H12") {
2025 a = 5 * 12.0107 + 12 * 1.00794;
2026 z = 5 * 6 + 12;
2027 } else if (gasname == "N2") {
2028 a = 2 * 14.0067;
2029 z = 2 * 7;
2030 } else if (gasname == "GeH4") {
2031 a = 72.64 + 4 * 1.00794;
2032 z = 32 + 4;
2033 } else if (gasname == "SiH4") {
2034 a = 28.0855 + 4 * 1.00794;
2035 z = 14 + 4;
2036 } else {
2037 a = 0.;
2038 z = 0.;
2039 return false;
2040 }
2041
2042 return true;
2043}

Referenced by LoadGasFile(), MediumGas(), and SetComposition().

◆ GetGasName() [1/2]

bool Garfield::MediumGas::GetGasName ( const int  gasnumber,
const int  version,
std::string &  gasname 
)
protected

Definition at line 2045 of file MediumGas.cc.

2046 {
2047
2048 switch (gasnumber) {
2049 case 1:
2050 gasname = "CF4";
2051 break;
2052 case 2:
2053 gasname = "Ar";
2054 break;
2055 case 3:
2056 gasname = "He";
2057 break;
2058 case 4:
2059 gasname = "He-3";
2060 break;
2061 case 5:
2062 gasname = "Ne";
2063 break;
2064 case 6:
2065 gasname = "Kr";
2066 break;
2067 case 7:
2068 gasname = "Xe";
2069 break;
2070 case 8:
2071 gasname = "CH4";
2072 break;
2073 case 9:
2074 gasname = "C2H6";
2075 break;
2076 case 10:
2077 gasname = "C3H8";
2078 break;
2079 case 11:
2080 gasname = "iC4H10";
2081 break;
2082 case 12:
2083 gasname = "CO2";
2084 break;
2085 case 13:
2086 gasname = "neoC5H12";
2087 break;
2088 case 14:
2089 gasname = "H2O";
2090 break;
2091 case 15:
2092 gasname = "O2";
2093 break;
2094 case 16:
2095 gasname = "N2";
2096 break;
2097 case 17:
2098 gasname = "NO";
2099 break;
2100 case 18:
2101 gasname = "N2O";
2102 break;
2103 case 19:
2104 gasname = "C2H4";
2105 break;
2106 case 20:
2107 gasname = "C2H2";
2108 break;
2109 case 21:
2110 gasname = "H2";
2111 break;
2112 case 22:
2113 gasname = "D2";
2114 break;
2115 case 23:
2116 gasname = "CO";
2117 break;
2118 case 24:
2119 gasname = "Methylal";
2120 break;
2121 case 25:
2122 gasname = "DME";
2123 break;
2124 case 26:
2125 gasname = "Reid-Step";
2126 break;
2127 case 27:
2128 gasname = "Maxwell-Model";
2129 break;
2130 case 28:
2131 gasname = "Reid-Ramp";
2132 break;
2133 case 29:
2134 gasname = "C2F6";
2135 break;
2136 case 30:
2137 gasname = "SF6";
2138 break;
2139 case 31:
2140 gasname = "NH3";
2141 break;
2142 case 32:
2143 gasname = "C3H6";
2144 break;
2145 case 33:
2146 gasname = "cC3H6";
2147 break;
2148 case 34:
2149 gasname = "CH3OH";
2150 break;
2151 case 35:
2152 gasname = "C2H5OH";
2153 break;
2154 case 36:
2155 gasname = "C3H7OH";
2156 break;
2157 case 37:
2158 gasname = "Cs";
2159 break;
2160 case 38:
2161 gasname = "F2";
2162 break;
2163 case 39:
2164 gasname = "CS2";
2165 break;
2166 case 40:
2167 gasname = "COS";
2168 break;
2169 case 41:
2170 gasname = "CD4";
2171 break;
2172 case 42:
2173 gasname = "BF3";
2174 break;
2175 case 43:
2176 gasname = "C2H2F4";
2177 break;
2178 case 44:
2179 if (version <= 11) {
2180 gasname = "He-3";
2181 } else {
2182 gasname = "TMA";
2183 }
2184 break;
2185 case 45:
2186 gasname = "He";
2187 break;
2188 case 46:
2189 gasname = "Ne";
2190 break;
2191 case 47:
2192 gasname = "Ar";
2193 break;
2194 case 48:
2195 gasname = "Kr";
2196 break;
2197 case 49:
2198 gasname = "Xe";
2199 break;
2200 case 50:
2201 gasname = "CHF3";
2202 break;
2203 case 51:
2204 gasname = "CF3Br";
2205 break;
2206 case 52:
2207 gasname = "C3F8";
2208 break;
2209 case 53:
2210 gasname = "O3";
2211 break;
2212 case 54:
2213 gasname = "Hg";
2214 break;
2215 case 55:
2216 gasname = "H2S";
2217 break;
2218 case 56:
2219 gasname = "nC4H10";
2220 break;
2221 case 57:
2222 gasname = "nC5H12";
2223 break;
2224 case 58:
2225 gasname = "N2";
2226 break;
2227 case 59:
2228 gasname = "GeH4";
2229 break;
2230 case 60:
2231 gasname = "SiH4";
2232 break;
2233 default:
2234 gasname = "";
2235 return false;
2236 break;
2237 }
2238 return true;
2239}

Referenced by Garfield::MediumMagboltz::DisablePenningTransfer(), Garfield::MediumMagboltz::EnablePenningTransfer(), LoadGasFile(), SetComposition(), and Garfield::MediumMagboltz::SetExcitationScalingFactor().

◆ GetGasName() [2/2]

bool Garfield::MediumGas::GetGasName ( std::string  input,
std::string &  gasname 
) const
protected

Definition at line 2241 of file MediumGas.cc.

2241 {
2242
2243 // Convert to upper-case
2244 for (unsigned int i = 0; i < input.length(); ++i) {
2245 input[i] = toupper(input[i]);
2246 }
2247
2248 gasname = "";
2249
2250 if (input == "") return false;
2251
2252 // CF4
2253 if (input == "CF4" || input == "FREON" || input == "FREON-14" ||
2254 input == "TETRAFLUOROMETHANE") {
2255 gasname = "CF4";
2256 return true;
2257 }
2258 // Argon
2259 if (input == "AR" || input == "ARGON") {
2260 gasname = "Ar";
2261 return true;
2262 }
2263 // Helium 4
2264 if (input == "HE" || input == "HELIUM" || input == "HE-4" ||
2265 input == "HE 4" || input == "HE4" || input == "4-HE" || input == "4 HE" ||
2266 input == "4HE" || input == "HELIUM-4" || input == "HELIUM 4" ||
2267 input == "HELIUM4") {
2268 gasname = "He";
2269 return true;
2270 }
2271 // Helium 3
2272 if (input == "HE-3" || input == "HE3" || input == "HELIUM-3" ||
2273 input == "HELIUM 3" || input == "HELIUM3") {
2274 gasname = "He-3";
2275 return true;
2276 }
2277 // Neon
2278 if (input == "NE" || input == "NEON") {
2279 gasname = "Ne";
2280 return true;
2281 }
2282 // Krypton
2283 if (input == "KR" || input == "KRYPTON") {
2284 gasname = "Kr";
2285 return true;
2286 }
2287 // Xenon
2288 if (input == "XE" || input == "XENON") {
2289 gasname = "Xe";
2290 return true;
2291 }
2292 // Methane
2293 if (input == "CH4" || input == "METHANE") {
2294 gasname = "CH4";
2295 return true;
2296 }
2297 // Ethane
2298 if (input == "C2H6" || input == "ETHANE") {
2299 gasname = "C2H6";
2300 return true;
2301 }
2302 // Propane
2303 if (input == "C3H8" || input == "PROPANE") {
2304 gasname = "C3H8";
2305 return true;
2306 }
2307 // Isobutane
2308 if (input == "C4H10" || input == "ISOBUTANE" || input == "ISO" ||
2309 input == "IC4H10" || input == "ISO-C4H10" || input == "ISOC4H10") {
2310 gasname = "iC4H10";
2311 return true;
2312 }
2313 // Carbon dioxide (CO2)
2314 if (input == "CO2" || input == "CARBON-DIOXIDE" ||
2315 input == "CARBON DIOXIDE" || input == "CARBONDIOXIDE") {
2316 gasname = "CO2";
2317 return true;
2318 }
2319 // Neopentane
2320 if (input == "NEOPENTANE" || input == "NEO-PENTANE" || input == "NEO-C5H12" ||
2321 input == "NEOC5H12" || input == "DIMETHYLPROPANE" || input == "C5H12") {
2322 gasname = "neoC5H12";
2323 return true;
2324 }
2325 // Water
2326 if (input == "H2O" || input == "WATER" || input == "WATER-VAPOUR" ||
2327 input == "WATER VAPOUR") {
2328 gasname = "H2O";
2329 return true;
2330 }
2331 // Oxygen
2332 if (input == "O2" || input == "OXYGEN") {
2333 gasname = "O2";
2334 return true;
2335 }
2336 // Nitrogen
2337 if (input == "NI" || input == "NITRO" || input == "N2" ||
2338 input == "NITROGEN") {
2339 gasname = "N2";
2340 return true;
2341 }
2342 // Nitric oxide (NO)
2343 if (input == "NO" || input == "NITRIC-OXIDE" || input == "NITRIC OXIDE" ||
2344 input == "NITROGEN-MONOXIDE" || input == "NITROGEN MONOXIDE") {
2345 gasname = "NO";
2346 return true;
2347 }
2348 // Nitrous oxide (N2O)
2349 if (input == "N2O" || input == "NITROUS-OXIDE" || input == "NITROUS OXIDE" ||
2350 input == "DINITROGEN-MONOXIDE" || input == "LAUGHING-GAS") {
2351 gasname = "N2O";
2352 return true;
2353 }
2354 // Ethene (C2H4)
2355 if (input == "C2H4" || input == "ETHENE" || input == "ETHYLENE") {
2356 gasname = "C2H4";
2357 return true;
2358 }
2359 // Acetylene (C2H2)
2360 if (input == "C2H2" || input == "ACETYL" || input == "ACETYLENE" ||
2361 input == "ETHYNE") {
2362 gasname = "C2H2";
2363 return true;
2364 }
2365 // Hydrogen
2366 if (input == "H2" || input == "HYDROGEN") {
2367 gasname = "H2";
2368 return true;
2369 }
2370 // Deuterium
2371 if (input == "D2" || input == "DEUTERIUM") {
2372 gasname = "D2";
2373 return true;
2374 }
2375 // Carbon monoxide (CO)
2376 if (input == "CO" || input == "CARBON-MONOXIDE" ||
2377 input == "CARBON MONOXIDE") {
2378 gasname = "CO";
2379 return true;
2380 }
2381 // Methylal (dimethoxymethane, CH3-O-CH2-O-CH3, "hot" version)
2382 if (input == "METHYLAL" || input == "METHYLAL-HOT" || input == "DMM" ||
2383 input == "DIMETHOXYMETHANE" || input == "FORMAL" || input == "C3H8O2") {
2384 gasname = "Methylal";
2385 return true;
2386 }
2387 // DME
2388 if (input == "DME" || input == "DIMETHYL-ETHER" || input == "DIMETHYLETHER" ||
2389 input == "DIMETHYL ETHER" || input == "METHYL ETHER" ||
2390 input == "METHYL-ETHER" || input == "METHYLETHER" ||
2391 input == "WOOD-ETHER" || input == "WOODETHER" || input == "WOOD ETHER" ||
2392 input == "DIMETHYL OXIDE" || input == "DIMETHYL-OXIDE" ||
2393 input == "DEMEON" || input == "METHOXYMETHANE" || input == "C4H10O2") {
2394 gasname = "DME";
2395 return true;
2396 }
2397 // Reid step
2398 if (input == "REID-STEP") {
2399 gasname = "Reid-Step";
2400 return true;
2401 }
2402 // Maxwell model
2403 if (input == "MAXWELL-MODEL") {
2404 gasname = "Maxwell-Model";
2405 return true;
2406 }
2407 // Reid ramp
2408 if (input == "REID-RAMP") {
2409 gasname = "Reid-Ramp";
2410 return true;
2411 }
2412 // C2F6
2413 if (input == "C2F6" || input == "FREON-116" || input == "ZYRON-116" ||
2414 input == "ZYRON-116-N5" || input == "HEXAFLUOROETHANE") {
2415 gasname = "C2F6";
2416 return true;
2417 }
2418 // SF6
2419 if (input == "SF6" || input == "SULPHUR-HEXAFLUORIDE" ||
2420 input == "SULFUR-HEXAFLUORIDE" || input == "SULPHUR HEXAFLUORIDE" ||
2421 input == "SULFUR HEXAFLUORIDE") {
2422 gasname = "SF6";
2423 return true;
2424 }
2425 // NH3
2426 if (input == "NH3" || input == "AMMONIA") {
2427 gasname = "NH3";
2428 return true;
2429 }
2430 // Propene
2431 if (input == "C3H6" || input == "PROPENE" || input == "PROPYLENE") {
2432 gasname = "C3H6";
2433 return true;
2434 }
2435 // Cyclopropane
2436 if (input == "C-PROPANE" || input == "CYCLO-PROPANE" ||
2437 input == "CYCLO PROPANE" || input == "CYCLOPROPANE" ||
2438 input == "C-C3H6" || input == "CC3H6" || input == "CYCLO-C3H6") {
2439 gasname = "cC3H6";
2440 return true;
2441 }
2442 // Methanol
2443 if (input == "METHANOL" || input == "METHYL-ALCOHOL" ||
2444 input == "METHYL ALCOHOL" || input == "WOOD ALCOHOL" ||
2445 input == "WOOD-ALCOHOL" || input == "CH3OH") {
2446 gasname = "CH3OH";
2447 return true;
2448 }
2449 // Ethanol
2450 if (input == "ETHANOL" || input == "ETHYL-ALCOHOL" ||
2451 input == "ETHYL ALCOHOL" || input == "GRAIN ALCOHOL" ||
2452 input == "GRAIN-ALCOHOL" || input == "C2H5OH") {
2453 gasname = "C2H5OH";
2454 return true;
2455 }
2456 // Propanol
2457 if (input == "PROPANOL" || input == "2-PROPANOL" || input == "ISOPROPYL" ||
2458 input == "ISO-PROPANOL" || input == "ISOPROPANOL" ||
2459 input == "ISOPROPYL ALCOHOL" || input == "ISOPROPYL-ALCOHOL" ||
2460 input == "C3H7OH") {
2461 gasname = "C3H7OH";
2462 return true;
2463 }
2464 // Cesium / Caesium.
2465 if (input == "CS" || input == "CESIUM" || input == "CAESIUM") {
2466 gasname = "Cs";
2467 return true;
2468 }
2469 // Fluorine
2470 if (input == "F2" || input == "FLUOR" || input == "FLUORINE") {
2471 gasname = "F2";
2472 return true;
2473 }
2474 // CS2
2475 if (input == "CS2" || input == "CARBON-DISULPHIDE" ||
2476 input == "CARBON-DISULFIDE" || input == "CARBON DISULPHIDE" ||
2477 input == "CARBON DISULFIDE") {
2478 gasname = "CS2";
2479 return true;
2480 }
2481 // COS
2482 if (input == "COS" || input == "CARBONYL-SULPHIDE" ||
2483 input == "CARBONYL-SULFIDE" || input == "CARBONYL SULFIDE") {
2484 gasname = "COS";
2485 return true;
2486 }
2487 // Deuterated methane
2488 if (input == "DEUT-METHANE" || input == "DEUTERIUM-METHANE" ||
2489 input == "DEUTERATED-METHANE" || input == "DEUTERATED METHANE" ||
2490 input == "DEUTERIUM METHANE" || input == "CD4") {
2491 gasname = "CD4";
2492 return true;
2493 }
2494 // BF3
2495 if (input == "BF3" || input == "BORON-TRIFLUORIDE" ||
2496 input == "BORON TRIFLUORIDE") {
2497 gasname = "BF3";
2498 return true;
2499 }
2500 // C2H2F4 (and C2HF5).
2501 if (input == "C2HF5" || input == "C2H2F4" || input == "C2F5H" ||
2502 input == "C2F4H2" || input == "FREON 134" || input == "FREON 134A" ||
2503 input == "FREON-134" || input == "FREON-134-A" || input == "FREON 125" ||
2504 input == "ZYRON 125" || input == "FREON-125" || input == "ZYRON-125" ||
2505 input == "TETRAFLUOROETHANE" || input == "PENTAFLUOROETHANE") {
2506 gasname = "C2H2F4";
2507 return true;
2508 }
2509 // TMA
2510 if (input == "TMA" || input == "TRIMETHYLAMINE" || input == "N(CH3)3" ||
2511 input == "N-(CH3)3") {
2512 gasname = "TMA";
2513 return true;
2514 }
2515 // CHF3
2516 if (input == "CHF3" || input == "FREON-23" || input == "TRIFLUOROMETHANE" ||
2517 input == "FLUOROFORM") {
2518 gasname = "CHF3";
2519 return true;
2520 }
2521 // CF3Br
2522 if (input == "CF3BR" || input == "TRIFLUOROBROMOMETHANE" ||
2523 input == "BROMOTRIFLUOROMETHANE" || input == "HALON-1301" ||
2524 input == "HALON 1301" || input == "FREON-13B1" || input == "FREON 13BI") {
2525 gasname = "CF3Br";
2526 return true;
2527 }
2528 // C3F8
2529 if (input == "C3F8" || input == "OCTAFLUOROPROPANE" || input == "R218" ||
2530 input == "R-218" || input == "FREON 218" || input == "FREON-218" ||
2531 input == "PERFLUOROPROPANE" || input == "RC 218" || input == "PFC 218" ||
2532 input == "RC-218" || input == "PFC-218" || input == "FLUTEC PP30" ||
2533 input == "GENETRON 218") {
2534 gasname = "C3F8";
2535 return true;
2536 }
2537 // Ozone
2538 if (input == "OZONE" || input == "O3") {
2539 gasname = "O3";
2540 return true;
2541 }
2542 // Mercury
2543 if (input == "MERCURY" || input == "HG" || input == "HG2") {
2544 gasname = "Hg";
2545 return true;
2546 }
2547 // H2S
2548 if (input == "H2S" || input == "HYDROGEN SULPHIDE" || input == "SEWER GAS" ||
2549 input == "HYDROGEN-SULPHIDE" || input == "SEWER-GAS" ||
2550 input == "HYDROGEN SULFIDE" || input == "HEPATIC ACID" ||
2551 input == "HYDROGEN-SULFIDE" || input == "HEPATIC-ACID" ||
2552 input == "SULFUR HYDRIDE" || input == "DIHYDROGEN MONOSULFIDE" ||
2553 input == "SULFUR-HYDRIDE" || input == "DIHYDROGEN-MONOSULFIDE" ||
2554 input == "DIHYDROGEN MONOSULPHIDE" || input == "SULPHUR HYDRIDE" ||
2555 input == "DIHYDROGEN-MONOSULPHIDE" || input == "SULPHUR-HYDRIDE" ||
2556 input == "STINK DAMP" || input == "SULFURATED HYDROGEN" ||
2557 input == "STINK-DAMP" || input == "SULFURATED-HYDROGEN") {
2558 gasname = "H2S";
2559 return true;
2560 }
2561 // n-Butane
2562 if (input == "N-BUTANE" || input == "N-C4H10" || input == "NBUTANE" ||
2563 input == "NC4H10") {
2564 gasname = "nC4H10";
2565 return true;
2566 }
2567 // n-Pentane
2568 if (input == "N-PENTANE" || input == "N-C5H12" || input == "NPENTANE" ||
2569 input == "NC5H12") {
2570 gasname = "nC5H12";
2571 return true;
2572 }
2573 // Nitrogen
2574 if (input == "NI-PHELPS" || input == "NI PHELPS" ||
2575 input == "NITROGEN-PHELPS" || input == "NITROGEN PHELPHS" ||
2576 input == "N2-PHELPS" || input == "N2 PHELPS" || input == "N2 (PHELPS)") {
2577 gasname = "N2 (Phelps)";
2578 return true;
2579 }
2580 // Germane, GeH4
2581 if (input == "GERMANE" || input == "GERM" || input == "GERMANIUM-HYDRIDE" ||
2582 input == "GERMANIUM HYDRIDE" || input == "GERMANIUM TETRAHYDRIDE" ||
2583 input == "GERMANIUM-TETRAHYDRIDE" || input == "GERMANOMETHANE" ||
2584 input == "MONOGERMANE" || input == "GEH4") {
2585 gasname = "GeH4";
2586 return true;
2587 }
2588 // Silane, SiH4
2589 if (input == "SILANE" || input == "SIL" || input == "SILICON-HYDRIDE" ||
2590 input == "SILICON HYDRIDE" || input == "SILICON-TETRAHYDRIDE" ||
2591 input == "SILICANE" || input == "MONOSILANE" || input == "SIH4") {
2592 gasname = "SiH4";
2593 return true;
2594 }
2595
2596 std::cerr << m_className << "::GetGasName:\n";
2597 std::cerr << " Gas " << input << " is not defined.\n";
2598 return false;
2599}

◆ GetGasNumberGasFile()

bool Garfield::MediumGas::GetGasNumberGasFile ( const std::string  input,
int &  number 
) const
protected

Definition at line 2601 of file MediumGas.cc.

2602 {
2603
2604 if (input == "") {
2605 number = 0;
2606 return false;
2607 }
2608
2609 // CF4
2610 if (input == "CF4") {
2611 number = 1;
2612 return true;
2613 }
2614 // Argon
2615 if (input == "Ar") {
2616 number = 2;
2617 return true;
2618 }
2619 // Helium 4
2620 if (input == "He" || input == "He-4") {
2621 number = 3;
2622 return true;
2623 }
2624 // Helium 3
2625 if (input == "He-3") {
2626 number = 4;
2627 return true;
2628 }
2629 // Neon
2630 if (input == "Ne") {
2631 number = 5;
2632 return true;
2633 }
2634 // Krypton
2635 if (input == "Kr") {
2636 number = 6;
2637 return true;
2638 }
2639 // Xenon
2640 if (input == "Xe") {
2641 number = 7;
2642 return true;
2643 }
2644 // Methane
2645 if (input == "CH4") {
2646 number = 8;
2647 return true;
2648 }
2649 // Ethane
2650 if (input == "C2H6") {
2651 number = 9;
2652 return true;
2653 }
2654 // Propane
2655 if (input == "C3H8") {
2656 number = 10;
2657 return true;
2658 }
2659 // Isobutane
2660 if (input == "iC4H10") {
2661 number = 11;
2662 return true;
2663 }
2664 // Carbon dioxide (CO2)
2665 if (input == "CO2") {
2666 number = 12;
2667 return true;
2668 }
2669 // Neopentane
2670 if (input == "neoC5H12") {
2671 number = 13;
2672 return true;
2673 }
2674 // Water
2675 if (input == "H2O") {
2676 number = 14;
2677 return true;
2678 }
2679 // Oxygen
2680 if (input == "O2") {
2681 number = 15;
2682 return true;
2683 }
2684 // Nitrogen
2685 if (input == "N2") {
2686 number = 16;
2687 return true;
2688 }
2689 // Nitric oxide (NO)
2690 if (input == "NO") {
2691 number = 17;
2692 return true;
2693 }
2694 // Nitrous oxide (N2O)
2695 if (input == "N2O") {
2696 number = 18;
2697 return true;
2698 }
2699 // Ethene (C2H4)
2700 if (input == "C2H4") {
2701 number = 19;
2702 return true;
2703 }
2704 // Acetylene (C2H2)
2705 if (input == "C2H2") {
2706 number = 20;
2707 return true;
2708 }
2709 // Hydrogen
2710 if (input == "H2") {
2711 number = 21;
2712 return true;
2713 }
2714 // Deuterium
2715 if (input == "D2") {
2716 number = 22;
2717 return true;
2718 }
2719 // Carbon monoxide (CO)
2720 if (input == "CO") {
2721 number = 23;
2722 return true;
2723 }
2724 // Methylal (dimethoxymethane, CH3-O-CH2-O-CH3, "hot" version)
2725 if (input == "Methylal") {
2726 number = 24;
2727 return true;
2728 }
2729 // DME
2730 if (input == "DME") {
2731 number = 25;
2732 return true;
2733 }
2734 // Reid step
2735 if (input == "Reid-Step") {
2736 number = 26;
2737 return true;
2738 }
2739 // Maxwell model
2740 if (input == "Maxwell-Model") {
2741 number = 27;
2742 return true;
2743 }
2744 // Reid ramp
2745 if (input == "Reid-Ramp") {
2746 number = 28;
2747 return true;
2748 }
2749 // C2F6
2750 if (input == "C2F6") {
2751 number = 29;
2752 return true;
2753 }
2754 // SF6
2755 if (input == "SF6") {
2756 number = 30;
2757 return true;
2758 }
2759 // NH3
2760 if (input == "NH3") {
2761 number = 31;
2762 return true;
2763 }
2764 // Propene
2765 if (input == "C3H6") {
2766 number = 32;
2767 return true;
2768 }
2769 // Cyclopropane
2770 if (input == "cC3H6") {
2771 number = 33;
2772 return true;
2773 }
2774 // Methanol
2775 if (input == "CH3OH") {
2776 number = 34;
2777 return true;
2778 }
2779 // Ethanol
2780 if (input == "C2H5OH") {
2781 number = 35;
2782 return true;
2783 }
2784 // Propanol
2785 if (input == "C3H7OH") {
2786 number = 36;
2787 return true;
2788 }
2789 // Cesium / Caesium.
2790 if (input == "Cs") {
2791 number = 37;
2792 return true;
2793 }
2794 // Fluorine
2795 if (input == "F2") {
2796 number = 38;
2797 return true;
2798 }
2799 // CS2
2800 if (input == "CS2") {
2801 number = 39;
2802 return true;
2803 }
2804 // COS
2805 if (input == "COS") {
2806 number = 40;
2807 return true;
2808 }
2809 // Deuterated methane
2810 if (input == "CD4") {
2811 number = 41;
2812 return true;
2813 }
2814 // BF3
2815 if (input == "BF3") {
2816 number = 42;
2817 return true;
2818 }
2819 // C2HF5 and C2H2F4.
2820 if (input == "C2HF5" || input == "C2H2F4") {
2821 number = 43;
2822 return true;
2823 }
2824 // TMA
2825 if (input == "TMA") {
2826 number = 44;
2827 return true;
2828 }
2829 // CHF3
2830 if (input == "CHF3") {
2831 number = 50;
2832 return true;
2833 }
2834 // CF3Br
2835 if (input == "CF3Br") {
2836 number = 51;
2837 return true;
2838 }
2839 // C3F8
2840 if (input == "C3F8") {
2841 number = 52;
2842 return true;
2843 }
2844 // Ozone
2845 if (input == "O3") {
2846 number = 53;
2847 return true;
2848 }
2849 // Mercury
2850 if (input == "Hg") {
2851 number = 54;
2852 return true;
2853 }
2854 // H2S
2855 if (input == "H2S") {
2856 number = 55;
2857 return true;
2858 }
2859 // n-butane
2860 if (input == "nC4H10") {
2861 number = 56;
2862 return true;
2863 }
2864 // n-pentane
2865 if (input == "nC5H12") {
2866 number = 57;
2867 return true;
2868 }
2869 // Nitrogen
2870 if (input == "N2 (Phelps)") {
2871 number = 58;
2872 return true;
2873 }
2874 // Germane, GeH4
2875 if (input == "GeH4") {
2876 number = 59;
2877 return true;
2878 }
2879 // Silane, SiH4
2880 if (input == "SiH4") {
2881 number = 60;
2882 return true;
2883 }
2884
2885 std::cerr << m_className << "::GetGasNumberGasFile:\n";
2886 std::cerr << " Gas " << input << " is not defined.\n";
2887 return false;
2888}

Referenced by WriteGasFile().

◆ GetMassDensity()

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

Reimplemented from Garfield::Medium.

Definition at line 268 of file MediumGas.cc.

268 {
269
270 return GetNumberDensity() * GetAtomicWeight() * AtomicMassUnit;
271}
double GetAtomicWeight() const
Definition: MediumGas.cc:251
double GetNumberDensity() const
Definition: MediumGas.cc:261

◆ GetNumberDensity()

double Garfield::MediumGas::GetNumberDensity ( ) const
virtual

Reimplemented from Garfield::Medium.

Definition at line 261 of file MediumGas.cc.

261 {
262
263 // Ideal gas law.
264 return LoschmidtNumber * (m_pressure / AtmosphericPressure) *
265 (ZeroCelsius / m_temperature);
266}

Referenced by GetMassDensity(), and LoadIonMobility().

◆ GetPhotoabsorptionCrossSection()

bool Garfield::MediumGas::GetPhotoabsorptionCrossSection ( const double &  e,
double &  sigma,
const unsigned int &  i 
)

Definition at line 2890 of file MediumGas.cc.

2891 {
2892
2893 if (i >= m_nMaxGases) {
2894 std::cerr << m_className << "::GetPhotoabsorptionCrossSection:\n";
2895 std::cerr << " Index (" << i << ") out of range.\n";
2896 return false;
2897 }
2898
2899 OpticalData optData;
2900 if (!optData.IsAvailable(gas[i])) return false;
2901 double eta = 0.;
2902 return optData.GetPhotoabsorptionCrossSection(gas[i], e, sigma, eta);
2903}

◆ IsGas()

bool Garfield::MediumGas::IsGas ( ) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 19 of file MediumGas.hh.

19{ return true; }

◆ LoadGasFile()

bool Garfield::MediumGas::LoadGasFile ( const std::string &  filename)

Definition at line 283 of file MediumGas.cc.

283 {
284
285 std::ifstream gasfile;
286 // Open the file.
287 gasfile.open(filename.c_str());
288 // Make sure the file could be opened.
289 if (!gasfile.is_open()) {
290 std::cerr << m_className << "::LoadGasFile:\n";
291 std::cerr << " Gas file could not be opened.\n";
292 return false;
293 }
294
295 char line[256];
296 char* token;
297
298 // GASOK bits
299 std::string gasBits = "";
300
301 // Gas composition
302 const int nMagboltzGases = 60;
303 std::vector<double> mixture(nMagboltzGases);
304 for (int i = nMagboltzGases; i--;) mixture[i] = 0.;
305
306 int excCount = 0;
307 int ionCount = 0;
308
309 int eFieldRes = 1;
310 int bFieldRes = 1;
311 int angRes = 1;
312
313 int versionNumber = 12;
314
315 // Start reading the data.
316 bool atTables = false;
317 while (!atTables) {
318 gasfile.getline(line, 256);
319 if (strncmp(line, " The gas tables follow:", 8) == 0 ||
320 strncmp(line, "The gas tables follow:", 7) == 0) {
321 atTables = true;
322 if (m_debug) {
323 std::cout << " Entering tables.\n";
324 getchar();
325 }
326 }
327 if (m_debug) {
328 std::cout << " Line: " << line << "\n";
329 getchar();
330 }
331 if (!atTables) {
332 token = strtok(line, " :,%");
333 while (token != NULL) {
334 if (m_debug) std::cout << " Token: " << token << "\n";
335 if (strcmp(token, "Version") == 0) {
336 token = strtok(NULL, " :,%");
337 versionNumber = atoi(token);
338 // Check the version number.
339 if (versionNumber != 10 && versionNumber != 11 &&
340 versionNumber != 12) {
341 std::cerr << m_className << "::LoadGasFile:\n";
342 std::cerr << " The file has version number " << versionNumber
343 << ".\n";
344 std::cerr << " Files written in this format cannot be read.\n";
345 gasfile.close();
346 return false;
347 } else {
348 std::cout << m_className << "::LoadGasFile:\n";
349 std::cout << " Version: " << versionNumber << "\n";
350 }
351 } else if (strcmp(token, "GASOK") == 0) {
352 // Get the GASOK bits indicating if a parameter
353 // is present in the table (T) or not (F).
354 token = strtok(NULL, " :,%\t");
355 token = strtok(NULL, " :,%\t");
356 gasBits += token;
357 } else if (strcmp(token, "Identifier") == 0) {
358 // Get the identification string.
359 std::string identifier = "";
360 token = strtok(NULL, "\n");
361 if (token != NULL) identifier += token;
362 if (m_debug) {
363 std::cout << m_className << "::LoadGasFile:\n";
364 std::cout << " Identifier:\n";
365 std::cout << " " << token << "\n";
366 }
367 } else if (strcmp(token, "Dimension") == 0) {
368 token = strtok(NULL, " :,%\t");
369 if (strcmp(token, "F") == 0) {
370 m_map2d = false;
371 } else {
372 m_map2d = true;
373 }
374 token = strtok(NULL, " :,%\t");
375 eFieldRes = atoi(token);
376 // Check the number of E points.
377 if (eFieldRes <= 0) {
378 std::cerr << m_className << "::LoadGasFile:\n";
379 std::cerr << " Number of E fields out of range.\n";
380 gasfile.close();
381 return false;
382 }
383 token = strtok(NULL, " :,%\t");
384 angRes = atoi(token);
385 // Check the number of angles.
386 if (m_map2d && angRes <= 0) {
387 std::cerr << m_className << "::LoadGasFile:\n";
388 std::cerr << " Number of E-B angles out of range.\n";
389 gasfile.close();
390 return false;
391 }
392
393 token = strtok(NULL, " :,%\t");
394 bFieldRes = atoi(token);
395 // Check the number of B points.
396 if (m_map2d && bFieldRes <= 0) {
397 std::cerr << m_className << "::LoadGasFile:\n";
398 std::cerr << " Number of B fields out of range.\n";
399 gasfile.close();
400 return false;
401 }
402
403 eFields.resize(eFieldRes);
404 m_nEfields = eFieldRes;
405 bFields.resize(bFieldRes);
406 m_nBfields = bFieldRes;
407 bAngles.resize(angRes);
408 m_nAngles = angRes;
409
410 // Fill in the excitation/ionisation structs
411 // Excitation
412 token = strtok(NULL, " :,%\t");
413 if (m_debug) {
414 std::cout << " " << token << "\n";
415 }
416 int nexc = atoi(token);
417 if (nexc >= 0) nExcListElements = nexc;
418 if (nExcListElements > 0) {
420 }
421 // Ionization
422 token = strtok(NULL, " :,%\t");
423 int nion = atoi(token);
424 if (nion >= 0) nIonListElements = nion;
425 if (nIonListElements > 0) {
427 }
428 if (m_debug) {
429 std::cout << " Finished initializing excitation/ionisation"
430 << " structs.\n";
431 }
432 } else if (strcmp(token, "E") == 0) {
433 token = strtok(NULL, " :,%");
434 if (strcmp(token, "fields") == 0) {
435 for (int i = 0; i < eFieldRes; i++) {
436 gasfile >> eFields[i];
437 }
438 }
439 } else if (strcmp(token, "E-B") == 0) {
440 token = strtok(NULL, " :,%");
441 if (strcmp(token, "angles") == 0) {
442 for (int i = 0; i < angRes; i++) {
443 gasfile >> bAngles[i];
444 }
445 }
446 } else if (strcmp(token, "B") == 0) {
447 token = strtok(NULL, " :,%");
448 if (strcmp(token, "fields") == 0) {
449 double bstore = 0.;
450 for (int i = 0; i < bFieldRes; i++) {
451 // B fields are stored in hGauss (to be checked!).
452 gasfile >> bstore;
453 bFields[i] = bstore / 100.;
454 }
455 }
456 } else if (strcmp(token, "Mixture") == 0) {
457 for (int i = 0; i < nMagboltzGases; ++i) {
458 gasfile >> mixture[i];
459 }
460 } else if (strcmp(token, "Excitation") == 0) {
461 // Skip number.
462 token = strtok(NULL, " :,%");
463 // Get label.
464 token = strtok(NULL, " :,%");
465 excitationList[excCount].label += token;
466 // Get energy.
467 token = strtok(NULL, " :,%");
468 excitationList[excCount].energy = atof(token);
469 // Get Penning probability.
470 token = strtok(NULL, " :,%");
471 excitationList[excCount].prob = atof(token);
472 if (versionNumber >= 11) {
473 // Get Penning rms distance.
474 token = strtok(NULL, " :,%");
475 excitationList[excCount].rms = atof(token);
476 // Get decay time.
477 token = strtok(NULL, " :,%");
478 excitationList[excCount].dt = atof(token);
479 } else {
480 excitationList[excCount].rms = 0.;
481 excitationList[excCount].dt = 0.;
482 }
483 // Increase counter.
484 excCount++;
485 } else if (strcmp(token, "Ionisation") == 0) {
486 // Skip number.
487 token = strtok(NULL, " :,%");
488 // Get label.
489 token = strtok(NULL, " :,%");
490 ionisationList[ionCount].label += token;
491 // Get energy.
492 token = strtok(NULL, " :,%");
493 ionisationList[ionCount].energy = atof(token);
494 // Increase counter.
495 ionCount++;
496 }
497 token = strtok(NULL, " :,%");
498 }
499 }
500 }
501
502 // Decode the GASOK bits.
503 // GASOK(I) : .TRUE. if present
504 // (1) electron drift velocity || E
505 // (2) ion mobility,
506 // (3) longitudinal diffusion || E
507 // (4) Townsend coefficient,
508 // (5) cluster size distribution.
509 // (6) attachment coefficient,
510 // (7) Lorentz angle,
511 // (8) transverse diffusion || ExB and Bt
512 // (9) electron drift velocity || Bt
513 // (10) electron drift velocity || ExB
514 // (11) diffusion tensor
515 // (12) ion dissociation
516 // (13) allocated for SRIM data (not used)
517 // (14) allocated for HEED data (not used)
518 // (15) excitation rates
519 // (16) ionisation rates
520
521 if (m_debug) {
522 std::cout << m_className << "::LoadGasFile:\n";
523 std::cout << " GASOK bits: " << gasBits << "\n";
524 }
525
526 if (gasBits[0] == 'T') {
528 InitParamArrays(eFieldRes, bFieldRes, angRes, tabElectronVelocityE, 0.);
529 } else {
531 tabElectronVelocityE.clear();
532 }
533 if (gasBits[1] == 'T') {
534 m_hasIonMobility = true;
535 InitParamArrays(eFieldRes, bFieldRes, angRes, tabIonMobility, 0.);
536 } else {
537 m_hasIonMobility = false;
538 tabIonMobility.clear();
539 }
540 if (gasBits[2] == 'T') {
542 InitParamArrays(eFieldRes, bFieldRes, angRes, tabElectronDiffLong, 0.);
543 } else {
544 m_hasElectronDiffLong = false;
545 tabElectronDiffLong.clear();
546 }
547 if (gasBits[3] == 'T') {
549 InitParamArrays(eFieldRes, bFieldRes, angRes, tabElectronTownsend, -30.);
550 InitParamArrays(eFieldRes, bFieldRes, angRes, tabTownsendNoPenning, -30.);
551 } else {
552 m_hasElectronTownsend = false;
553 tabElectronTownsend.clear();
554 tabTownsendNoPenning.clear();
555 }
556 // gasBits[4]: cluster size distribution; skipped
557 if (gasBits[5] == 'T') {
559 InitParamArrays(eFieldRes, bFieldRes, angRes, tabElectronAttachment, -30.);
560 } else {
562 tabElectronAttachment.clear();
563 }
564 // gasBits[6]: Lorentz angle; skipped
565 if (gasBits[7] == 'T') {
567 InitParamArrays(eFieldRes, bFieldRes, angRes, tabElectronDiffTrans, 0.);
568 } else {
570 tabElectronDiffTrans.clear();
571 }
572 if (gasBits[8] == 'T') {
574 InitParamArrays(eFieldRes, bFieldRes, angRes, tabElectronVelocityB, 0.);
575 } else {
577 tabElectronVelocityB.clear();
578 }
579 if (gasBits[9] == 'T') {
581 InitParamArrays(eFieldRes, bFieldRes, angRes, tabElectronVelocityExB, 0.);
582 } else {
585 }
586 if (gasBits[10] == 'T') {
588 InitParamTensor(eFieldRes, bFieldRes, angRes, 6, tabElectronDiffTens, 0.);
589 } else {
590 m_hasElectronDiffTens = false;
591 tabElectronDiffTens.clear();
592 }
593 if (gasBits[11] == 'T') {
595 InitParamArrays(eFieldRes, bFieldRes, angRes, tabIonDissociation, -30.);
596 } else {
597 m_hasIonDissociation = false;
598 tabIonDissociation.clear();
599 }
600 // gasBits[12]: SRIM; skipped
601 // gasBits[13]: HEED; skipped
602 if (gasBits[14] == 'T') {
603 m_hasExcRates = true;
604 InitParamTensor(eFieldRes, bFieldRes, angRes, nExcListElements, tabExcRates,
605 0.);
606 } else {
607 m_hasExcRates = false;
608 tabExcRates.clear();
609 }
610 if (gasBits[15] == 'T') {
611 m_hasIonRates = true;
612 InitParamTensor(eFieldRes, bFieldRes, angRes, nIonListElements, tabIonRates,
613 0.);
614 } else {
615 m_hasIonRates = false;
616 tabIonRates.clear();
617 }
618
619 // Check the gas mixture.
620 std::vector<std::string> gasnames;
621 gasnames.clear();
622 std::vector<double> percentages;
623 percentages.clear();
624 bool gasMixOk = true;
625 unsigned int gasCount = 0;
626 for (int i = 0; i < nMagboltzGases; ++i) {
627 if (mixture[i] > 0.) {
628 std::string gasname = "";
629 if (!GetGasName(i + 1, versionNumber, gasname)) {
630 std::cerr << m_className << "::LoadGasFile:\n";
631 std::cerr << " Unknown gas (gas number ";
632 std::cerr << i + 1 << ")\n";
633 gasMixOk = false;
634 break;
635 }
636 gasnames.push_back(gasname);
637 percentages.push_back(mixture[i]);
638 ++gasCount;
639 }
640 }
641 if (gasCount > m_nMaxGases) {
642 std::cerr << m_className << "::LoadGasFile:\n";
643 std::cerr << " Gas mixture has " << gasCount << " components.\n";
644 std::cerr << " Number of gases is limited to " << m_nMaxGases << ".\n";
645 gasMixOk = false;
646 } else if (gasCount == 0) {
647 std::cerr << m_className << "::LoadGasFile:\n";
648 std::cerr << " Gas mixture is not defined (zero components).\n";
649 gasMixOk = false;
650 }
651 double sum = 0.;
652 for (unsigned int i = 0; i < gasCount; ++i) sum += percentages[i];
653 if (gasMixOk && sum != 100.) {
654 std::cerr << m_className << "::LoadGasFile:\n";
655 std::cerr << " Percentages are not normalized to 100.\n";
656 for (unsigned int i = 0; i < gasCount; ++i) percentages[i] *= 100. / sum;
657 }
658
659 // Force re-initialisation of collision rates etc.
660 m_isChanged = true;
661
662 if (gasMixOk) {
663 m_name = "";
664 m_nComponents = gasCount;
665 for (unsigned int i = 0; i < m_nComponents; ++i) {
666 if (i > 0) m_name += "/";
667 m_name += gasnames[i];
668 gas[i] = gasnames[i];
669 fraction[i] = percentages[i] / 100.;
670 GetGasInfo(gas[i], atWeight[i], atNum[i]);
671 }
672 std::cout << m_className << "::LoadGasFile:\n";
673 std::cout << " Gas composition set to " << m_name;
674 if (m_nComponents > 1) {
675 std::cout << " (" << fraction[0] * 100;
676 for (unsigned int i = 1; i < m_nComponents; ++i) {
677 std::cout << "/" << fraction[i] * 100;
678 }
679 std::cout << ")";
680 }
681 std::cout << "\n";
682 } else {
683 std::cerr << m_className << "::LoadGasFile:\n";
684 std::cerr << " Gas composition could not be established.\n";
685 }
686
687 if (m_debug) {
688 std::cout << m_className << "::LoadGasFile:\n";
689 std::cout << " Reading gas tables.\n";
690 }
691
692 // Temporary variables
693 // Velocities
694 double ve = 0., vb = 0., vexb = 0.;
695 // Lorentz angle
696 double lor = 0.;
697 // Diffusion coefficients
698 double dl = 0., dt = 0.;
699 // Towsend and attachment coefficients
700 double alpha = 0., alpha0 = 0., eta = 0.;
701 // Ion mobility and dissociation coefficient
702 double mu = 0., diss = 0.;
703 double ionDiffLong = 0., ionDiffTrans = 0.;
704 double diff = 0.;
705 double rate = 0.;
706 double waste = 0.;
707
708 if (m_map2d) {
709 if (m_debug) {
710 std::cout << m_className << "::LoadGasFile:\n";
711 std::cout << " Gas table is 3D.\n";
712 }
713 for (int i = 0; i < eFieldRes; i++) {
714 for (int j = 0; j < angRes; j++) {
715 for (int k = 0; k < bFieldRes; k++) {
716 // Drift velocity along E, Bt and ExB
717 gasfile >> ve >> vb >> vexb;
718 // Convert from cm / us to cm / ns
719 ve *= 1.e-3;
720 vb *= 1.e-3;
721 vexb *= 1.e-3;
725 // Longitudinal and transverse diffusion coefficient
726 gasfile >> dl >> dt;
729 // Townsend and attachment coefficient
730 gasfile >> alpha >> alpha0 >> eta;
732 tabElectronTownsend[j][k][i] = alpha;
733 tabTownsendNoPenning[j][k][i] = alpha0;
734 }
736 tabElectronAttachment[j][k][i] = eta;
737 }
738 // Ion mobility
739 gasfile >> mu;
740 // Convert from cm2 / (V us) to cm2 / (V ns)
741 mu *= 1.e-3;
742 if (m_hasIonMobility) tabIonMobility[j][k][i] = mu;
743 // Lorentz angle (unused)
744 gasfile >> lor;
745 // Ion dissociation
746 gasfile >> diss;
747 if (m_hasIonDissociation) tabIonDissociation[j][k][i] = diss;
748 // Diffusion tensor
749 for (int l = 0; l < 6; l++) {
750 gasfile >> diff;
751 if (m_hasElectronDiffTens) tabElectronDiffTens[l][j][k][i] = diff;
752 }
753 // Excitation rates
754 for (int l = 0; l < nExcListElements; l++) {
755 gasfile >> rate;
756 if (m_hasExcRates) tabExcRates[l][j][k][i] = rate;
757 }
758 // Ionization rates
759 for (int l = 0; l < nIonListElements; l++) {
760 gasfile >> rate;
761 if (m_hasIonRates) tabIonRates[l][j][k][i] = rate;
762 }
763 }
764 }
765 }
766 } else {
767 if (m_debug) {
768 std::cout << m_className << "::LoadGasFile:\n";
769 std::cout << " Gas table is 1D.\n";
770 }
771 for (int i = 0; i < eFieldRes; i++) {
772 if (m_debug) std::cout << " Done table: " << i << "\n";
773 // Drift velocity along E, Bt, ExB
774 gasfile >> ve >> waste >> vb >> waste >> vexb >> waste;
775 ve *= 1.e-3;
776 vb *= 1.e-3;
777 vexb *= 1.e-3;
781 // Longitudinal and transferse diffusion coefficients
782 gasfile >> dl >> waste >> dt >> waste;
785 // Townsend and attachment coefficients
786 gasfile >> alpha >> waste >> alpha0 >> eta >> waste;
788 tabElectronTownsend[0][0][i] = alpha;
789 tabTownsendNoPenning[0][0][i] = alpha0;
790 }
792 tabElectronAttachment[0][0][i] = eta;
793 }
794 // Ion mobility
795 gasfile >> mu >> waste;
796 mu *= 1.e-3;
797 if (m_hasIonMobility) tabIonMobility[0][0][i] = mu;
798 // Lorentz angle (unused)
799 gasfile >> lor >> waste;
800 // Ion dissociation
801 gasfile >> diss >> waste;
802 if (m_hasIonDissociation) tabIonDissociation[0][0][i] = diss;
803 // Diffusion tensor
804 for (int j = 0; j < 6; j++) {
805 gasfile >> diff >> waste;
806 if (m_hasElectronDiffTens) tabElectronDiffTens[j][0][0][i] = diff;
807 }
808 // Excitation rates
809 for (int j = 0; j < nExcListElements; j++) {
810 gasfile >> rate >> waste;
811 if (m_hasExcRates) tabExcRates[j][0][0][i] = rate;
812 }
813 // Ionization rates
814 for (int j = 0; j < nIonListElements; j++) {
815 gasfile >> rate >> waste;
816 if (m_hasIonRates) tabIonRates[j][0][0][i] = rate;
817 }
818 }
819 }
820 if (m_debug) std::cout << " Done with gas tables.\n";
821
822 // Extrapolation methods
823 int hExtrap[13], lExtrap[13];
824 // Interpolation methods
825 int interpMeth[13];
826
827 // Moving on to the file footer
828 bool done = false;
829 while (!done) {
830 gasfile.getline(line, 256);
831 token = strtok(line, " :,%=\t");
832 while (token != NULL) {
833 if (strcmp(token, "H") == 0) {
834 token = strtok(NULL, " :,%=\t");
835 for (int i = 0; i < 13; i++) {
836 token = strtok(NULL, " :,%=\t");
837 if (token != NULL) hExtrap[i] = atoi(token);
838 }
839 } else if (strcmp(token, "L") == 0) {
840 token = strtok(NULL, " :,%=\t");
841 for (int i = 0; i < 13; i++) {
842 token = strtok(NULL, " :,%=\t");
843 if (token != NULL) lExtrap[i] = atoi(token);
844 }
845 } else if (strcmp(token, "Thresholds") == 0) {
846 token = strtok(NULL, " :,%=\t");
847 if (token != NULL) thrElectronTownsend = atoi(token);
848 token = strtok(NULL, " :,%=\t");
849 if (token != NULL) thrElectronAttachment = atoi(token);
850 token = strtok(NULL, " :,%=\t");
851 if (token != NULL) thrIonDissociation = atoi(token);
852 } else if (strcmp(token, "Interp") == 0) {
853 for (int i = 0; i < 13; i++) {
854 token = strtok(NULL, " :,%=\t");
855 if (token != NULL) interpMeth[i] = atoi(token);
856 }
857 } else if (strcmp(token, "A") == 0) {
858 token = strtok(NULL, " :,%=\t");
859 // Parameter for energy loss distribution, currently not used
860 // double a;
861 // if (token != NULL) a = atof(token);
862 } else if (strcmp(token, "Z") == 0) {
863 // Parameter for energy loss distribution, currently not used
864 token = strtok(NULL, " :,%=\t");
865 // double z;
866 // if (token != NULL) z = atof(token);
867 } else if (strcmp(token, "EMPROB") == 0) {
868 // Parameter for energy loss distribution, currently not used
869 token = strtok(NULL, " :,%=\t");
870 // double emprob;
871 // if (token != NULL) emprob = atof(token);
872 } else if (strcmp(token, "EPAIR") == 0) {
873 // Parameter for energy loss distribution, currently not used
874 token = strtok(NULL, " :,%=\t");
875 // double epair;
876 // if (token != NULL) epair = atof(token);
877 } else if (strcmp(token, "Ion") == 0) {
878 // Ion diffusion coefficients
879 token = strtok(NULL, " :,%=\t");
880 token = strtok(NULL, " :,%=\t");
881 if (token != NULL) ionDiffLong = atof(token);
882 token = strtok(NULL, " :,%=\t");
883 if (token != NULL) ionDiffTrans = atof(token);
884 } else if (strcmp(token, "CMEAN") == 0) {
885 // Cluster parameter, currently not used
886 token = strtok(NULL, " :,%=\t");
887 // double clsPerCm;
888 // if (token != NULL) clsPerCm = atof(token);
889 } else if (strcmp(token, "RHO") == 0) {
890 // Parameter for energy loss distribution, currently not used
891 token = strtok(NULL, " :,%=\t");
892 // double rho;
893 // if (token != NULL) rho = atof(token);
894 } else if (strcmp(token, "PGAS") == 0) {
895 token = strtok(NULL, " :,%=\t");
896 double pTorr = 760.;
897 if (token != NULL) pTorr = atof(token);
898 if (pTorr > 0.) m_pressure = pTorr;
899 } else if (strcmp(token, "TGAS") == 0) {
900 token = strtok(NULL, " :,%=\t");
901 double tKelvin = 293.15;
902 if (token != NULL) tKelvin = atof(token);
903 if (tKelvin > 0.) m_temperature = tKelvin;
904 done = true;
905 break;
906 } else {
907 done = true;
908 break;
909 }
910 token = strtok(NULL, " :,%=\t");
911 }
912 }
913
914 gasfile.close();
915
916 // Set the reference pressure and temperature.
919
920 // Multiply the E/p values by the pressure.
921 for (int i = eFieldRes; i--;) {
923 }
924 // Scale the parameters.
925 const double sqrtPressure = sqrt(pressureTable);
926 const double logPressure = log(pressureTable);
927 for (int i = eFieldRes; i--;) {
928 for (int j = angRes; j--;) {
929 for (int k = bFieldRes; k--;) {
931 tabElectronDiffLong[j][k][i] /= sqrtPressure;
932 }
934 tabElectronDiffTrans[j][k][i] /= sqrtPressure;
935 }
937 for (int l = 6; l--;) {
938 tabElectronDiffTens[l][j][k][i] /= pressureTable;
939 }
940 }
942 tabElectronTownsend[j][k][i] += logPressure;
943 }
945 tabElectronAttachment[j][k][i] += logPressure;
946 }
948 tabIonDissociation[j][k][i] += logPressure;
949 }
950 }
951 }
952 }
953
954 // Decode the extrapolation and interpolation tables.
955 m_extrHighVelocity = hExtrap[0];
956 m_extrLowVelocity = lExtrap[0];
957 m_intpVelocity = interpMeth[0];
958 // Indices 1 and 2 correspond to velocities along Bt and ExB.
959 m_extrHighDiffusion = hExtrap[3];
960 m_extrLowDiffusion = lExtrap[3];
961 m_intpDiffusion = interpMeth[3];
962 m_extrHighTownsend = hExtrap[4];
963 m_extrLowTownsend = lExtrap[4];
964 m_intpTownsend = interpMeth[4];
965 m_extrHighAttachment = hExtrap[5];
966 m_extrLowAttachment = lExtrap[5];
967 m_intpAttachment = interpMeth[5];
968 m_extrHighMobility = hExtrap[6];
969 m_extrLowMobility = lExtrap[6];
970 m_intpMobility = interpMeth[6];
971 // Index 7, 8: Lorentz angle, transv. diff.
972 m_extrHighDissociation = hExtrap[9];
973 m_extrLowDissociation = lExtrap[9];
974 m_intpDissociation = interpMeth[9];
975 // Index 10: diff. tensor
976 m_extrHighExcRates = hExtrap[11];
977 m_extrLowExcRates = lExtrap[11];
978 m_intpExcRates = interpMeth[11];
979 m_extrHighIonRates = hExtrap[12];
980 m_extrLowIonRates = lExtrap[12];
981 m_intpIonRates = interpMeth[12];
982
983 // Ion diffusion
984 if (ionDiffLong > 0.) {
985 m_hasIonDiffLong = true;
986 InitParamArrays(eFieldRes, bFieldRes, angRes, tabIonDiffLong, 0.);
987 for (int i = eFieldRes; i--;) {
988 for (int j = angRes; j--;) {
989 for (int k = bFieldRes; k--;) {
990 tabIonDiffLong[j][k][i] = ionDiffLong;
991 }
992 }
993 }
994 } else {
995 m_hasIonDiffLong = false;
996 tabIonDiffLong.clear();
997 }
998 if (ionDiffTrans > 0.) {
999 m_hasIonDiffTrans = true;
1000 InitParamArrays(eFieldRes, bFieldRes, angRes, tabIonDiffTrans, 0.);
1001 for (int i = eFieldRes; i--;) {
1002 for (int j = angRes; j--;) {
1003 for (int k = bFieldRes; k--;) {
1004 tabIonDiffTrans[j][k][i] = ionDiffTrans;
1005 }
1006 }
1007 }
1008 } else {
1009 m_hasIonDiffTrans = false;
1010 tabIonDiffTrans.clear();
1011 }
1012
1013 if (m_debug) {
1014 std::cout << m_className << "::LoadGasFile:\n";
1015 std::cout << " Gas file sucessfully read.\n";
1016 }
1017
1018 return true;
1019}
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:313
std::vector< std::vector< std::vector< std::vector< double > > > > tabIonRates
Definition: MediumGas.hh:110
bool GetGasName(const int gasnumber, const int version, std::string &gasname)
Definition: MediumGas.cc:2045
std::vector< std::vector< std::vector< double > > > tabTownsendNoPenning
Definition: MediumGas.hh:105
std::vector< std::vector< std::vector< std::vector< double > > > > tabExcRates
Definition: MediumGas.hh:109
bool m_hasElectronDiffTens
Definition: Medium.hh:334
unsigned int m_nBfields
Definition: Medium.hh:323
unsigned int m_extrLowMobility
Definition: Medium.hh:383
unsigned int m_intpDiffusion
Definition: Medium.hh:388
unsigned int m_extrLowDiffusion
Definition: Medium.hh:380
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)
Definition: Medium.cc:2792
bool m_hasElectronDiffTrans
Definition: Medium.hh:334
bool m_hasElectronTownsend
Definition: Medium.hh:335
std::vector< double > bFields
Definition: Medium.hh:327
std::vector< std::vector< std::vector< double > > > tabElectronVelocityB
Definition: Medium.hh:338
unsigned int m_extrHighAttachment
Definition: Medium.hh:382
unsigned int m_nAngles
Definition: Medium.hh:324
unsigned int m_intpMobility
Definition: Medium.hh:391
unsigned int m_intpAttachment
Definition: Medium.hh:390
bool m_hasElectronVelocityB
Definition: Medium.hh:333
unsigned int m_extrHighDiffusion
Definition: Medium.hh:380
std::vector< double > eFields
Definition: Medium.hh:326
std::vector< std::vector< std::vector< double > > > tabIonDissociation
Definition: Medium.hh:368
unsigned int m_extrHighMobility
Definition: Medium.hh:383
std::vector< std::vector< std::vector< double > > > tabElectronVelocityE
Definition: Medium.hh:336
unsigned int m_intpTownsend
Definition: Medium.hh:389
unsigned int m_extrLowAttachment
Definition: Medium.hh:382
std::vector< double > bAngles
Definition: Medium.hh:328
unsigned int m_extrHighDissociation
Definition: Medium.hh:384
unsigned int m_intpVelocity
Definition: Medium.hh:387
std::vector< std::vector< std::vector< double > > > tabElectronDiffLong
Definition: Medium.hh:339
std::vector< std::vector< std::vector< double > > > tabElectronAttachment
Definition: Medium.hh:342
int thrIonDissociation
Definition: Medium.hh:376
unsigned int m_extrLowVelocity
Definition: Medium.hh:379
unsigned int m_extrHighTownsend
Definition: Medium.hh:381
unsigned int m_extrHighVelocity
Definition: Medium.hh:379
unsigned int m_extrLowTownsend
Definition: Medium.hh:381
std::vector< std::vector< std::vector< double > > > tabIonDiffTrans
Definition: Medium.hh:367
unsigned int m_extrLowDissociation
Definition: Medium.hh:384
std::vector< std::vector< std::vector< std::vector< double > > > > tabElectronDiffTens
Definition: Medium.hh:345
unsigned int m_nEfields
Definition: Medium.hh:322
bool m_hasElectronVelocityExB
Definition: Medium.hh:333
std::vector< std::vector< std::vector< double > > > tabIonMobility
Definition: Medium.hh:365
bool m_hasIonDiffTrans
Definition: Medium.hh:363
int thrElectronAttachment
Definition: Medium.hh:372
bool m_hasIonMobility
Definition: Medium.hh:362
bool m_hasElectronDiffLong
Definition: Medium.hh:334
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:2772
std::vector< std::vector< std::vector< double > > > tabIonDiffLong
Definition: Medium.hh:366
bool m_hasIonDissociation
Definition: Medium.hh:364
bool m_hasIonDiffLong
Definition: Medium.hh:363
unsigned int m_intpDissociation
Definition: Medium.hh:392
bool m_hasElectronVelocityE
Definition: Medium.hh:333
int thrElectronTownsend
Definition: Medium.hh:371
bool m_hasElectronAttachment
Definition: Medium.hh:335
std::vector< std::vector< std::vector< double > > > tabElectronDiffTrans
Definition: Medium.hh:340
std::vector< std::vector< std::vector< double > > > tabElectronVelocityExB
Definition: Medium.hh:337

Referenced by GarfieldPhysics::InitializePhysics().

◆ LoadIonMobility()

bool Garfield::MediumGas::LoadIonMobility ( const std::string &  filename)

Definition at line 1717 of file MediumGas.cc.

1717 {
1718
1719 // Open the file.
1720 std::ifstream infile;
1721 infile.open(filename.c_str(), std::ios::in);
1722 // Make sure the file could actually be opened.
1723 if (!infile) {
1724 std::cerr << m_className << "::LoadIonMobility:\n";
1725 std::cerr << " Error opening file\n";
1726 std::cerr << " " << filename << ".\n";
1727 return false;
1728 }
1729
1730 double field = -1., mu = -1.;
1731 double lastField = field;
1732 std::vector<double> efields;
1733 std::vector<double> mobilities;
1734 efields.clear();
1735 mobilities.clear();
1736
1737 // Read the file line by line.
1738 char line[100];
1739
1740 int i = 0;
1741 while (!infile.eof()) {
1742 ++i;
1743 // Read the next line.
1744 infile.getline(line, 100);
1745
1746 char* token = NULL;
1747 token = strtok(line, " ,\t");
1748 if (!token) {
1749 break;
1750 } else if (strcmp(token, "#") == 0 || strcmp(token, "*") == 0 ||
1751 strcmp(token, "//") == 0) {
1752 continue;
1753 } else {
1754 field = atof(token);
1755 token = strtok(NULL, " ,\t");
1756 if (!token) {
1757 std::cerr << m_className << "::LoadIonMobility:\n";
1758 std::cerr << " Found E/N but no mobility before the end-of-line\n";
1759 std::cerr << " " << filename << " (line " << i << ").\n";
1760 continue;
1761 }
1762 mu = atof(token);
1763 }
1764 token = strtok(NULL, " ,\t");
1765 if (token && strcmp(token, "//") != 0) {
1766 std::cerr << m_className << "::LoadIonMobility:\n";
1767 std::cerr << " Unexpected non-comment characters after the mobility; "
1768 << " line skipped.\n";
1769 std::cerr << " " << filename << " (line " << i << ").\n";
1770 continue;
1771 }
1772 if (m_debug) {
1773 std::cout << " E/N = " << field << " Td: mu = " << mu << " cm2/(V.s)\n";
1774 }
1775 // Check if the data has been read correctly.
1776 if (infile.fail() && !infile.eof()) {
1777 std::cerr << m_className << "::LoadIonMobility:\n";
1778 std::cerr << " Error reading file\n";
1779 std::cerr << " " << filename << " (line " << i << ").\n";
1780 return false;
1781 }
1782 // Make sure the values make sense.
1783 // Negative field values are not allowed.
1784 if (field < 0.) {
1785 std::cerr << m_className << "::LoadIonMobility:\n";
1786 std::cerr << " Negative electric field (line " << i << ").\n";
1787 return false;
1788 }
1789 // The table has to be in ascending order.
1790 if (field <= lastField) {
1791 std::cerr << m_className << "::LoadIonMobility:\n";
1792 std::cerr << " Table is not in ascending order (line " << i << ").\n";
1793 return false;
1794 }
1795 // Add the values to the list.
1796 efields.push_back(field);
1797 mobilities.push_back(mu);
1798 lastField = field;
1799 }
1800
1801 const int ne = efields.size();
1802 if (ne <= 0) {
1803 std::cerr << m_className << "::LoadIonMobilities:\n";
1804 std::cerr << " No valid data found.\n";
1805 return false;
1806 }
1807
1808 // The E/N values in the file are supposed to be in Td (10^-17 V cm2).
1809 const double scaleField = 1.e-17 * GetNumberDensity();
1810 // The reduced mobilities in the file are supposed to be in cm2/(V s).
1811 const double scaleMobility =
1812 1.e-9 * (AtmosphericPressure / m_pressure) * (m_temperature / ZeroCelsius);
1813 for (int j = ne; j--;) {
1814 // Scale the fields and mobilities.
1815 efields[j] *= scaleField;
1816 mobilities[j] *= scaleMobility;
1817 }
1818
1819 std::cout << m_className << "::LoadIonMobility:\n";
1820 std::cout << " Read " << ne << " values from file " << filename << "\n";
1821
1822 return SetIonMobility(efields, mobilities);
1823}
bool SetIonMobility(const unsigned int &ie, const unsigned int &ib, const unsigned int &ia, const double &mu)
Definition: Medium.cc:2454

◆ PrintGas()

void Garfield::MediumGas::PrintGas ( )

Definition at line 1437 of file MediumGas.cc.

1437 {
1438
1439 // Print a summary.
1440 std::cout << m_className << "::PrintGas:\n";
1441 std::cout << " Gas composition: " << m_name;
1442 if (m_nComponents > 1) {
1443 std::cout << " (" << fraction[0] * 100;
1444 for (unsigned int i = 1; i < m_nComponents; ++i) {
1445 std::cout << "/" << fraction[i] * 100;
1446 }
1447 std::cout << ")";
1448 }
1449 std::cout << "\n";
1450 std::cout << " Pressure: " << m_pressure << " Torr\n";
1451 std::cout << " Temperature: " << m_temperature << " K\n";
1452 std::cout << " Gas file:\n";
1453 std::cout << " Pressure: " << pressureTable << " Torr\n";
1454 std::cout << " Temperature: " << temperatureTable << " K\n";
1455 if (m_nEfields > 1) {
1456 std::cout << " Electric field range: " << eFields[0] << " - "
1457 << eFields[m_nEfields - 1] << " V/cm in " << m_nEfields - 1
1458 << " steps.\n";
1459 } else if (m_nEfields == 1) {
1460 std::cout << " Electric field: " << eFields[0] << " V/cm\n";
1461 } else {
1462 std::cout << " Electric field range: not set\n";
1463 }
1464 if (m_nBfields > 1) {
1465 std::cout << " Magnetic field range: " << bFields[0] << " - "
1466 << bFields[m_nBfields - 1] << " T in " << m_nBfields - 1
1467 << " steps.\n";
1468 } else if (m_nBfields == 1) {
1469 std::cout << " Magnetic field: " << bFields[0] << "\n";
1470 } else {
1471 std::cout << " Magnetic field range: not set\n";
1472 }
1473 if (m_nAngles > 1) {
1474 std::cout << " Angular range: " << bAngles[0] << " - "
1475 << bAngles[m_nAngles - 1] << " in " << m_nAngles - 1
1476 << " steps.\n";
1477 } else if (m_nAngles == 1) {
1478 std::cout << " Angle between E and B: " << bAngles[0] << "\n";
1479 } else {
1480 std::cout << " Angular range: not set\n";
1481 }
1482
1483 std::cout << " Available electron transport data:\n";
1485 std::cout << " Velocity along E\n";
1486 }
1488 std::cout << " Velocity along Bt\n";
1489 }
1491 std::cout << " Velocity along ExB\n";
1492 }
1494 std::cout << " Low field extrapolation: ";
1495 if (m_extrLowVelocity == 0)
1496 std::cout << " constant\n";
1497 else if (m_extrLowVelocity == 1)
1498 std::cout << " linear\n";
1499 else if (m_extrLowVelocity == 2)
1500 std::cout << " exponential\n";
1501 else
1502 std::cout << " unknown\n";
1503 std::cout << " High field extrapolation: ";
1504 if (m_extrHighVelocity == 0)
1505 std::cout << " constant\n";
1506 else if (m_extrHighVelocity == 1)
1507 std::cout << " linear\n";
1508 else if (m_extrHighVelocity == 2)
1509 std::cout << " exponential\n";
1510 else
1511 std::cout << " unknown\n";
1512 std::cout << " Interpolation order: " << m_intpVelocity << "\n";
1513 }
1515 std::cout << " Longitudinal diffusion coefficient\n";
1516 }
1518 std::cout << " Transverse diffusion coefficient\n";
1519 }
1521 std::cout << " Diffusion tensor\n";
1522 }
1524 std::cout << " Low field extrapolation: ";
1525 if (m_extrLowDiffusion == 0)
1526 std::cout << " constant\n";
1527 else if (m_extrLowDiffusion == 1)
1528 std::cout << " linear\n";
1529 else if (m_extrLowDiffusion == 2)
1530 std::cout << " exponential\n";
1531 else
1532 std::cout << " unknown\n";
1533 std::cout << " High field extrapolation: ";
1534 if (m_extrHighDiffusion == 0)
1535 std::cout << " constant\n";
1536 else if (m_extrHighDiffusion == 1)
1537 std::cout << " linear\n";
1538 else if (m_extrHighDiffusion == 2)
1539 std::cout << " exponential\n";
1540 else
1541 std::cout << " unknown\n";
1542 std::cout << " Interpolation order: " << m_intpDiffusion << "\n";
1543 }
1545 std::cout << " Townsend coefficient\n";
1546 std::cout << " Low field extrapolation: ";
1547 if (m_extrLowTownsend == 0)
1548 std::cout << " constant\n";
1549 else if (m_extrLowTownsend == 1)
1550 std::cout << " linear\n";
1551 else if (m_extrLowTownsend == 2)
1552 std::cout << " exponential\n";
1553 else
1554 std::cout << " unknown\n";
1555 std::cout << " High field extrapolation: ";
1556 if (m_extrHighTownsend == 0)
1557 std::cout << " constant\n";
1558 else if (m_extrHighTownsend == 1)
1559 std::cout << " linear\n";
1560 else if (m_extrHighTownsend == 2)
1561 std::cout << " exponential\n";
1562 else
1563 std::cout << " unknown\n";
1564 std::cout << " Interpolation order: " << m_intpTownsend << "\n";
1565 }
1567 std::cout << " Attachment coefficient\n";
1568 std::cout << " Low field extrapolation: ";
1569 if (m_extrLowAttachment == 0)
1570 std::cout << " constant\n";
1571 else if (m_extrLowAttachment == 1)
1572 std::cout << " linear\n";
1573 else if (m_extrLowAttachment == 2)
1574 std::cout << " exponential\n";
1575 else
1576 std::cout << " unknown\n";
1577 std::cout << " High field extrapolation: ";
1578 if (m_extrHighAttachment == 0)
1579 std::cout << " constant\n";
1580 else if (m_extrHighAttachment == 1)
1581 std::cout << " linear\n";
1582 else if (m_extrHighAttachment == 2)
1583 std::cout << " exponential\n";
1584 else
1585 std::cout << " unknown\n";
1586 std::cout << " Interpolation order: " << m_intpAttachment << "\n";
1587 }
1588 if (m_hasExcRates) {
1589 std::cout << " Excitation rates\n";
1590 std::cout << " Low field extrapolation: ";
1591 if (m_extrLowExcRates == 0)
1592 std::cout << " constant\n";
1593 else if (m_extrLowExcRates == 1)
1594 std::cout << " linear\n";
1595 else if (m_extrLowExcRates == 2)
1596 std::cout << " exponential\n";
1597 else
1598 std::cout << " unknown\n";
1599 std::cout << " High field extrapolation: ";
1600 if (m_extrHighExcRates == 0)
1601 std::cout << " constant\n";
1602 else if (m_extrHighExcRates == 1)
1603 std::cout << " linear\n";
1604 else if (m_extrHighExcRates == 2)
1605 std::cout << " exponential\n";
1606 else
1607 std::cout << " unknown\n";
1608 std::cout << " Interpolation order: " << m_intpExcRates << "\n";
1609 }
1610 if (m_hasIonRates) {
1611 std::cout << " Ionisation rates\n";
1612 std::cout << " Low field extrapolation: ";
1613 if (m_extrLowIonRates == 0)
1614 std::cout << " constant\n";
1615 else if (m_extrLowIonRates == 1)
1616 std::cout << " linear\n";
1617 else if (m_extrLowIonRates == 2)
1618 std::cout << " exponential\n";
1619 else
1620 std::cout << " unknown\n";
1621 std::cout << " High field extrapolation: ";
1622 if (m_extrHighIonRates == 0)
1623 std::cout << " constant\n";
1624 else if (m_extrHighIonRates == 1)
1625 std::cout << " linear\n";
1626 else if (m_extrHighIonRates == 2)
1627 std::cout << " exponential\n";
1628 else
1629 std::cout << " unknown\n";
1630 std::cout << " Interpolation order: " << m_intpIonRates << "\n";
1631 }
1636 std::cout << " none\n";
1637 }
1638
1639 std::cout << " Available ion transport data:\n";
1640 if (m_hasIonMobility) {
1641 std::cout << " Mobility\n";
1642 std::cout << " Low field extrapolation: ";
1643 if (m_extrLowMobility == 0)
1644 std::cout << " constant\n";
1645 else if (m_extrLowMobility == 1)
1646 std::cout << " linear\n";
1647 else if (m_extrLowMobility == 2)
1648 std::cout << " exponential\n";
1649 else
1650 std::cout << " unknown\n";
1651 std::cout << " High field extrapolation: ";
1652 if (m_extrHighMobility == 0)
1653 std::cout << " constant\n";
1654 else if (m_extrHighMobility == 1)
1655 std::cout << " linear\n";
1656 else if (m_extrHighMobility == 2)
1657 std::cout << " exponential\n";
1658 else
1659 std::cout << " unknown\n";
1660 std::cout << " Interpolation order: " << m_intpMobility << "\n";
1661 }
1662 if (m_hasIonDiffLong) {
1663 std::cout << " Longitudinal diffusion coefficient\n";
1664 }
1665 if (m_hasIonDiffTrans) {
1666 std::cout << " Transverse diffusion coefficient\n";
1667 }
1669 std::cout << " Low field extrapolation: ";
1670 if (m_extrLowDiffusion == 0)
1671 std::cout << " constant\n";
1672 else if (m_extrLowDiffusion == 1)
1673 std::cout << " linear\n";
1674 else if (m_extrLowDiffusion == 2)
1675 std::cout << " exponential\n";
1676 else
1677 std::cout << " unknown\n";
1678 std::cout << " High field extrapolation: ";
1679 if (m_extrHighDiffusion == 0)
1680 std::cout << " constant\n";
1681 else if (m_extrHighDiffusion == 1)
1682 std::cout << " linear\n";
1683 else if (m_extrHighDiffusion == 2)
1684 std::cout << " exponential\n";
1685 else
1686 std::cout << " unknown\n";
1687 std::cout << " Interpolation order: " << m_intpDiffusion << "\n";
1688 }
1690 std::cout << " Dissociation coefficient\n";
1691 std::cout << " Low field extrapolation: ";
1692 if (m_extrLowDissociation == 0)
1693 std::cout << " constant\n";
1694 else if (m_extrLowDissociation == 1)
1695 std::cout << " linear\n";
1696 else if (m_extrLowDissociation == 2)
1697 std::cout << " exponential\n";
1698 else
1699 std::cout << " unknown\n";
1700 std::cout << " High field extrapolation: ";
1701 if (m_extrHighDissociation == 0)
1702 std::cout << " constant\n";
1703 else if (m_extrHighDissociation == 1)
1704 std::cout << " linear\n";
1705 else if (m_extrHighDissociation == 2)
1706 std::cout << " exponential\n";
1707 else
1708 std::cout << " unknown\n";
1709 std::cout << " Interpolation order: " << m_intpDissociation << "\n";
1710 }
1713 std::cout << " none\n";
1714 }
1715}

Referenced by Garfield::MediumMagboltz::PrintGas().

◆ ScaleAttachment()

double Garfield::MediumGas::ScaleAttachment ( const double &  eta) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 74 of file MediumGas.hh.

74 {
75 return eta * m_pressure / pressureTable;
76 }

◆ ScaleDiffusion()

double Garfield::MediumGas::ScaleDiffusion ( const double &  d) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 65 of file MediumGas.hh.

65 {
66 return d * sqrt(pressureTable / m_pressure);
67 }

◆ ScaleDiffusionTensor()

double Garfield::MediumGas::ScaleDiffusionTensor ( const double &  d) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 68 of file MediumGas.hh.

68 {
69 return d * pressureTable / m_pressure;
70 }

◆ ScaleElectricField()

double Garfield::MediumGas::ScaleElectricField ( const double &  e) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 59 of file MediumGas.hh.

59 {
60 return e * pressureTable / m_pressure;
61 }

◆ ScaleTownsend()

double Garfield::MediumGas::ScaleTownsend ( const double &  alpha) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 71 of file MediumGas.hh.

71 {
73 }

◆ SetAtomicNumber()

void Garfield::MediumGas::SetAtomicNumber ( const double &  z)
virtual

Reimplemented from Garfield::Medium.

Definition at line 220 of file MediumGas.cc.

220 {
221
222 std::cerr << m_className << "::SetAtomicNumber:\n";
223 std::cerr << " Effective Z cannot be changed"
224 << " directly to " << z << ".\n";
225 std::cerr << " Use SetComposition to define the gas mixture.\n";
226}

◆ SetAtomicWeight()

void Garfield::MediumGas::SetAtomicWeight ( const double &  a)
virtual

Reimplemented from Garfield::Medium.

Definition at line 228 of file MediumGas.cc.

228 {
229
230 std::cerr << m_className << "::SetAtomicWeight:\n";
231 std::cerr << " Effective A cannot be changed"
232 << " directly to " << a << ".\n";
233 std::cerr << " Use SetComposition to define the gas mixture.\n";
234}

◆ SetComposition()

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

Definition at line 64 of file MediumGas.cc.

69 {
70
71 // Make a backup copy of the gas composition.
72 std::string gasOld[m_nMaxGases];
73 for (unsigned int i = 0; i< m_nMaxGases; ++i) {
74 gasOld[i] = gas[i];
75 }
76 const unsigned int nComponentsOld = m_nComponents;
77 m_nComponents = 0;
78
79 // Find the gas name corresponding to the input string.
80 std::string gasname = "";
81 if (f1 > 0. && GetGasName(gas1, gasname)) {
82 gas[m_nComponents] = gasname;
85 }
86 if (f2 > 0. && GetGasName(gas2, gasname)) {
87 gas[m_nComponents] = gasname;
90 }
91 if (f3 > 0. && GetGasName(gas3, gasname)) {
92 gas[m_nComponents] = gasname;
95 }
96 if (f4 > 0. && GetGasName(gas4, gasname)) {
97 gas[m_nComponents] = gasname;
100 }
101 if (f5 > 0. && GetGasName(gas5, gasname)) {
102 gas[m_nComponents] = gasname;
105 }
106 if (f6 > 0. && GetGasName(gas6, gasname)) {
107 gas[m_nComponents] = gasname;
110 }
111
112 // Check if at least one valid ingredient was specified.
113 if (m_nComponents == 0) {
114 std::cerr << m_className << "::SetComposition:\n";
115 std::cerr << " Error setting the composition.\n";
116 std::cerr << " No valid ingredients were specified.\n";
117 return false;
118 }
119
120 // Establish the name.
121 m_name = "";
122 double sum = 0.;
123 for (unsigned int i = 0; i < m_nComponents; ++i) {
124 if (i > 0) m_name += "/";
125 m_name += gas[i];
126 sum += fraction[i];
127 }
128 // Normalise the fractions to one.
129 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
130 if (i < m_nComponents) {
131 fraction[i] /= sum;
132 } else {
133 fraction[i] = 0.;
134 }
135 }
136
137 // Set the atomic weight and number
138 for (unsigned int i = 0; i < m_nComponents; ++i) {
139 atWeight[i] = 0.;
140 atNum[i] = 0.;
141 GetGasInfo(gas[i], atWeight[i], atNum[i]);
142 }
143
144 // Print the composition.
145 std::cout << m_className << "::SetComposition:\n";
146 std::cout << " " << m_name;
147 if (m_nComponents > 1) {
148 std::cout << " (" << fraction[0] * 100;
149 for (unsigned int i = 1; i < m_nComponents; ++i) {
150 std::cout << "/" << fraction[i] * 100;
151 }
152 std::cout << ")";
153 }
154 std::cout << "\n";
155
156 // Force a recalculation of the collision rates.
157 m_isChanged = true;
158
159 // Copy the previous Penning transfer parameters.
160 double rPenningGasOld[m_nMaxGases];
161 double lambdaPenningGasOld[m_nMaxGases];
162 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
163 rPenningGasOld[i] = rPenningGas[i];
164 lambdaPenningGasOld[i] = lambdaPenningGas[i];
165 rPenningGas[i] = 0.;
166 lambdaPenningGas[i] = 0.;
167 }
168 for (unsigned int i = 0; i < m_nComponents; ++i) {
169 for (unsigned int j = 0; j < nComponentsOld; ++j) {
170 if (gas[i] == gasOld[j]) {
171 if (rPenningGasOld[j] > 0.) {
172 rPenningGas[i] = rPenningGasOld[j];
173 lambdaPenningGas[i] = lambdaPenningGasOld[i];
174 std::cout << m_className << "::SetComposition:\n";
175 std::cout << " Adopting Penning transfer parameters for " << gas[i]
176 << " from previous mixture.\n";
177 std::cout << " r = " << rPenningGas[i] << "\n";
178 std::cout << " lambda = " << lambdaPenningGas[i] << " cm\n";
179 }
180 }
181 }
182 }
183 return true;
184}

Referenced by GarfieldPhysics::InitializePhysics().

◆ SetExtrapolationMethodExcitationRates()

void Garfield::MediumGas::SetExtrapolationMethodExcitationRates ( const std::string  extrLow,
const std::string  extrHigh 
)

Definition at line 1825 of file MediumGas.cc.

1826 {
1827
1828 unsigned int iExtr = 0;
1829 if (GetExtrapolationIndex(extrLow, iExtr)) {
1830 m_extrLowExcRates = iExtr;
1831 } else {
1832 std::cerr << m_className << "::SetExtrapolationMethodExcitationRates:\n";
1833 std::cerr << " Unknown extrapolation method (" << extrLow << ")\n";
1834 }
1835 if (GetExtrapolationIndex(extrHigh, iExtr)) {
1836 m_extrHighExcRates = iExtr;
1837 } else {
1838 std::cerr << m_className << "::SetExtrapolationMethodExcitationRates:\n";
1839 std::cerr << " Unknown extrapolation method (" << extrHigh << ")\n";
1840 }
1841}
bool GetExtrapolationIndex(std::string extrStr, unsigned int &extrNb)
Definition: Medium.cc:2635

◆ SetExtrapolationMethodIonisationRates()

void Garfield::MediumGas::SetExtrapolationMethodIonisationRates ( const std::string  extrLow,
const std::string  extrHigh 
)

Definition at line 1843 of file MediumGas.cc.

1844 {
1845
1846 unsigned int iExtr = 0;
1847 if (GetExtrapolationIndex(extrLow, iExtr)) {
1848 m_extrLowIonRates = iExtr;
1849 } else {
1850 std::cerr << m_className << "::SetExtrapolationMethodIonisationRates:\n";
1851 std::cerr << " Unknown extrapolation method (" << extrLow << ")\n";
1852 }
1853 if (GetExtrapolationIndex(extrHigh, iExtr)) {
1854 m_extrHighIonRates = iExtr;
1855 } else {
1856 std::cerr << m_className << "::SetExtrapolationMethodIonisationRates:\n";
1857 std::cerr << " Unknown extrapolation method (" << extrHigh << ")\n";
1858 }
1859}

◆ SetInterpolationMethodExcitationRates()

void Garfield::MediumGas::SetInterpolationMethodExcitationRates ( const int  intrp)

Definition at line 1861 of file MediumGas.cc.

1861 {
1862
1863 if (intrp > 0) {
1864 m_intpExcRates = intrp;
1865 }
1866}

◆ SetInterpolationMethodIonisationRates()

void Garfield::MediumGas::SetInterpolationMethodIonisationRates ( const int  intrp)

Definition at line 1868 of file MediumGas.cc.

1868 {
1869
1870 if (intrp > 0) {
1871 m_intpIonRates = intrp;
1872 }
1873}

◆ SetMassDensity()

void Garfield::MediumGas::SetMassDensity ( const double &  rho)
virtual

Reimplemented from Garfield::Medium.

Definition at line 243 of file MediumGas.cc.

243 {
244
245 std::cerr << m_className << "::SetMassDensity:\n";
246 std::cerr << " Density cannot directly be changed to " << rho << ".\n";
247 std::cerr << " Use SetTemperature, SetPressure"
248 << " and SetComposition.\n";
249}

◆ SetNumberDensity()

void Garfield::MediumGas::SetNumberDensity ( const double &  n)
virtual

Reimplemented from Garfield::Medium.

Definition at line 236 of file MediumGas.cc.

236 {
237
238 std::cerr << m_className << "::SetNumberDensity:\n";
239 std::cerr << " Density cannot directly be changed to " << n << ".\n";
240 std::cerr << " Use SetTemperature and SetPressure.\n";
241}

◆ UnScaleElectricField()

double Garfield::MediumGas::UnScaleElectricField ( const double &  e) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 62 of file MediumGas.hh.

62 {
63 return e * m_pressure / pressureTable;
64 }

◆ WriteGasFile()

bool Garfield::MediumGas::WriteGasFile ( const std::string &  filename)

Definition at line 1021 of file MediumGas.cc.

1021 {
1022
1023 const int nMagboltzGases = 60;
1024 std::vector<double> mixture(nMagboltzGases);
1025 for (int i = nMagboltzGases; i--;) mixture[i] = 0.;
1026 // Set the gas mixture.
1027 for (unsigned int i = 0; i < m_nComponents; ++i) {
1028 int ng = 0;
1029 if (!GetGasNumberGasFile(gas[i], ng)) {
1030 std::cerr << m_className << "::WriteGasFile:\n";
1031 std::cerr << " Error retrieving gas number for gas " << gas[i]
1032 << ".\n";
1033 } else {
1034 mixture[ng - 1] = fraction[i] * 100.;
1035 }
1036 }
1037
1038 const int eFieldRes = m_nEfields;
1039 const int bFieldRes = m_nBfields;
1040 const int angRes = m_nAngles;
1041
1042 if (m_debug) {
1043 std::cout << m_className << "::WriteGasFile:\n";
1044 std::cout << " Writing gas tables to file: " << filename << "\n";
1045 }
1046
1047 std::ofstream outFile;
1048 outFile.open(filename.c_str(), std::ios::out);
1049 if (!outFile.is_open()) {
1050 std::cerr << m_className << "::WriteGasFile:\n";
1051 std::cerr << " File could not be opened.\n";
1052 outFile.close();
1053 return false;
1054 }
1055
1056 // Assemble the GASOK bits.
1057 std::string gasBits = "FFFFFFFFFFFFFFFFFFFF";
1058 if (m_hasElectronVelocityE) gasBits[0] = 'T';
1059 if (m_hasIonMobility) gasBits[1] = 'T';
1060 if (m_hasElectronDiffLong) gasBits[2] = 'T';
1061 if (m_hasElectronTownsend) gasBits[3] = 'T';
1062 // Custer size distribution; skippped
1063 if (m_hasElectronAttachment) gasBits[5] = 'T';
1064 // Lorentz angle; skipped
1065 if (m_hasElectronDiffTrans) gasBits[7] = 'T';
1066 if (m_hasElectronVelocityB) gasBits[8] = 'T';
1067 if (m_hasElectronVelocityExB) gasBits[9] = 'T';
1068 if (m_hasElectronDiffTens) gasBits[10] = 'T';
1069 if (m_hasIonDissociation) gasBits[11] = 'T';
1070 // SRIM, HEED; skipped
1071 if (m_hasExcRates) gasBits[14] = 'T';
1072 if (m_hasIonRates) gasBits[15] = 'T';
1073
1074 // Get the current time.
1075 time_t rawtime = time(0);
1076 tm timeinfo = *localtime(&rawtime);
1077 char datebuf[80] = {0};
1078 char timebuf[80] = {0};
1079 // Format date and time.
1080 strftime(datebuf, sizeof(datebuf) - 1, "%d/%m/%y", &timeinfo);
1081 strftime(timebuf, sizeof(timebuf) - 1, "%H.%M.%S", &timeinfo);
1082 // Set the member name.
1083 std::string member = "< none >";
1084 // Write the header.
1085 outFile << "*----.----1----.----2----.----3----.----4----.----"
1086 << "5----.----6----.----7----.----8----.----9----.---"
1087 << "10----.---11----.---12----.---13--\n";
1088 outFile << "% Created " << datebuf << " at " << timebuf << " ";
1089 outFile << member << " GAS ";
1090 // Add remark.
1091 std::string buffer;
1092 buffer = std::string(25, ' ');
1093 outFile << "\"none" << buffer << "\"\n";
1094 const int versionNumber = 12;
1095 outFile << " Version : " << versionNumber << "\n";
1096 outFile << " GASOK bits: " << gasBits << "\n";
1097 std::stringstream idStream;
1098 idStream.str("");
1099 idStream << m_name << ", p = " << pressureTable / AtmosphericPressure
1100 << " atm, T = " << temperatureTable << " K";
1101 std::string idString = idStream.str();
1102 outFile << " Identifier: " << std::setw(80) << std::left << idString << "\n";
1103 outFile << std::right;
1104 buffer = std::string(80, ' ');
1105 outFile << " Clusters : " << buffer << "\n";
1106 outFile << " Dimension : ";
1107 if (m_map2d) {
1108 outFile << "T ";
1109 } else {
1110 outFile << "F ";
1111 }
1112 outFile << std::setw(9) << eFieldRes << " " << std::setw(9) << angRes << " "
1113 << std::setw(9) << bFieldRes << " " << std::setw(9)
1114 << nExcListElements << " " << std::setw(9) << nIonListElements
1115 << "\n";
1116 outFile << " E fields \n";
1117 outFile << std::scientific << std::setw(15) << std::setprecision(8);
1118 for (int i = 0; i < eFieldRes; i++) {
1119 // List 5 values, then new line.
1120 outFile << std::setw(15) << eFields[i] / m_pressure;
1121 if ((i + 1) % 5 == 0) outFile << "\n";
1122 }
1123 if (eFieldRes % 5 != 0) outFile << "\n";
1124 outFile << " E-B angles \n";
1125 for (int i = 0; i < angRes; i++) {
1126 // List 5 values, then new line.
1127 outFile << std::setw(15) << bAngles[i];
1128 if ((i + 1) % 5 == 0) outFile << "\n";
1129 }
1130 if (angRes % 5 != 0) outFile << "\n";
1131 outFile << " B fields \n";
1132 for (int i = 0; i < bFieldRes; i++) {
1133 // List 5 values, then new line.
1134 // B fields are stored in hGauss (to be checked!).
1135 outFile << std::setw(15) << bFields[i] * 100.;
1136 if ((i + 1) % 5 == 0) outFile << "\n";
1137 }
1138 if (bFieldRes % 5 != 0) outFile << "\n";
1139 outFile << " Mixture: \n";
1140 for (int i = 0; i < nMagboltzGases; i++) {
1141 // List 5 values, then new line.
1142 outFile << std::setw(15) << mixture[i];
1143 if ((i + 1) % 5 == 0) outFile << "\n";
1144 }
1145 if (nMagboltzGases % 5 != 0) outFile << "\n";
1146 for (int i = 0; i < nExcListElements; i++) {
1147 outFile << " Excitation " << std::setw(5) << i + 1 << ": " << std::setw(45)
1148 << std::left << excitationList[i].label << " " << std::setw(15)
1149 << std::right << excitationList[i].energy << std::setw(15)
1150 << excitationList[i].prob << std::setw(15) << excitationList[i].rms
1151 << std::setw(15) << excitationList[i].dt << "\n";
1152 }
1153 for (int i = 0; i < nIonListElements; i++) {
1154 outFile << " Ionisation " << std::setw(5) << i + 1 << ": " << std::setw(45)
1155 << std::left << ionisationList[i].label << " " << std::setw(15)
1156 << std::right << ionisationList[i].energy << "\n";
1157 }
1158
1159 const double sqrtPressure = sqrt(pressureTable);
1160 const double logPressure = log(pressureTable);
1161
1162 outFile << " The gas tables follow:\n";
1163 int cnt = 0;
1164 for (int i = 0; i < eFieldRes; i++) {
1165 for (int j = 0; j < angRes; j++) {
1166 for (int k = 0; k < bFieldRes; k++) {
1167 double ve = 0., vb = 0., vexb = 0.;
1171 // Convert from cm / ns to cm / us.
1172 ve *= 1.e3;
1173 vb *= 1.e3;
1174 vexb *= 1.e3;
1175 double dl = 0., dt = 0.;
1176 if (m_hasElectronDiffLong) dl = tabElectronDiffLong[j][k][i];
1178 dl *= sqrtPressure;
1179 dt *= sqrtPressure;
1180 double alpha = -30., alpha0 = -30., eta = -30.;
1182 alpha = tabElectronTownsend[j][k][i];
1183 alpha0 = tabTownsendNoPenning[j][k][i];
1184 alpha -= logPressure;
1185 alpha0 -= logPressure;
1186 }
1188 eta = tabElectronAttachment[j][k][i];
1189 eta -= logPressure;
1190 }
1191 // Ion mobility
1192 double mu = 0.;
1193 if (m_hasIonMobility) mu = tabIonMobility[j][k][i];
1194 // Convert from cm2 / (V ns) to cm2 / (V us).
1195 mu *= 1.e3;
1196 // Lorentz angle
1197 double lor = 0.;
1198 // Dissociation coefficient
1199 double diss = -30.;
1201 diss = tabIonDissociation[j][k][i];
1202 diss -= logPressure;
1203 }
1204 // Set spline coefficient to dummy value.
1205 const double spl = 0.;
1206 // Write the values to file.
1207 outFile << std::setw(15);
1208 outFile << ve;
1209 ++cnt;
1210 if (cnt % 8 == 0) outFile << "\n";
1211 if (!m_map2d) {
1212 outFile << std::setw(15);
1213 outFile << spl;
1214 ++cnt;
1215 if (cnt % 8 == 0) outFile << "\n";
1216 }
1217 outFile << std::setw(15);
1218 outFile << vb;
1219 ++cnt;
1220 if (cnt % 8 == 0) outFile << "\n";
1221 if (!m_map2d) {
1222 outFile << std::setw(15);
1223 outFile << spl;
1224 ++cnt;
1225 if (cnt % 8 == 0) outFile << "\n";
1226 }
1227 outFile << std::setw(15);
1228 outFile << vexb;
1229 ++cnt;
1230 if (cnt % 8 == 0) outFile << "\n";
1231 if (!m_map2d) {
1232 outFile << std::setw(15);
1233 outFile << spl;
1234 ++cnt;
1235 if (cnt % 8 == 0) outFile << "\n";
1236 }
1237 outFile << std::setw(15);
1238 outFile << dl;
1239 ++cnt;
1240 if (cnt % 8 == 0) outFile << "\n";
1241 if (!m_map2d) {
1242 outFile << std::setw(15);
1243 outFile << spl;
1244 ++cnt;
1245 if (cnt % 8 == 0) outFile << "\n";
1246 }
1247 outFile << std::setw(15);
1248 outFile << dt;
1249 ++cnt;
1250 if (cnt % 8 == 0) outFile << "\n";
1251 if (!m_map2d) {
1252 outFile << std::setw(15);
1253 outFile << spl;
1254 ++cnt;
1255 if (cnt % 8 == 0) outFile << "\n";
1256 }
1257 outFile << std::setw(15);
1258 outFile << alpha;
1259 ++cnt;
1260 if (cnt % 8 == 0) outFile << "\n";
1261 if (!m_map2d) {
1262 outFile << std::setw(15);
1263 outFile << spl;
1264 ++cnt;
1265 if (cnt % 8 == 0) outFile << "\n";
1266 }
1267 outFile << std::setw(15);
1268 outFile << alpha0;
1269 ++cnt;
1270 if (cnt % 8 == 0) outFile << "\n";
1271 outFile << std::setw(15);
1272 outFile << eta;
1273 ++cnt;
1274 if (cnt % 8 == 0) outFile << "\n";
1275 outFile << std::setw(15);
1276 if (!m_map2d) {
1277 outFile << spl;
1278 ++cnt;
1279 if (cnt % 8 == 0) outFile << "\n";
1280 outFile << std::setw(15);
1281 }
1282 outFile << mu;
1283 ++cnt;
1284 if (cnt % 8 == 0) outFile << "\n";
1285 if (!m_map2d) {
1286 outFile << std::setw(15);
1287 outFile << spl;
1288 ++cnt;
1289 if (cnt % 8 == 0) outFile << "\n";
1290 }
1291 outFile << std::setw(15);
1292 outFile << lor;
1293 ++cnt;
1294 if (cnt % 8 == 0) outFile << "\n";
1295 if (!m_map2d) {
1296 outFile << std::setw(15);
1297 outFile << spl;
1298 ++cnt;
1299 if (cnt % 8 == 0) outFile << "\n";
1300 }
1301 outFile << std::setw(15);
1302 outFile << diss;
1303 ++cnt;
1304 if (cnt % 8 == 0) outFile << "\n";
1305 if (!m_map2d) {
1306 outFile << std::setw(15);
1307 outFile << spl;
1308 ++cnt;
1309 if (cnt % 8 == 0) outFile << "\n";
1310 }
1311 outFile << std::setw(15);
1312 for (int l = 0; l < 6; ++l) {
1313 double diff = 0.;
1315 diff = tabElectronDiffTens[l][j][k][i];
1316 diff *= pressureTable;
1317 }
1318 outFile << std::setw(15);
1319 outFile << diff;
1320 ++cnt;
1321 if (cnt % 8 == 0) outFile << "\n";
1322 if (!m_map2d) {
1323 outFile << std::setw(15) << spl;
1324 ++cnt;
1325 if (cnt % 8 == 0) outFile << "\n";
1326 }
1327 }
1328 if (m_hasExcRates && nExcListElements > 0) {
1329 for (int l = 0; l < nExcListElements; l++) {
1330 outFile << std::setw(15);
1331 outFile << tabExcRates[l][j][k][i];
1332 ++cnt;
1333 if (cnt % 8 == 0) outFile << "\n";
1334 if (!m_map2d) {
1335 outFile << std::setw(15) << spl;
1336 ++cnt;
1337 if (cnt % 8 == 0) outFile << "\n";
1338 }
1339 }
1340 }
1341 if (m_hasIonRates && nIonListElements > 0) {
1342 for (int l = 0; l < nIonListElements; l++) {
1343 outFile << std::setw(15);
1344 outFile << tabIonRates[l][j][k][i];
1345 ++cnt;
1346 if (cnt % 8 == 0) outFile << "\n";
1347 if (!m_map2d) {
1348 outFile << std::setw(15) << spl;
1349 ++cnt;
1350 if (cnt % 8 == 0) outFile << "\n";
1351 }
1352 }
1353 }
1354 }
1355 }
1356 if (cnt % 8 != 0) outFile << "\n";
1357 cnt = 0;
1358 }
1359
1360 // Extrapolation methods
1361 int hExtrap[13], lExtrap[13];
1362 int interpMeth[13];
1363
1364 hExtrap[0] = hExtrap[1] = hExtrap[2] = m_extrHighVelocity;
1365 lExtrap[0] = lExtrap[1] = lExtrap[2] = m_extrLowVelocity;
1366 interpMeth[0] = interpMeth[1] = interpMeth[2] = m_intpVelocity;
1367 hExtrap[3] = hExtrap[8] = hExtrap[10] = m_extrHighDiffusion;
1368 lExtrap[3] = lExtrap[8] = lExtrap[10] = m_extrLowDiffusion;
1369 interpMeth[3] = interpMeth[8] = interpMeth[10] = m_intpDiffusion;
1370 hExtrap[4] = m_extrHighTownsend;
1371 lExtrap[4] = m_extrLowTownsend;
1372 interpMeth[4] = m_intpTownsend;
1373 hExtrap[5] = m_extrHighAttachment;
1374 lExtrap[5] = m_extrLowAttachment;
1375 interpMeth[5] = m_intpAttachment;
1376 hExtrap[6] = m_extrHighMobility;
1377 lExtrap[6] = m_extrLowMobility;
1378 interpMeth[6] = m_intpMobility;
1379 // Lorentz angle
1380 hExtrap[7] = 1;
1381 lExtrap[7] = 0;
1382 interpMeth[7] = 2;
1383 hExtrap[9] = m_extrHighDissociation;
1384 lExtrap[9] = m_extrLowDissociation;
1385 interpMeth[9] = m_intpDissociation;
1386 hExtrap[11] = m_extrHighExcRates;
1387 lExtrap[11] = m_extrLowExcRates;
1388 interpMeth[11] = m_intpExcRates;
1389 hExtrap[12] = m_extrHighIonRates;
1390 lExtrap[12] = m_extrLowIonRates;
1391 interpMeth[12] = m_intpIonRates;
1392
1393 outFile << " H Extr: ";
1394 for (int i = 0; i < 13; i++) {
1395 outFile << std::setw(5) << hExtrap[i];
1396 }
1397 outFile << "\n";
1398 outFile << " L Extr: ";
1399 for (int i = 0; i < 13; i++) {
1400 outFile << std::setw(5) << lExtrap[i];
1401 }
1402 outFile << "\n";
1403 outFile << " Thresholds: " << std::setw(10) << thrElectronTownsend
1404 << std::setw(10) << thrElectronAttachment << std::setw(10)
1405 << thrIonDissociation << "\n";
1406 outFile << " Interp: ";
1407 for (int i = 0; i < 13; i++) {
1408 outFile << std::setw(5) << interpMeth[i];
1409 }
1410 outFile << "\n";
1411 outFile << " A =" << std::setw(15) << 0. << ","
1412 << " Z =" << std::setw(15) << 0. << ","
1413 << " EMPROB=" << std::setw(15) << 0. << ","
1414 << " EPAIR =" << std::setw(15) << 0. << "\n";
1415 double ionDiffLong = 0., ionDiffTrans = 0.;
1416 if (m_hasIonDiffLong) ionDiffLong = tabIonDiffLong[0][0][0];
1417 if (m_hasIonDiffTrans) ionDiffTrans = tabIonDiffTrans[0][0][0];
1418 outFile << " Ion diffusion: " << std::setw(15) << ionDiffLong << std::setw(15)
1419 << ionDiffTrans << "\n";
1420 outFile << " CMEAN =" << std::setw(15) << 0. << ","
1421 << " RHO =" << std::setw(15) << 0. << ","
1422 << " PGAS =" << std::setw(15) << pressureTable << ","
1423 << " TGAS =" << std::setw(15) << temperatureTable << "\n";
1424 outFile << " CLSTYP : NOT SET \n";
1425 buffer = std::string(80, ' ');
1426 outFile << " FCNCLS : " << buffer << "\n";
1427 outFile << " NCLS : " << std::setw(10) << 0 << "\n";
1428 outFile << " Average : " << std::setw(25) << std::setprecision(18) << 0.
1429 << "\n";
1430 outFile << " Heed initialisation done: F\n";
1431 outFile << " SRIM initialisation done: F\n";
1432 outFile.close();
1433
1434 return true;
1435}
bool GetGasNumberGasFile(const std::string input, int &number) const
Definition: MediumGas.cc:2601

Member Data Documentation

◆ atNum

double Garfield::MediumGas::atNum[m_nMaxGases]
protected

Definition at line 88 of file MediumGas.hh.

Referenced by GetAtomicNumber(), LoadGasFile(), MediumGas(), and SetComposition().

◆ atWeight

double Garfield::MediumGas::atWeight[m_nMaxGases]
protected

Definition at line 87 of file MediumGas.hh.

Referenced by GetAtomicWeight(), LoadGasFile(), MediumGas(), and SetComposition().

◆ excitationList

std::vector<excListElement> Garfield::MediumGas::excitationList
protected

◆ fraction

◆ gas

◆ ionisationList

std::vector<ionListElement> Garfield::MediumGas::ionisationList
protected

◆ lambdaPenningGas

double Garfield::MediumGas::lambdaPenningGas[m_nMaxGases]
protected

◆ lambdaPenningGlobal

double Garfield::MediumGas::lambdaPenningGlobal
protected

◆ m_extrHighExcRates

unsigned int Garfield::MediumGas::m_extrHighExcRates
protected

◆ m_extrHighIonRates

unsigned int Garfield::MediumGas::m_extrHighIonRates
protected

◆ m_extrLowExcRates

unsigned int Garfield::MediumGas::m_extrLowExcRates
protected

◆ m_extrLowIonRates

unsigned int Garfield::MediumGas::m_extrLowIonRates
protected

◆ m_hasExcRates

bool Garfield::MediumGas::m_hasExcRates
protected

◆ m_hasIonRates

bool Garfield::MediumGas::m_hasIonRates
protected

◆ m_intpExcRates

unsigned int Garfield::MediumGas::m_intpExcRates
protected

◆ m_intpIonRates

unsigned int Garfield::MediumGas::m_intpIonRates
protected

◆ m_nMaxGases

const unsigned int Garfield::MediumGas::m_nMaxGases = 6
staticprotected

◆ nExcListElements

int Garfield::MediumGas::nExcListElements
protected

◆ nIonListElements

int Garfield::MediumGas::nIonListElements
protected

◆ pressureTable

◆ rPenningGas

double Garfield::MediumGas::rPenningGas[m_nMaxGases]
protected

◆ rPenningGlobal

double Garfield::MediumGas::rPenningGlobal
protected

◆ tabExcRates

std::vector<std::vector<std::vector<std::vector<double> > > > Garfield::MediumGas::tabExcRates
protected

◆ tabIonRates

std::vector<std::vector<std::vector<std::vector<double> > > > Garfield::MediumGas::tabIonRates
protected

◆ tabTownsendNoPenning

std::vector<std::vector<std::vector<double> > > Garfield::MediumGas::tabTownsendNoPenning
protected

◆ temperatureTable

double Garfield::MediumGas::temperatureTable
protected

◆ usePenning


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