Garfield++ v2r0
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

Base class for gas media. More...

#include <MediumGas.hh>

+ Inheritance diagram for Garfield::MediumGas:

Classes

struct  excListElement
 
struct  ionListElement
 

Public Member Functions

 MediumGas ()
 
virtual ~MediumGas ()
 
bool IsGas () const
 
bool SetComposition (const std::string &gas1, const double f1=1., const std::string &gas2="", const double f2=0., const std::string &gas3="", const double f3=0., const std::string &gas4="", const double f4=0., const std::string &gas5="", const double f5=0., const std::string &gas6="", const double f6=0.)
 
void GetComposition (std::string &gas1, double &f1, std::string &gas2, double &f2, std::string &gas3, double &f3, std::string &gas4, double &f4, std::string &gas5, double &f5, std::string &gas6, double &f6)
 
void GetComponent (const unsigned int i, std::string &label, double &f)
 
void SetAtomicNumber (const double z)
 
double GetAtomicNumber () const
 
void SetAtomicWeight (const double a)
 
double GetAtomicWeight () const
 
void SetNumberDensity (const double n)
 
double GetNumberDensity () const
 
void SetMassDensity (const double rho)
 
double GetMassDensity () const
 
bool LoadGasFile (const std::string &filename)
 
bool WriteGasFile (const std::string &filename)
 
void PrintGas ()
 
bool LoadIonMobility (const std::string &filename)
 
void SetExtrapolationMethodExcitationRates (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodIonisationRates (const std::string &extrLow, const std::string &extrHigh)
 
void SetInterpolationMethodExcitationRates (const int intrp)
 
void SetInterpolationMethodIonisationRates (const int intrp)
 
double ScaleElectricField (const double e) const
 
double UnScaleElectricField (const double e) const
 
double ScaleDiffusion (const double d) const
 
double ScaleDiffusionTensor (const double d) const
 
double ScaleTownsend (const double alpha) const
 
double ScaleAttachment (const double eta) const
 
double ScaleLorentzAngle (const double lor) const
 
bool GetPhotoabsorptionCrossSection (const double &e, double &sigma, const unsigned int &i)
 
- Public Member Functions inherited from Garfield::Medium
 Medium ()
 
virtual ~Medium ()
 
int GetId () const
 
const std::string & GetName () const
 
virtual bool IsGas () const
 
virtual bool IsSemiconductor () const
 
void SetTemperature (const double t)
 
double GetTemperature () const
 
void SetPressure (const double p)
 
double GetPressure () const
 
void SetDielectricConstant (const double eps)
 
double GetDielectricConstant () const
 
unsigned int GetNumberOfComponents () const
 
virtual void GetComponent (const unsigned int i, std::string &label, double &f)
 
virtual void SetAtomicNumber (const double z)
 
virtual double GetAtomicNumber () const
 
virtual void SetAtomicWeight (const double a)
 
virtual double GetAtomicWeight () const
 
virtual void SetNumberDensity (const double n)
 
virtual double GetNumberDensity () const
 
virtual void SetMassDensity (const double rho)
 
virtual double GetMassDensity () const
 
virtual void EnableDrift ()
 
void DisableDrift ()
 
virtual void EnablePrimaryIonisation ()
 
void DisablePrimaryIonisation ()
 
bool IsDriftable () const
 
bool IsMicroscopic () const
 
bool IsIonisable () const
 
void SetW (const double w)
 
double GetW ()
 
void SetFanoFactor (const double f)
 
double GetFanoFactor ()
 
virtual bool ElectronVelocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
 
virtual bool ElectronDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
 
virtual bool ElectronDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double cov[3][3])
 
virtual bool ElectronTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
virtual bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
virtual bool ElectronLorentzAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
 
virtual double GetElectronEnergy (const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
 
virtual void GetElectronMomentum (const double e, double &px, double &py, double &pz, int &band)
 
virtual double GetElectronNullCollisionRate (const int band=0)
 
virtual double GetElectronCollisionRate (const double e, const int band=0)
 
virtual bool GetElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, int &nion, int &ndxc, int &band)
 
virtual unsigned int GetNumberOfIonisationProducts () const
 
virtual bool GetIonisationProduct (const unsigned int i, int &type, double &energy) const
 
virtual unsigned int GetNumberOfDeexcitationProducts () const
 
virtual bool GetDeexcitationProduct (const unsigned int i, double &t, double &s, int &type, double &energy) const
 
virtual bool HoleVelocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
 
virtual bool HoleDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
 
virtual bool HoleDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double cov[3][3])
 
virtual bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 
virtual bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 
virtual bool IonVelocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
 
virtual bool IonDiffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
 
virtual bool IonDissociation (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &diss)
 
void SetFieldGrid (double emin, double emax, int ne, bool logE, double bmin=0., double bmax=0., int nb=1, double amin=0., double amax=0., int na=1)
 
void SetFieldGrid (const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles)
 
void GetFieldGrid (std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles)
 
bool GetElectronVelocityE (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetElectronVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetElectronVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetElectronLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 
bool GetElectronTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 
bool GetElectronTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &alpha)
 
bool GetElectronAttachment (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &eta)
 
bool GetElectronLorentzAngle (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &lor)
 
bool GetHoleVelocityE (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetHoleVelocityExB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetHoleVelocityB (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &v)
 
bool GetHoleLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 
bool GetHoleTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 
bool GetHoleTownsend (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &alpha)
 
bool GetHoleAttachment (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &eta)
 
bool GetIonMobility (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &mu)
 
bool GetIonLongitudinalDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dl)
 
bool GetIonTransverseDiffusion (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &dt)
 
bool GetIonDissociation (const unsigned int ie, const unsigned int ib, const unsigned int ia, double &diss)
 
void ResetElectronVelocity ()
 
void ResetElectronDiffusion ()
 
void ResetElectronTownsend ()
 
void ResetElectronAttachment ()
 
void ResetElectronLorentzAngle ()
 
void ResetHoleVelocity ()
 
void ResetHoleDiffusion ()
 
void ResetHoleTownsend ()
 
void ResetHoleAttachment ()
 
void ResetIonMobility ()
 
void ResetIonDiffusion ()
 
void ResetIonDissociation ()
 
bool SetIonMobility (const unsigned int ie, const unsigned int ib, const unsigned int ia, const double mu)
 
bool SetIonMobility (const std::vector< double > &fields, const std::vector< double > &mobilities)
 
void SetExtrapolationMethodVelocity (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodDiffusion (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodTownsend (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodAttachment (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodIonMobility (const std::string &extrLow, const std::string &extrHigh)
 
void SetExtrapolationMethodIonDissociation (const std::string &extrLow, const std::string &extrHigh)
 
void SetInterpolationMethodVelocity (const unsigned int intrp)
 
void SetInterpolationMethodDiffusion (const unsigned int intrp)
 
void SetInterpolationMethodTownsend (const unsigned int intrp)
 
void SetInterpolationMethodAttachment (const unsigned int intrp)
 
void SetInterpolationMethodIonMobility (const unsigned int intrp)
 
void SetInterpolationMethodIonDissociation (const unsigned int intrp)
 
virtual double ScaleElectricField (const double e) const
 
virtual double UnScaleElectricField (const double e) const
 
virtual double ScaleVelocity (const double v) const
 
virtual double ScaleDiffusion (const double d) const
 
virtual double ScaleDiffusionTensor (const double d) const
 
virtual double ScaleTownsend (const double alpha) const
 
virtual double ScaleAttachment (const double eta) const
 
virtual double ScaleLorentzAngle (const double lor) const
 
virtual double ScaleDissociation (const double diss) const
 
virtual bool GetOpticalDataRange (double &emin, double &emax, const unsigned int i=0)
 
virtual bool GetDielectricFunction (const double e, double &eps1, double &eps2, const unsigned int i=0)
 
virtual bool GetPhotoAbsorptionCrossSection (const double e, double &sigma, const unsigned int i=0)
 
virtual double GetPhotonCollisionRate (const double e)
 
virtual bool GetPhotonCollision (const double e, int &type, int &level, double &e1, double &ctheta, int &nsec, double &esec)
 
void EnableDebugging ()
 
void DisableDebugging ()
 

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 GetAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const double e, const double b) const
 
double Interpolate1D (const double e, const std::vector< double > &table, const std::vector< double > &fields, const unsigned int intpMeth, const int jExtr, const int iExtr)
 
bool GetExtrapolationIndex (std::string extrStr, unsigned int &extrNb)
 
void CloneTable (std::vector< std::vector< std::vector< double > > > &tab, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int intp, const unsigned int extrLow, const unsigned int extrHigh, const double init, const std::string &label)
 
void CloneTensor (std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const unsigned int n, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int intp, const unsigned int extrLow, const unsigned int extrHigh, const double init, const std::string &label)
 
void InitParamArrays (const unsigned int eRes, const unsigned int bRes, const unsigned int aRes, std::vector< std::vector< std::vector< double > > > &tab, const double val)
 
void InitParamTensor (const unsigned int eRes, const unsigned int bRes, const unsigned int aRes, const unsigned int tRes, std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const double val)
 

Protected Attributes

std::string m_gas [m_nMaxGases]
 
double m_fraction [m_nMaxGases]
 
double m_atWeight [m_nMaxGases]
 
double m_atNum [m_nMaxGases]
 
bool m_usePenning
 
double m_rPenningGlobal
 
double m_rPenningGas [m_nMaxGases]
 
double m_lambdaPenningGlobal
 
double m_lambdaPenningGas [m_nMaxGases]
 
double m_pressureTable
 
double m_temperatureTable
 
std::vector< std::vector< std::vector< double > > > m_tabTownsendNoPenning
 
bool m_hasExcRates
 
bool m_hasIonRates
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_tabExcRates
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_tabIonRates
 
std::vector< excListElementm_excitationList
 
std::vector< ionListElementm_ionisationList
 
unsigned int m_extrLowExcRates
 
unsigned int m_extrHighExcRates
 
unsigned int m_extrLowIonRates
 
unsigned int m_extrHighIonRates
 
unsigned int m_intpExcRates
 
unsigned int m_intpIonRates
 
- Protected Attributes inherited from Garfield::Medium
std::string m_className
 
int m_id
 
std::string m_name
 
double m_temperature
 
double m_pressure
 
double m_epsilon
 
unsigned int m_nComponents
 
double m_z
 
double m_a
 
double m_density
 
bool m_driftable
 
bool m_microscopic
 
bool m_ionisable
 
double m_w
 
double m_fano
 
bool m_isChanged
 
bool m_debug
 
std::vector< double > m_eFields
 
std::vector< double > m_bFields
 
std::vector< double > m_bAngles
 
bool m_map2d
 
bool m_hasElectronVelocityE
 
bool m_hasElectronVelocityB
 
bool m_hasElectronVelocityExB
 
bool m_hasElectronDiffLong
 
bool m_hasElectronDiffTrans
 
bool m_hasElectronDiffTens
 
bool m_hasElectronAttachment
 
bool m_hasElectronLorentzAngle
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityE
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityExB
 
std::vector< std::vector< std::vector< double > > > tabElectronVelocityB
 
std::vector< std::vector< std::vector< double > > > tabElectronDiffLong
 
std::vector< std::vector< std::vector< double > > > tabElectronDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabElectronTownsend
 
std::vector< std::vector< std::vector< double > > > tabElectronAttachment
 
std::vector< std::vector< std::vector< double > > > tabElectronLorentzAngle
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabElectronDiffTens
 
bool m_hasHoleVelocityE
 
bool m_hasHoleVelocityB
 
bool m_hasHoleVelocityExB
 
bool m_hasHoleDiffLong
 
bool m_hasHoleDiffTrans
 
bool m_hasHoleDiffTens
 
bool m_hasHoleTownsend
 
bool m_hasHoleAttachment
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityE
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityExB
 
std::vector< std::vector< std::vector< double > > > tabHoleVelocityB
 
std::vector< std::vector< std::vector< double > > > tabHoleDiffLong
 
std::vector< std::vector< std::vector< double > > > tabHoleDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabHoleTownsend
 
std::vector< std::vector< std::vector< double > > > tabHoleAttachment
 
std::vector< std::vector< std::vector< std::vector< double > > > > tabHoleDiffTens
 
bool m_hasIonMobility
 
bool m_hasIonDiffLong
 
bool m_hasIonDiffTrans
 
bool m_hasIonDissociation
 
std::vector< std::vector< std::vector< double > > > tabIonMobility
 
std::vector< std::vector< std::vector< double > > > tabIonDiffLong
 
std::vector< std::vector< std::vector< double > > > tabIonDiffTrans
 
std::vector< std::vector< std::vector< double > > > tabIonDissociation
 
int thrElectronTownsend
 
int thrElectronAttachment
 
int thrHoleTownsend
 
int thrHoleAttachment
 
int thrIonDissociation
 
unsigned int m_extrLowVelocity
 
unsigned int m_extrHighVelocity
 
unsigned int m_extrLowDiffusion
 
unsigned int m_extrHighDiffusion
 
unsigned int m_extrLowTownsend
 
unsigned int m_extrHighTownsend
 
unsigned int m_extrLowAttachment
 
unsigned int m_extrHighAttachment
 
unsigned int m_extrLowLorentzAngle
 
unsigned int m_extrHighLorentzAngle
 
unsigned int m_extrLowMobility
 
unsigned int m_extrHighMobility
 
unsigned int m_extrLowDissociation
 
unsigned int m_extrHighDissociation
 
unsigned int m_intpVelocity
 
unsigned int m_intpDiffusion
 
unsigned int m_intpTownsend
 
unsigned int m_intpAttachment
 
unsigned int m_intpLorentzAngle
 
unsigned int m_intpMobility
 
unsigned int m_intpDissociation
 

Static Protected Attributes

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

Detailed Description

Base class for gas media.

Definition at line 13 of file MediumGas.hh.

Constructor & Destructor Documentation

◆ MediumGas()

Garfield::MediumGas::MediumGas ( )

Definition at line 16 of file MediumGas.cc.

17 : Medium(),
18 m_usePenning(false),
23 m_hasExcRates(false),
24 m_hasIonRates(false),
31
32 m_className = "MediumGas";
33
34 // Default gas mixture: pure argon
35 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
36 m_fraction[i] = 0.;
37 m_gas[i] = "";
38 m_atWeight[i] = 0.;
39 m_atNum[i] = 0.;
40 }
41 m_gas[0] = "Ar";
42 m_fraction[0] = 1.;
43 m_name = m_gas[0];
45
46 m_isChanged = true;
47
50
51 // Initialise Penning parameters
52 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
53 m_rPenningGas[i] = 0.;
54 m_lambdaPenningGas[i] = 0.;
55 }
56
57}
double m_rPenningGlobal
Definition: MediumGas.hh:99
std::string m_gas[m_nMaxGases]
Definition: MediumGas.hh:90
bool GetGasInfo(const std::string &gasname, double &a, double &z) const
Definition: MediumGas.cc:1861
double m_lambdaPenningGlobal
Definition: MediumGas.hh:102
double m_rPenningGas[m_nMaxGases]
Definition: MediumGas.hh:100
double m_atNum[m_nMaxGases]
Definition: MediumGas.hh:93
double m_temperatureTable
Definition: MediumGas.hh:108
unsigned int m_intpExcRates
Definition: MediumGas.hh:137
unsigned int m_extrHighExcRates
Definition: MediumGas.hh:135
double m_atWeight[m_nMaxGases]
Definition: MediumGas.hh:92
static const unsigned int m_nMaxGases
Definition: MediumGas.hh:87
double m_fraction[m_nMaxGases]
Definition: MediumGas.hh:91
unsigned int m_intpIonRates
Definition: MediumGas.hh:138
unsigned int m_extrHighIonRates
Definition: MediumGas.hh:136
unsigned int m_extrLowIonRates
Definition: MediumGas.hh:136
unsigned int m_extrLowExcRates
Definition: MediumGas.hh:135
double m_lambdaPenningGas[m_nMaxGases]
Definition: MediumGas.hh:103
double m_pressure
Definition: Medium.hh:305
std::string m_name
Definition: Medium.hh:301
virtual void EnableDrift()
Definition: Medium.hh:52
virtual void EnablePrimaryIonisation()
Definition: Medium.hh:54
std::string m_className
Definition: Medium.hh:294
bool m_isChanged
Definition: Medium.hh:326
double m_temperature
Definition: Medium.hh:303

◆ ~MediumGas()

virtual Garfield::MediumGas::~MediumGas ( )
inlinevirtual

Definition at line 19 of file MediumGas.hh.

19{}

Member Function Documentation

◆ GetAtomicNumber()

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

Reimplemented from Garfield::Medium.

Definition at line 241 of file MediumGas.cc.

241 {
242
243 // Effective Z, weighted by the fractions of the components.
244 double z = 0.;
245 for (unsigned int i = 0; i < m_nComponents; ++i) {
246 z += m_atNum[i] * m_fraction[i];
247 }
248 return z;
249}
unsigned int m_nComponents
Definition: Medium.hh:309

◆ GetAtomicWeight()

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

Reimplemented from Garfield::Medium.

Definition at line 219 of file MediumGas.cc.

219 {
220
221 // Effective A, weighted by the fractions of the components.
222 double a = 0.;
223 for (unsigned int i = 0; i < m_nComponents; ++i) {
224 a += m_atWeight[i] * m_fraction[i];
225 }
226 return a;
227}

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 177 of file MediumGas.cc.

178 {
179
180 if (i >= m_nComponents) {
181 std::cerr << m_className << "::GetComponent:\n Index out of range.\n";
182 label = "";
183 f = 0.;
184 return;
185 }
186
187 label = m_gas[i];
188 f = m_fraction[i];
189}

◆ 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 158 of file MediumGas.cc.

161 {
162
163 gas1 = m_gas[0];
164 gas2 = m_gas[1];
165 gas3 = m_gas[2];
166 gas4 = m_gas[3];
167 gas5 = m_gas[4];
168 gas6 = m_gas[5];
169 f1 = m_fraction[0];
170 f2 = m_fraction[1];
171 f3 = m_fraction[2];
172 f4 = m_fraction[3];
173 f5 = m_fraction[4];
174 f6 = m_fraction[5];
175}

◆ GetGasInfo()

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

Definition at line 1861 of file MediumGas.cc.

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

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 2031 of file MediumGas.cc.

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

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 2227 of file MediumGas.cc.

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

◆ GetGasNumberGasFile()

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

Definition at line 2587 of file MediumGas.cc.

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

Referenced by WriteGasFile().

◆ GetMassDensity()

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

Reimplemented from Garfield::Medium.

Definition at line 236 of file MediumGas.cc.

236 {
237
238 return GetNumberDensity() * GetAtomicWeight() * AtomicMassUnit;
239}
double GetAtomicWeight() const
Definition: MediumGas.cc:219
double GetNumberDensity() const
Definition: MediumGas.cc:229

◆ GetNumberDensity()

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

Reimplemented from Garfield::Medium.

Definition at line 229 of file MediumGas.cc.

229 {
230
231 // Ideal gas law.
232 return LoschmidtNumber * (m_pressure / AtmosphericPressure) *
233 (ZeroCelsius / m_temperature);
234}

Referenced by GetMassDensity(), and LoadIonMobility().

◆ GetPhotoabsorptionCrossSection()

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

Definition at line 2876 of file MediumGas.cc.

2877 {
2878
2879 if (i >= m_nMaxGases) {
2880 std::cerr << m_className << "::GetPhotoabsorptionCrossSection:\n";
2881 std::cerr << " Index (" << i << ") out of range.\n";
2882 return false;
2883 }
2884
2885 OpticalData optData;
2886 if (!optData.IsAvailable(m_gas[i])) return false;
2887 double eta = 0.;
2888 return optData.GetPhotoabsorptionCrossSection(m_gas[i], e, sigma, eta);
2889}

◆ IsGas()

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

Reimplemented from Garfield::Medium.

Definition at line 21 of file MediumGas.hh.

21{ return true; }

◆ LoadGasFile()

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

Definition at line 251 of file MediumGas.cc.

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

Referenced by GarfieldPhysics::InitializePhysics().

◆ LoadIonMobility()

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

Definition at line 1707 of file MediumGas.cc.

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

◆ PrintGas()

void Garfield::MediumGas::PrintGas ( )

Definition at line 1405 of file MediumGas.cc.

1405 {
1406
1407 // Print a summary.
1408 std::cout << m_className << "::PrintGas:\n";
1409 std::cout << " Gas composition: " << m_name;
1410 if (m_nComponents > 1) {
1411 std::cout << " (" << m_fraction[0] * 100;
1412 for (unsigned int i = 1; i < m_nComponents; ++i) {
1413 std::cout << "/" << m_fraction[i] * 100;
1414 }
1415 std::cout << ")";
1416 }
1417 std::cout << "\n";
1418 std::cout << " Pressure: " << m_pressure << " Torr\n";
1419 std::cout << " Temperature: " << m_temperature << " K\n";
1420 std::cout << " Gas file:\n";
1421 std::cout << " Pressure: " << m_pressureTable << " Torr\n";
1422 std::cout << " Temperature: " << m_temperatureTable << " K\n";
1423 if (m_eFields.size() > 1) {
1424 std::cout << " Electric field range: " << m_eFields[0] << " - "
1425 << m_eFields.back() << " V/cm in " << m_eFields.size() - 1
1426 << " steps.\n";
1427 } else if (m_eFields.size() == 1) {
1428 std::cout << " Electric field: " << m_eFields[0] << " V/cm\n";
1429 } else {
1430 std::cout << " Electric field range: not set\n";
1431 }
1432 if (m_bFields.size() > 1) {
1433 std::cout << " Magnetic field range: " << m_bFields[0] << " - "
1434 << m_bFields.back() << " T in " << m_bFields.size() - 1
1435 << " steps.\n";
1436 } else if (m_bFields.size() == 1) {
1437 std::cout << " Magnetic field: " << m_bFields[0] << "\n";
1438 } else {
1439 std::cout << " Magnetic field range: not set\n";
1440 }
1441 if (m_bAngles.size() > 1) {
1442 std::cout << " Angular range: " << m_bAngles[0] << " - "
1443 << m_bAngles.back() << " in " << m_bAngles.size() - 1
1444 << " steps.\n";
1445 } else if (m_bAngles.size() == 1) {
1446 std::cout << " Angle between E and B: " << m_bAngles[0] << "\n";
1447 } else {
1448 std::cout << " Angular range: not set\n";
1449 }
1450
1451 std::cout << " Available electron transport data:\n";
1453 std::cout << " Velocity along E\n";
1454 }
1456 std::cout << " Velocity along Bt\n";
1457 }
1459 std::cout << " Velocity along ExB\n";
1460 }
1462 std::cout << " Low field extrapolation: ";
1463 if (m_extrLowVelocity == 0)
1464 std::cout << " constant\n";
1465 else if (m_extrLowVelocity == 1)
1466 std::cout << " linear\n";
1467 else if (m_extrLowVelocity == 2)
1468 std::cout << " exponential\n";
1469 else
1470 std::cout << " unknown\n";
1471 std::cout << " High field extrapolation: ";
1472 if (m_extrHighVelocity == 0)
1473 std::cout << " constant\n";
1474 else if (m_extrHighVelocity == 1)
1475 std::cout << " linear\n";
1476 else if (m_extrHighVelocity == 2)
1477 std::cout << " exponential\n";
1478 else
1479 std::cout << " unknown\n";
1480 std::cout << " Interpolation order: " << m_intpVelocity << "\n";
1481 }
1483 std::cout << " Longitudinal diffusion coefficient\n";
1484 }
1486 std::cout << " Transverse diffusion coefficient\n";
1487 }
1489 std::cout << " Diffusion tensor\n";
1490 }
1492 std::cout << " Low field extrapolation: ";
1493 if (m_extrLowDiffusion == 0)
1494 std::cout << " constant\n";
1495 else if (m_extrLowDiffusion == 1)
1496 std::cout << " linear\n";
1497 else if (m_extrLowDiffusion == 2)
1498 std::cout << " exponential\n";
1499 else
1500 std::cout << " unknown\n";
1501 std::cout << " High field extrapolation: ";
1502 if (m_extrHighDiffusion == 0)
1503 std::cout << " constant\n";
1504 else if (m_extrHighDiffusion == 1)
1505 std::cout << " linear\n";
1506 else if (m_extrHighDiffusion == 2)
1507 std::cout << " exponential\n";
1508 else
1509 std::cout << " unknown\n";
1510 std::cout << " Interpolation order: " << m_intpDiffusion << "\n";
1511 }
1512 if (!tabElectronTownsend.empty()) {
1513 std::cout << " Townsend coefficient\n";
1514 std::cout << " Low field extrapolation: ";
1515 if (m_extrLowTownsend == 0)
1516 std::cout << " constant\n";
1517 else if (m_extrLowTownsend == 1)
1518 std::cout << " linear\n";
1519 else if (m_extrLowTownsend == 2)
1520 std::cout << " exponential\n";
1521 else
1522 std::cout << " unknown\n";
1523 std::cout << " High field extrapolation: ";
1524 if (m_extrHighTownsend == 0)
1525 std::cout << " constant\n";
1526 else if (m_extrHighTownsend == 1)
1527 std::cout << " linear\n";
1528 else if (m_extrHighTownsend == 2)
1529 std::cout << " exponential\n";
1530 else
1531 std::cout << " unknown\n";
1532 std::cout << " Interpolation order: " << m_intpTownsend << "\n";
1533 }
1535 std::cout << " Attachment coefficient\n";
1536 std::cout << " Low field extrapolation: ";
1537 if (m_extrLowAttachment == 0)
1538 std::cout << " constant\n";
1539 else if (m_extrLowAttachment == 1)
1540 std::cout << " linear\n";
1541 else if (m_extrLowAttachment == 2)
1542 std::cout << " exponential\n";
1543 else
1544 std::cout << " unknown\n";
1545 std::cout << " High field extrapolation: ";
1546 if (m_extrHighAttachment == 0)
1547 std::cout << " constant\n";
1548 else if (m_extrHighAttachment == 1)
1549 std::cout << " linear\n";
1550 else if (m_extrHighAttachment == 2)
1551 std::cout << " exponential\n";
1552 else
1553 std::cout << " unknown\n";
1554 std::cout << " Interpolation order: " << m_intpAttachment << "\n";
1555 }
1557 std::cout << " Lorentz Angle\n";
1558 std::cout << " Low field extrapolation: ";
1559 if (m_extrLowLorentzAngle == 0)
1560 std::cout << " constant\n";
1561 else if (m_extrLowLorentzAngle == 1)
1562 std::cout << " linear\n";
1563 else if (m_extrLowLorentzAngle == 2)
1564 std::cout << " exponential\n";
1565 else
1566 std::cout << " unknown\n";
1567 std::cout << " High field extrapolation: ";
1568 if (m_extrHighLorentzAngle == 0)
1569 std::cout << " constant\n";
1570 else if (m_extrHighLorentzAngle == 1)
1571 std::cout << " linear\n";
1572 else if (m_extrHighLorentzAngle == 2)
1573 std::cout << " exponential\n";
1574 else
1575 std::cout << " unknown\n";
1576 std::cout << " Interpolation order: " << m_intpLorentzAngle << "\n";
1577 }
1578 if (m_hasExcRates) {
1579 std::cout << " Excitation rates\n";
1580 std::cout << " Low field extrapolation: ";
1581 if (m_extrLowExcRates == 0)
1582 std::cout << " constant\n";
1583 else if (m_extrLowExcRates == 1)
1584 std::cout << " linear\n";
1585 else if (m_extrLowExcRates == 2)
1586 std::cout << " exponential\n";
1587 else
1588 std::cout << " unknown\n";
1589 std::cout << " High field extrapolation: ";
1590 if (m_extrHighExcRates == 0)
1591 std::cout << " constant\n";
1592 else if (m_extrHighExcRates == 1)
1593 std::cout << " linear\n";
1594 else if (m_extrHighExcRates == 2)
1595 std::cout << " exponential\n";
1596 else
1597 std::cout << " unknown\n";
1598 std::cout << " Interpolation order: " << m_intpExcRates << "\n";
1599 }
1600 if (m_hasIonRates) {
1601 std::cout << " Ionisation rates\n";
1602 std::cout << " Low field extrapolation: ";
1603 if (m_extrLowIonRates == 0)
1604 std::cout << " constant\n";
1605 else if (m_extrLowIonRates == 1)
1606 std::cout << " linear\n";
1607 else if (m_extrLowIonRates == 2)
1608 std::cout << " exponential\n";
1609 else
1610 std::cout << " unknown\n";
1611 std::cout << " High field extrapolation: ";
1612 if (m_extrHighIonRates == 0)
1613 std::cout << " constant\n";
1614 else if (m_extrHighIonRates == 1)
1615 std::cout << " linear\n";
1616 else if (m_extrHighIonRates == 2)
1617 std::cout << " exponential\n";
1618 else
1619 std::cout << " unknown\n";
1620 std::cout << " Interpolation order: " << m_intpIonRates << "\n";
1621 }
1626 std::cout << " none\n";
1627 }
1628
1629 std::cout << " Available ion transport data:\n";
1630 if (m_hasIonMobility) {
1631 std::cout << " Mobility\n";
1632 std::cout << " Low field extrapolation: ";
1633 if (m_extrLowMobility == 0)
1634 std::cout << " constant\n";
1635 else if (m_extrLowMobility == 1)
1636 std::cout << " linear\n";
1637 else if (m_extrLowMobility == 2)
1638 std::cout << " exponential\n";
1639 else
1640 std::cout << " unknown\n";
1641 std::cout << " High field extrapolation: ";
1642 if (m_extrHighMobility == 0)
1643 std::cout << " constant\n";
1644 else if (m_extrHighMobility == 1)
1645 std::cout << " linear\n";
1646 else if (m_extrHighMobility == 2)
1647 std::cout << " exponential\n";
1648 else
1649 std::cout << " unknown\n";
1650 std::cout << " Interpolation order: " << m_intpMobility << "\n";
1651 }
1652 if (m_hasIonDiffLong) {
1653 std::cout << " Longitudinal diffusion coefficient\n";
1654 }
1655 if (m_hasIonDiffTrans) {
1656 std::cout << " Transverse diffusion coefficient\n";
1657 }
1659 std::cout << " Low field extrapolation: ";
1660 if (m_extrLowDiffusion == 0)
1661 std::cout << " constant\n";
1662 else if (m_extrLowDiffusion == 1)
1663 std::cout << " linear\n";
1664 else if (m_extrLowDiffusion == 2)
1665 std::cout << " exponential\n";
1666 else
1667 std::cout << " unknown\n";
1668 std::cout << " High field extrapolation: ";
1669 if (m_extrHighDiffusion == 0)
1670 std::cout << " constant\n";
1671 else if (m_extrHighDiffusion == 1)
1672 std::cout << " linear\n";
1673 else if (m_extrHighDiffusion == 2)
1674 std::cout << " exponential\n";
1675 else
1676 std::cout << " unknown\n";
1677 std::cout << " Interpolation order: " << m_intpDiffusion << "\n";
1678 }
1680 std::cout << " Dissociation coefficient\n";
1681 std::cout << " Low field extrapolation: ";
1682 if (m_extrLowDissociation == 0)
1683 std::cout << " constant\n";
1684 else if (m_extrLowDissociation == 1)
1685 std::cout << " linear\n";
1686 else if (m_extrLowDissociation == 2)
1687 std::cout << " exponential\n";
1688 else
1689 std::cout << " unknown\n";
1690 std::cout << " High field extrapolation: ";
1691 if (m_extrHighDissociation == 0)
1692 std::cout << " constant\n";
1693 else if (m_extrHighDissociation == 1)
1694 std::cout << " linear\n";
1695 else if (m_extrHighDissociation == 2)
1696 std::cout << " exponential\n";
1697 else
1698 std::cout << " unknown\n";
1699 std::cout << " Interpolation order: " << m_intpDissociation << "\n";
1700 }
1703 std::cout << " none\n";
1704 }
1705}

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

◆ ScaleAttachment()

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

Reimplemented from Garfield::Medium.

Definition at line 76 of file MediumGas.hh.

76 {
77 return eta * m_pressure / m_pressureTable;
78 }

◆ ScaleDiffusion()

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

Reimplemented from Garfield::Medium.

Definition at line 67 of file MediumGas.hh.

67 {
68 return d * sqrt(m_pressureTable / m_pressure);
69 }

◆ ScaleDiffusionTensor()

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

Reimplemented from Garfield::Medium.

Definition at line 70 of file MediumGas.hh.

70 {
71 return d * m_pressureTable / m_pressure;
72 }

◆ ScaleElectricField()

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

Reimplemented from Garfield::Medium.

Definition at line 61 of file MediumGas.hh.

61 {
62 return e * m_pressureTable / m_pressure;
63 }

◆ ScaleLorentzAngle()

double Garfield::MediumGas::ScaleLorentzAngle ( const double  lor) const
inlinevirtual

Reimplemented from Garfield::Medium.

Definition at line 79 of file MediumGas.hh.

79 {
80 return lor * m_pressure / m_pressureTable;
81 }

◆ ScaleTownsend()

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

Reimplemented from Garfield::Medium.

Definition at line 73 of file MediumGas.hh.

73 {
75 }

◆ SetAtomicNumber()

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

Reimplemented from Garfield::Medium.

Definition at line 191 of file MediumGas.cc.

191 {
192
193 std::cerr << m_className << "::SetAtomicNumber:\n"
194 << " Effective Z cannot be changed directly to " << z << ".\n"
195 << " Use SetComposition to define the gas mixture.\n";
196}

◆ SetAtomicWeight()

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

Reimplemented from Garfield::Medium.

Definition at line 198 of file MediumGas.cc.

198 {
199
200 std::cerr << m_className << "::SetAtomicWeight:\n"
201 << " Effective A cannot be changed directly to " << a << ".\n"
202 << " Use SetComposition to define the gas mixture.\n";
203}

◆ 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 59 of file MediumGas.cc.

64 {
65
66 // Make a backup copy of the gas composition.
67 std::string gasOld[m_nMaxGases];
68 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
69 gasOld[i] = m_gas[i];
70 }
71 const unsigned int nComponentsOld = m_nComponents;
72 m_nComponents = 0;
73
74 std::string gases[6] = {gas1, gas2, gas3, gas4, gas5, gas6};
75 double fractions[6] = {f1, f2, f3, f4, f5, f6};
76 for (unsigned int i = 0; i < 6; ++i) {
77 // Find the gas name corresponding to the input string.
78 std::string gasname = "";
79 if (fractions[i] > 0. && GetGasName(gases[i], gasname)) {
80 m_gas[m_nComponents] = gasname;
81 m_fraction[m_nComponents] = fractions[i];
83 }
84 }
85
86 // Check if at least one valid ingredient was specified.
87 if (m_nComponents == 0) {
88 std::cerr << m_className << "::SetComposition:\n"
89 << " Error setting the composition.\n"
90 << " No valid ingredients were specified.\n";
91 return false;
92 }
93
94 // Establish the name.
95 m_name = "";
96 double sum = 0.;
97 for (unsigned int i = 0; i < m_nComponents; ++i) {
98 if (i > 0) m_name += "/";
99 m_name += m_gas[i];
100 sum += m_fraction[i];
101 }
102 // Normalise the fractions to one.
103 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
104 if (i < m_nComponents) {
105 m_fraction[i] /= sum;
106 } else {
107 m_fraction[i] = 0.;
108 }
109 }
110
111 // Set the atomic weight and number
112 for (unsigned int i = 0; i < m_nComponents; ++i) {
113 m_atWeight[i] = 0.;
114 m_atNum[i] = 0.;
116 }
117
118 // Print the composition.
119 std::cout << m_className << "::SetComposition:\n " << m_name;
120 if (m_nComponents > 1) {
121 std::cout << " (" << m_fraction[0] * 100;
122 for (unsigned int i = 1; i < m_nComponents; ++i) {
123 std::cout << "/" << m_fraction[i] * 100;
124 }
125 std::cout << ")";
126 }
127 std::cout << "\n";
128
129 // Force a recalculation of the collision rates.
130 m_isChanged = true;
131
132 // Copy the previous Penning transfer parameters.
133 double rPenningGasOld[m_nMaxGases];
134 double lambdaPenningGasOld[m_nMaxGases];
135 for (unsigned int i = 0; i < m_nMaxGases; ++i) {
136 rPenningGasOld[i] = m_rPenningGas[i];
137 lambdaPenningGasOld[i] = m_lambdaPenningGas[i];
138 m_rPenningGas[i] = 0.;
139 m_lambdaPenningGas[i] = 0.;
140 }
141 for (unsigned int i = 0; i < m_nComponents; ++i) {
142 for (unsigned int j = 0; j < nComponentsOld; ++j) {
143 if (m_gas[i] != gasOld[j]) continue;
144 if (rPenningGasOld[j] > 0.) {
145 m_rPenningGas[i] = rPenningGasOld[j];
146 m_lambdaPenningGas[i] = lambdaPenningGasOld[i];
147 std::cout << m_className << "::SetComposition:\n"
148 << " Using Penning transfer parameters for "
149 << m_gas[i] << " from previous mixture.\n"
150 << " r = " << m_rPenningGas[i] << "\n"
151 << " lambda = " << m_lambdaPenningGas[i] << " cm\n";
152 }
153 }
154 }
155 return true;
156}

Referenced by GarfieldPhysics::InitializePhysics().

◆ SetExtrapolationMethodExcitationRates()

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

Definition at line 1811 of file MediumGas.cc.

1812 {
1813
1814 unsigned int iExtr = 0;
1815 if (GetExtrapolationIndex(extrLow, iExtr)) {
1816 m_extrLowExcRates = iExtr;
1817 } else {
1818 std::cerr << m_className << "::SetExtrapolationMethodExcitationRates:\n";
1819 std::cerr << " Unknown extrapolation method (" << extrLow << ")\n";
1820 }
1821 if (GetExtrapolationIndex(extrHigh, iExtr)) {
1822 m_extrHighExcRates = iExtr;
1823 } else {
1824 std::cerr << m_className << "::SetExtrapolationMethodExcitationRates:\n";
1825 std::cerr << " Unknown extrapolation method (" << extrHigh << ")\n";
1826 }
1827}
bool GetExtrapolationIndex(std::string extrStr, unsigned int &extrNb)
Definition: Medium.cc:2453

◆ SetExtrapolationMethodIonisationRates()

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

Definition at line 1829 of file MediumGas.cc.

1830 {
1831
1832 unsigned int iExtr = 0;
1833 if (GetExtrapolationIndex(extrLow, iExtr)) {
1834 m_extrLowIonRates = iExtr;
1835 } else {
1836 std::cerr << m_className << "::SetExtrapolationMethodIonisationRates:\n";
1837 std::cerr << " Unknown extrapolation method (" << extrLow << ")\n";
1838 }
1839 if (GetExtrapolationIndex(extrHigh, iExtr)) {
1840 m_extrHighIonRates = iExtr;
1841 } else {
1842 std::cerr << m_className << "::SetExtrapolationMethodIonisationRates:\n";
1843 std::cerr << " Unknown extrapolation method (" << extrHigh << ")\n";
1844 }
1845}

◆ SetInterpolationMethodExcitationRates()

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

Definition at line 1847 of file MediumGas.cc.

1847 {
1848
1849 if (intrp > 0) {
1850 m_intpExcRates = intrp;
1851 }
1852}

◆ SetInterpolationMethodIonisationRates()

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

Definition at line 1854 of file MediumGas.cc.

1854 {
1855
1856 if (intrp > 0) {
1857 m_intpIonRates = intrp;
1858 }
1859}

◆ SetMassDensity()

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

Reimplemented from Garfield::Medium.

Definition at line 212 of file MediumGas.cc.

212 {
213
214 std::cerr << m_className << "::SetMassDensity:\n"
215 << " Density cannot directly be changed to " << rho << ".\n"
216 << " Use SetTemperature, SetPressure and SetComposition.\n";
217}

◆ SetNumberDensity()

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

Reimplemented from Garfield::Medium.

Definition at line 205 of file MediumGas.cc.

205 {
206
207 std::cerr << m_className << "::SetNumberDensity:\n"
208 << " Density cannot directly be changed to " << n << ".\n"
209 << " Use SetTemperature and SetPressure.\n";
210}

◆ UnScaleElectricField()

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

Reimplemented from Garfield::Medium.

Definition at line 64 of file MediumGas.hh.

64 {
65 return e * m_pressure / m_pressureTable;
66 }

◆ WriteGasFile()

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

Definition at line 986 of file MediumGas.cc.

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

Member Data Documentation

◆ m_atNum

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

Definition at line 93 of file MediumGas.hh.

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

◆ m_atWeight

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

Definition at line 92 of file MediumGas.hh.

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

◆ m_excitationList

std::vector<excListElement> Garfield::MediumGas::m_excitationList
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_fraction

◆ m_gas

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

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

◆ m_lambdaPenningGas

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

◆ m_lambdaPenningGlobal

double Garfield::MediumGas::m_lambdaPenningGlobal
protected

◆ m_nMaxGases

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

◆ m_pressureTable

◆ m_rPenningGas

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

◆ m_rPenningGlobal

double Garfield::MediumGas::m_rPenningGlobal
protected

◆ m_tabExcRates

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

◆ m_tabIonRates

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

◆ m_tabTownsendNoPenning

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

◆ m_temperatureTable

double Garfield::MediumGas::m_temperatureTable
protected

◆ m_usePenning


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