![]() |
Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
|
Interface to Magboltz (version 11). More...
#include <MediumMagboltz.hh>
Classes | |
struct | Deexcitation |
Public Member Functions | |
MediumMagboltz () | |
Default constructor. | |
MediumMagboltz (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.) | |
Constructor. | |
virtual | ~MediumMagboltz () |
Destructor. | |
bool | SetMaxElectronEnergy (const double e) |
Set the highest electron energy to be included in the table of scattering rates. | |
double | GetMaxElectronEnergy () const |
Get the highest electron energy in the table of scattering rates. | |
bool | SetMaxPhotonEnergy (const double e) |
Set the highest photon energy to be included in the table of scattering rates. | |
double | GetMaxPhotonEnergy () const |
Get the highest photon energy in the table of scattering rates. | |
void | EnableAnisotropicScattering (const bool on=true) |
Switch on/off anisotropic scattering (enabled by default) | |
void | SetSplittingFunctionOpalBeaty () |
Sample the secondary electron energy according to the Opal-Beaty model. | |
void | SetSplittingFunctionGreenSawada () |
Sample the secondary electron energy according to the Green-Sawada model. | |
void | SetSplittingFunctionFlat () |
Sample the secondary electron energy from a flat distribution. | |
void | EnableDeexcitation () |
Switch on (microscopic) de-excitation handling. | |
void | DisableDeexcitation () |
Switch off (microscopic) de-excitation handling. | |
void | EnableRadiationTrapping () |
Switch on discrete photoabsorption levels. | |
void | DisableRadiationTrapping () |
Switch off discrete photoabsorption levels. | |
void | SetLineWidth (const double n) |
Set the number of emission line widths within which to apply absorption by discrete lines. | |
bool | EnablePenningTransfer () override |
Switch on simulation of Penning transfers, using pre-implemented parameterisations of the transfer probability (if available). | |
bool | EnablePenningTransfer (const double r, const double lambda) override |
Switch on simulation of Penning transfers by means of transfer probabilities, for all excitation levels in the mixture. | |
bool | EnablePenningTransfer (const double r, const double lambda, std::string gasname) override |
Switch on simulation of Penning transfers by means of transfer probabilities, for all excitations of a given component. | |
void | DisablePenningTransfer () override |
Switch the simulation of Penning transfers off globally. | |
bool | DisablePenningTransfer (std::string gasname) override |
Switch the simulation of Penning transfers off for a given component. | |
void | EnableCrossSectionOutput (const bool on=true) |
Write the gas cross-section table to a file during the initialisation. | |
void | SetExcitationScaling (const double r, std::string gasname) |
Multiply all excitation cross-sections by a uniform scaling factor. | |
bool | Initialise (const bool verbose=false) |
Initialise the table of scattering rates (called internally when a collision rate is requested and the gas mixture or other parameters have changed). | |
void | PrintGas () override |
Print information about the present gas mixture and available data. | |
double | GetElectronNullCollisionRate (const int band) override |
double | GetElectronCollisionRate (const double e, const int band) override |
Get the (real) collision rate [ns-1] at a given electron energy e [eV]. | |
double | GetElectronCollisionRate (const double e, const unsigned int level, const int band) |
Get the collision rate [ns-1] for a specific level. | |
bool | ElectronCollision (const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< Secondary > &secondaries, int &band) override |
Sample the collision type. | |
void | ComputeDeexcitation (int iLevel, int &fLevel, std::vector< Secondary > &secondaries) |
double | GetPhotonCollisionRate (const double e) override |
bool | PhotonCollision (const double e, int &type, int &level, double &e1, double &ctheta, std::vector< Secondary > &secondaries) override |
void | ResetCollisionCounters () |
Reset the collision counters. | |
unsigned int | GetNumberOfElectronCollisions () const |
Get the total number of electron collisions. | |
unsigned int | GetNumberOfElectronCollisions (unsigned int &nElastic, unsigned int &nIonising, unsigned int &nAttachment, unsigned int &nInelastic, unsigned int &nExcitation, unsigned int &nSuperelastic) const |
Get the number of collisions broken down by cross-section type. | |
unsigned int | GetNumberOfLevels () |
Get the number of cross-section terms. | |
bool | GetLevel (const unsigned int i, int &ngas, int &type, std::string &descr, double &e) |
Get detailed information about a given cross-section term i. | |
bool | GetPenningTransfer (const unsigned int i, double &r, double &lambda) |
Get the Penning transfer probability and distance of a specific level. | |
unsigned int | GetNumberOfElectronCollisions (const unsigned int level) const |
Get the number of collisions for a specific cross-section term. | |
unsigned int | GetNumberOfPenningTransfers () const |
Get the number of Penning transfers that occured since the last reset. | |
unsigned int | GetNumberOfPhotonCollisions () const |
Get the total number of photon collisions. | |
unsigned int | GetNumberOfPhotonCollisions (unsigned int &nElastic, unsigned int &nIonising, unsigned int &nInelastic) const |
Get number of photon collisions by collision type. | |
void | EnableThermalMotion (const bool on=true) |
Take the thermal motion of the gas at the selected temperature into account in the calculations done by Magboltz. | |
void | EnableAutoEnergyLimit (const bool on=true) |
Let Magboltz determine the upper energy limit (this is the default) or use the energy limit specified using SetMaxElectronEnergy). | |
void | RunMagboltz (const double e, const double b, const double btheta, const int ncoll, bool verbose, double &vx, double &vy, double &vz, double &wv, double &wr, double &dl, double &dt, double &alpha, double &eta, double &riontof, double &ratttof, double &lor, double &vxerr, double &vyerr, double &vzerr, double &wverr, double &wrerr, double &dlerr, double &dterr, double &alphaerr, double &etaerr, double &riontoferr, double &ratttoferr, double &lorerr, double &alphatof, std::array< double, 6 > &difftens) |
Run Magboltz for a given electric field, magnetic field and angle. | |
void | GenerateGasTable (const int numCollisions=10, const bool verbose=true) |
Generate a new gas table (can later be saved to file) by running Magboltz for all electric fields, magnetic fields, and angles in the currently set grid. | |
void | PlotElectronCrossSections (const unsigned int i, TPad *pad) |
void | PlotElectronCollisionRates (TPad *pad) |
void | PlotElectronInverseMeanFreePath (TPad *pad) |
double | CreateGPUTransferObject (MediumGPU *&med_gpu) override |
Create and initialise GPU Transfer class. | |
Public Member Functions inherited from Garfield::MediumGas | |
MediumGas () | |
Constructor. | |
virtual | ~MediumGas () |
Destructor. | |
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.) |
Set the gas mixture. | |
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) const |
Retrieve the gas mixture. | |
bool | LoadGasFile (const std::string &filename, const bool quiet=false) |
Read table of gas properties (transport parameters) from file. | |
bool | WriteGasFile (const std::string &filename) |
Save the present table of gas properties (transport parameters) to a file. | |
bool | MergeGasFile (const std::string &filename, const bool replaceOld) |
Read table of gas properties from and merge with the existing dataset. | |
bool | GetPenningTransfer (const std::string &gasname, double &r, double &lambda) |
Retrieve the Penning transfer probability and distance for a specific component. | |
bool | LoadIonMobility (const std::string &filename, const bool quiet=false) |
Read a table of (positive) ion mobilities vs. electric field from file. | |
bool | LoadNegativeIonMobility (const std::string &filename, const bool quiet=false) |
Read a table of negative ion mobilities vs. electric field from file. | |
bool | AdjustTownsendCoefficient () |
Adjust the Townsend coefficient using the excitation and ionisation rates stored in the gas table and the Penning transfer probabilities. | |
size_t | GetNumberOfIonisationLevels () const |
Return the number of ionisation levels in the table. | |
size_t | GetNumberOfExcitationLevels () const |
Return the number of excitation levels in the table. | |
void | GetIonisationLevel (const size_t level, std::string &label, double &energy) const |
Return the identifier and threshold of an ionisation level. | |
void | GetExcitationLevel (const size_t level, std::string &label, double &energy) const |
Return the identifier and energy of an excitation level. | |
bool | GetElectronIonisationRate (const size_t level, const size_t ie, const size_t ib, const size_t ia, double &f) const |
Get an entry in the table of ionisation rates. | |
bool | GetElectronExcitationRate (const size_t level, const size_t ie, const size_t ib, const size_t ia, double &f) const |
Get an entry in the table of excitation rates. | |
bool | IsGas () const override |
void | GetComponent (const unsigned int i, std::string &label, double &f) override |
void | SetAtomicNumber (const double z) override |
double | GetAtomicNumber () const override |
void | SetAtomicWeight (const double a) override |
double | GetAtomicWeight () const override |
void | SetNumberDensity (const double n) override |
double | GetNumberDensity () const override |
void | SetMassDensity (const double rho) override |
double | GetMassDensity () const override |
void | ResetTables () override |
void | SetExtrapolationMethodExcitationRates (const std::string &low, const std::string &high) |
void | SetExtrapolationMethodIonisationRates (const std::string &low, const std::string &high) |
void | SetInterpolationMethodExcitationRates (const unsigned int intrp) |
void | SetInterpolationMethodIonisationRates (const unsigned int intrp) |
double | ScaleElectricField (const double e) const override |
double | UnScaleElectricField (const double e) const override |
double | ScaleDiffusion (const double d) const override |
double | ScaleDiffusionTensor (const double d) const override |
double | ScaleTownsend (const double alpha) const override |
double | ScaleAttachment (const double eta) const override |
double | ScaleLorentzAngle (const double lor) const override |
bool | GetPhotoAbsorptionCrossSection (const double e, double &sigma, const unsigned int i) override |
Static Public Member Functions | |
static int | GetGasNumberMagboltz (const std::string &input) |
Static Public Member Functions inherited from Garfield::MediumGas | |
static void | PrintGases () |
Print a list of all available gases. |
Private Member Functions | |
bool | Update (const bool verbose=false) |
bool | Mixer (const bool verbose=false) |
void | SetupGreenSawada () |
void | GetExcitationIonisationLevels () |
void | ComputeDeexcitationTable (const bool verbose) |
void | AddPenningDeexcitation (Deexcitation &dxc, const double rate, const double pPenning) |
double | RateConstantWK (const double energy, const double osc, const double pacs, const int igas1, const int igas2) const |
double | RateConstantHardSphere (const double r1, const double r2, const int igas1, const int igas2) const |
double | CalcDiscreteLineCf (const Deexcitation &dxc, double e, double cfOth) const |
void | ComputeDeexcitationInternal (int iLevel, int &fLevel, std::vector< Secondary > &secondaries) |
bool | ComputePhotonCollisionTable (const bool verbose) |
Private Attributes | |
std::mutex | m_mutex |
Mutex. | |
bool | m_useGasMotion = false |
Simulate thermal motion of the gas or not (when running Magboltz). | |
bool | m_autoEnergyLimit = true |
Automatic calculation of the energy limit by Magboltz or not. | |
double | m_eMax |
Max. electron energy in the collision rate tables. | |
double | m_eStep |
Energy spacing in the linear part of the collision rate tables. | |
double | m_eStepInv |
Inverse energy spacing. | |
double | m_eHigh |
double | m_eHighLog |
double | m_lnStep |
bool | m_useCsOutput = false |
Flag enabling/disabling output of cross-section table to file. | |
unsigned int | m_nTerms = 0 |
Number of different cross-section types in the current gas mixture. | |
std::array< double, m_nMaxGases > | m_mgas |
Mass. | |
std::array< double, m_nMaxGases > | m_rgas |
Recoil energy parameter. | |
std::array< double, m_nMaxGases > | m_s2 |
std::array< double, Magboltz::nMaxLevels > | m_wOpalBeaty |
Opal-Beaty-Peterson splitting parameter [eV]. | |
std::array< std::array< double, 5 >, m_nMaxGases > | m_parGreenSawada |
Green-Sawada splitting parameters [eV] (Γs, Γb, Ts, Ta, Tb). | |
std::array< bool, m_nMaxGases > | m_hasGreenSawada |
bool | m_useOpalBeaty = true |
Sample secondary electron energies using Opal-Beaty parameterisation. | |
bool | m_useGreenSawada = false |
Sample secondary electron energies using Green-Sawada parameterisation. | |
std::array< double, Magboltz::nMaxLevels > | m_energyLoss |
Energy loss. | |
std::array< int, Magboltz::nMaxLevels > | m_csType |
Cross-section type. | |
std::array< double, Magboltz::nMaxLevels > | m_yFluorescence |
Fluorescence yield. | |
std::array< unsigned int, Magboltz::nMaxLevels > | m_nAuger1 |
Number of Auger electrons produced in a collision. | |
std::array< unsigned int, Magboltz::nMaxLevels > | m_nAuger2 |
std::array< double, Magboltz::nMaxLevels > | m_eAuger1 |
Energy imparted to Auger electrons. | |
std::array< double, Magboltz::nMaxLevels > | m_eAuger2 |
std::array< unsigned int, Magboltz::nMaxLevels > | m_nFluorescence |
std::array< double, Magboltz::nMaxLevels > | m_eFluorescence |
std::vector< std::vector< double > > | m_scatPar |
std::vector< std::vector< double > > | m_scatCut |
std::vector< std::vector< double > > | m_scatParLog |
std::vector< std::vector< double > > | m_scatCutLog |
std::array< int, Magboltz::nMaxLevels > | m_scatModel |
std::vector< std::string > | m_description |
Level description. | |
std::vector< double > | m_cfTot |
std::vector< double > | m_cfTotLog |
std::vector< std::vector< double > > | m_cf |
std::vector< std::vector< double > > | m_cfLog |
bool | m_useAnisotropic = true |
double | m_cfNull = 0. |
Null-collision frequency. | |
std::array< unsigned int, nCsTypes > | m_nCollisions |
Collision counters 0: elastic 1: ionisation 2: attachment 3: inelastic 4: excitation 5: super-elastic. | |
std::vector< unsigned int > | m_nCollisionsDetailed |
Number of collisions for each cross-section term. | |
std::array< double, Magboltz::nMaxLevels > | m_rPenning |
Penning transfer probability (by level) | |
std::array< double, Magboltz::nMaxLevels > | m_lambdaPenning |
Mean distance of Penning ionisation (by level) | |
unsigned int | m_nPenning = 0 |
Number of Penning ionisations. | |
bool | m_useDeexcitation = false |
Flag enabling/disabling detailed simulation of de-excitation process. | |
bool | m_useRadTrap = true |
Flag enabling/disable radiation trapping (absorption of photons discrete excitation lines) | |
std::vector< Deexcitation > | m_deexcitations |
std::array< int, Magboltz::nMaxLevels > | m_iDeexcitation |
double | m_nAbsWidths = 1000. |
std::array< double, m_nMaxGases > | m_ionPot |
Ionisation potentials of each component. | |
double | m_minIonPot = -1. |
Minimum ionisation potential. | |
std::array< double, m_nMaxGases > | m_scaleExc |
double | m_eFinalGamma |
double | m_eStepGamma |
unsigned int | m_nPhotonTerms = 0 |
std::vector< double > | m_cfTotGamma |
std::vector< std::vector< double > > | m_cfGamma |
std::vector< int > | csTypeGamma |
std::array< unsigned int, nCsTypesGamma > | m_nPhotonCollisions |
Static Private Attributes | |
static constexpr int | nEnergyStepsLog = 1000 |
static constexpr int | nEnergyStepsGamma = 5000 |
static constexpr int | nCsTypes = 7 |
static constexpr int | nCsTypesGamma = 4 |
static const int | DxcTypeRad |
static const int | DxcTypeCollIon |
static const int | DxcTypeCollNonIon |
Additional Inherited Members | |
Protected Member Functions inherited from Garfield::MediumGas | |
bool | LoadMobility (const std::string &filename, const bool quiet, const bool negative) |
bool | ReadHeader (std::ifstream &gasfile, int &version, std::bitset< 20 > &gasok, bool &is3d, std::vector< double > &mixture, std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles, std::vector< ExcLevel > &excLevels, std::vector< IonLevel > &ionLevels) |
void | ReadFooter (std::ifstream &gasfile, std::array< unsigned int, 13 > &extrapH, std::array< unsigned int, 13 > &extrapL, std::array< unsigned int, 13 > &interp, unsigned int &thrAlp, unsigned int &thrAtt, unsigned int &thrDis, double &ionDiffL, double &ionDiffT, double &pgas, double &tgas) |
void | ReadRecord3D (std::ifstream &gasfile, double &ve, double &vb, double &vx, double &dl, double &dt, double &alpha, double &alpha0, double &eta, double &mu, double &lor, double &dis, std::array< double, 6 > &dif, std::vector< double > &rexc, std::vector< double > &rion, std::bitset< 20 > gasok) |
void | ReadRecord1D (std::ifstream &gasfile, double &ve, double &vb, double &vx, double &wv, double &wr, double &dl, double &dt, double &alpha, double &alpha0, double &eta, double &riontof, double &ratttof, double &mu, double &lor, double &dis, std::array< double, 6 > &dif, std::vector< double > &rexc, std::vector< double > &rion, std::bitset< 20 > gasok) |
void | InsertE (const int ie, const int ne, const int nb, const int na) |
void | InsertB (const int ib, const int ne, const int nb, const int na) |
void | InsertA (const int ia, const int ne, const int nb, const int na) |
void | ZeroRowE (const int ie, const int nb, const int na) |
void | ZeroRowB (const int ib, const int ne, const int na) |
void | ZeroRowA (const int ia, const int ne, const int nb) |
bool | GetMixture (const std::vector< double > &mixture, const int version, std::vector< std::string > &gasnames, std::vector< double > &percentages) const |
void | GetGasBits (std::bitset< 20 > &gasok) const |
Static Protected Member Functions inherited from Garfield::MediumGas | |
static bool | GetGasInfo (const std::string &gasname, double &a, double &z, double &w, double &f) |
static std::string | GetGasName (const int gasnumber, const int version) |
static std::string | GetGasName (std::string input) |
static int | GetGasNumberGasFile (const std::string &input) |
static const std::vector< std::string > | GetAliases (const std::string &gas) |
Protected Attributes inherited from Garfield::MediumGas | |
std::array< std::string, m_nMaxGases > | m_gas |
std::array< double, m_nMaxGases > | m_fraction |
std::array< double, m_nMaxGases > | m_atWeight |
std::array< double, m_nMaxGases > | m_atNum |
bool | m_usePenning = false |
double | m_rPenningGlobal = 0. |
double | m_lambdaPenningGlobal = 0. |
std::array< double, m_nMaxGases > | m_rPenningGas |
std::array< double, m_nMaxGases > | m_lambdaPenningGas |
double | m_pressureTable |
double | m_temperatureTable |
std::vector< std::vector< std::vector< double > > > | m_eAlp0 |
std::vector< std::vector< std::vector< std::vector< double > > > > | m_excRates |
std::vector< std::vector< std::vector< std::vector< double > > > > | m_ionRates |
std::vector< ExcLevel > | m_excLevels |
std::vector< IonLevel > | m_ionLevels |
std::pair< unsigned int, unsigned int > | m_extrExc = {0, 1} |
std::pair< unsigned int, unsigned int > | m_extrIon = {0, 1} |
unsigned int | m_intpExc = 2 |
unsigned int | m_intpIon = 2 |
Static Protected Attributes inherited from Garfield::MediumGas | |
static constexpr unsigned int | m_nMaxGases = 6 |
Interface to Magboltz (version 11).
Definition at line 22 of file MediumMagboltz.hh.
Garfield::MediumMagboltz::MediumMagboltz | ( | ) |
Default constructor.
Garfield::MediumMagboltz::MediumMagboltz | ( | 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. ) |
Constructor.
|
inlinevirtual |
|
inlineprivate |
Definition at line 449 of file MediumMagboltz.hh.
|
private |
void Garfield::MediumMagboltz::ComputeDeexcitation | ( | int | iLevel, |
int & | fLevel, | ||
std::vector< Secondary > & | secondaries ) |
|
private |
|
private |
|
private |
|
overridevirtual |
Create and initialise GPU Transfer class.
Reimplemented from Garfield::MediumGas.
|
inline |
Switch off (microscopic) de-excitation handling.
Definition at line 64 of file MediumMagboltz.hh.
|
overridevirtual |
Switch the simulation of Penning transfers off globally.
Reimplemented from Garfield::MediumGas.
|
overridevirtual |
Switch the simulation of Penning transfers off for a given component.
Reimplemented from Garfield::MediumGas.
|
inline |
Switch off discrete photoabsorption levels.
Definition at line 68 of file MediumMagboltz.hh.
|
override |
Sample the collision type.
|
inline |
Switch on/off anisotropic scattering (enabled by default)
Definition at line 49 of file MediumMagboltz.hh.
|
inline |
Let Magboltz determine the upper energy limit (this is the default) or use the energy limit specified using SetMaxElectronEnergy).
Definition at line 164 of file MediumMagboltz.hh.
|
inline |
Write the gas cross-section table to a file during the initialisation.
Definition at line 81 of file MediumMagboltz.hh.
void Garfield::MediumMagboltz::EnableDeexcitation | ( | ) |
Switch on (microscopic) de-excitation handling.
|
overridevirtual |
Switch on simulation of Penning transfers, using pre-implemented parameterisations of the transfer probability (if available).
Reimplemented from Garfield::MediumGas.
|
overridevirtual |
Switch on simulation of Penning transfers by means of transfer probabilities, for all excitation levels in the mixture.
r | transfer probability [0, 1] |
lambda | parameter for sampling the distance of the Penning electron with respect to the excitation. |
Reimplemented from Garfield::MediumGas.
|
overridevirtual |
Switch on simulation of Penning transfers by means of transfer probabilities, for all excitations of a given component.
Reimplemented from Garfield::MediumGas.
void Garfield::MediumMagboltz::EnableRadiationTrapping | ( | ) |
Switch on discrete photoabsorption levels.
|
inline |
Take the thermal motion of the gas at the selected temperature into account in the calculations done by Magboltz.
By the default, this feature is off (static gas at 0 K).
Definition at line 161 of file MediumMagboltz.hh.
void Garfield::MediumMagboltz::GenerateGasTable | ( | const int | numCollisions = 10, |
const bool | verbose = true ) |
Generate a new gas table (can later be saved to file) by running Magboltz for all electric fields, magnetic fields, and angles in the currently set grid.
|
override |
Get the (real) collision rate [ns-1] at a given electron energy e [eV].
double Garfield::MediumMagboltz::GetElectronCollisionRate | ( | const double | e, |
const unsigned int | level, | ||
const int | band ) |
Get the collision rate [ns-1] for a specific level.
|
override |
|
private |
|
static |
bool Garfield::MediumMagboltz::GetLevel | ( | const unsigned int | i, |
int & | ngas, | ||
int & | type, | ||
std::string & | descr, | ||
double & | e ) |
Get detailed information about a given cross-section term i.
|
inline |
Get the highest electron energy in the table of scattering rates.
Definition at line 40 of file MediumMagboltz.hh.
|
inline |
Get the highest photon energy in the table of scattering rates.
Definition at line 46 of file MediumMagboltz.hh.
unsigned int Garfield::MediumMagboltz::GetNumberOfElectronCollisions | ( | ) | const |
Get the total number of electron collisions.
unsigned int Garfield::MediumMagboltz::GetNumberOfElectronCollisions | ( | const unsigned int | level | ) | const |
Get the number of collisions for a specific cross-section term.
unsigned int Garfield::MediumMagboltz::GetNumberOfElectronCollisions | ( | unsigned int & | nElastic, |
unsigned int & | nIonising, | ||
unsigned int & | nAttachment, | ||
unsigned int & | nInelastic, | ||
unsigned int & | nExcitation, | ||
unsigned int & | nSuperelastic ) const |
Get the number of collisions broken down by cross-section type.
unsigned int Garfield::MediumMagboltz::GetNumberOfLevels | ( | ) |
Get the number of cross-section terms.
|
inline |
Get the number of Penning transfers that occured since the last reset.
Definition at line 149 of file MediumMagboltz.hh.
unsigned int Garfield::MediumMagboltz::GetNumberOfPhotonCollisions | ( | ) | const |
Get the total number of photon collisions.
unsigned int Garfield::MediumMagboltz::GetNumberOfPhotonCollisions | ( | unsigned int & | nElastic, |
unsigned int & | nIonising, | ||
unsigned int & | nInelastic ) const |
Get number of photon collisions by collision type.
bool Garfield::MediumMagboltz::GetPenningTransfer | ( | const unsigned int | i, |
double & | r, | ||
double & | lambda ) |
Get the Penning transfer probability and distance of a specific level.
|
override |
bool Garfield::MediumMagboltz::Initialise | ( | const bool | verbose = false | ) |
Initialise the table of scattering rates (called internally when a collision rate is requested and the gas mixture or other parameters have changed).
|
private |
|
override |
void Garfield::MediumMagboltz::PlotElectronCollisionRates | ( | TPad * | pad | ) |
void Garfield::MediumMagboltz::PlotElectronCrossSections | ( | const unsigned int | i, |
TPad * | pad ) |
void Garfield::MediumMagboltz::PlotElectronInverseMeanFreePath | ( | TPad * | pad | ) |
|
overridevirtual |
Print information about the present gas mixture and available data.
Reimplemented from Garfield::MediumGas.
|
private |
|
private |
void Garfield::MediumMagboltz::ResetCollisionCounters | ( | ) |
Reset the collision counters.
void Garfield::MediumMagboltz::RunMagboltz | ( | const double | e, |
const double | b, | ||
const double | btheta, | ||
const int | ncoll, | ||
bool | verbose, | ||
double & | vx, | ||
double & | vy, | ||
double & | vz, | ||
double & | wv, | ||
double & | wr, | ||
double & | dl, | ||
double & | dt, | ||
double & | alpha, | ||
double & | eta, | ||
double & | riontof, | ||
double & | ratttof, | ||
double & | lor, | ||
double & | vxerr, | ||
double & | vyerr, | ||
double & | vzerr, | ||
double & | wverr, | ||
double & | wrerr, | ||
double & | dlerr, | ||
double & | dterr, | ||
double & | alphaerr, | ||
double & | etaerr, | ||
double & | riontoferr, | ||
double & | ratttoferr, | ||
double & | lorerr, | ||
double & | alphatof, | ||
std::array< double, 6 > & | difftens ) |
Run Magboltz for a given electric field, magnetic field and angle.
[in] | e | electric field |
[in] | b | magnetic field |
[in] | btheta | angle between electric and magnetic field |
[in] | ncoll | number of collisions (in multiples of 107) to be simulated |
[in] | verbose | verbosity flag |
[out] | vx,vy,vz | drift velocity vector |
[out] | wv,wr | flux and bulk drift velocity |
[out] | dl,dt | diffusion cofficients |
[out] | alpha | Townsend cofficient |
[out] | eta | attachment cofficient |
[out] | riontof | TOF ionisation rate |
[out] | ratttof | TOF attachment rate |
[out] | lor | Lorentz angle |
[out] | vxerr,vyerr,vzerr | errors on drift velocity |
[out] | wverr,wrerr | errors on TOF drifts |
[out] | dlerr,dterr | errors on diffusion coefficients |
[out] | alphaerr,etaerr | errors on Townsend/attachment coefficients |
[out] | riontoferr,ratttoferr | errors on TOF rates |
[out] | lorerr | error on Lorentz angle |
[out] | alphatof | effective Townsend coefficient ![]() |
[out] | difftens | components of the diffusion tensor (zz, xx, yy, xz, yz, xy) |
void Garfield::MediumMagboltz::SetExcitationScaling | ( | const double | r, |
std::string | gasname ) |
Multiply all excitation cross-sections by a uniform scaling factor.
|
inline |
Set the number of emission line widths within which to apply absorption by discrete lines.
Definition at line 71 of file MediumMagboltz.hh.
bool Garfield::MediumMagboltz::SetMaxElectronEnergy | ( | const double | e | ) |
Set the highest electron energy to be included in the table of scattering rates.
bool Garfield::MediumMagboltz::SetMaxPhotonEnergy | ( | const double | e | ) |
Set the highest photon energy to be included in the table of scattering rates.
void Garfield::MediumMagboltz::SetSplittingFunctionFlat | ( | ) |
Sample the secondary electron energy from a flat distribution.
void Garfield::MediumMagboltz::SetSplittingFunctionGreenSawada | ( | ) |
Sample the secondary electron energy according to the Green-Sawada model.
void Garfield::MediumMagboltz::SetSplittingFunctionOpalBeaty | ( | ) |
Sample the secondary electron energy according to the Opal-Beaty model.
|
private |
|
private |
|
private |
Definition at line 434 of file MediumMagboltz.hh.
|
staticprivate |
Definition at line 225 of file MediumMagboltz.hh.
|
staticprivate |
Definition at line 226 of file MediumMagboltz.hh.
|
staticprivate |
Definition at line 224 of file MediumMagboltz.hh.
|
private |
Automatic calculation of the energy limit by Magboltz or not.
Definition at line 234 of file MediumMagboltz.hh.
|
private |
Definition at line 332 of file MediumMagboltz.hh.
|
private |
Definition at line 433 of file MediumMagboltz.hh.
|
private |
Definition at line 333 of file MediumMagboltz.hh.
|
private |
Null-collision frequency.
Definition at line 338 of file MediumMagboltz.hh.
|
private |
Definition at line 329 of file MediumMagboltz.hh.
|
private |
Definition at line 431 of file MediumMagboltz.hh.
|
private |
Definition at line 330 of file MediumMagboltz.hh.
|
private |
Cross-section type.
Definition at line 305 of file MediumMagboltz.hh.
|
private |
Definition at line 410 of file MediumMagboltz.hh.
|
private |
Level description.
Definition at line 326 of file MediumMagboltz.hh.
|
private |
Energy imparted to Auger electrons.
Definition at line 313 of file MediumMagboltz.hh.
|
private |
Definition at line 314 of file MediumMagboltz.hh.
|
private |
Definition at line 427 of file MediumMagboltz.hh.
|
private |
Definition at line 316 of file MediumMagboltz.hh.
|
private |
Definition at line 243 of file MediumMagboltz.hh.
|
private |
Definition at line 243 of file MediumMagboltz.hh.
|
private |
Max. electron energy in the collision rate tables.
Definition at line 238 of file MediumMagboltz.hh.
|
private |
Energy loss.
Definition at line 303 of file MediumMagboltz.hh.
|
private |
Energy spacing in the linear part of the collision rate tables.
Definition at line 240 of file MediumMagboltz.hh.
|
private |
Definition at line 427 of file MediumMagboltz.hh.
|
private |
Inverse energy spacing.
Definition at line 242 of file MediumMagboltz.hh.
|
private |
Definition at line 261 of file MediumMagboltz.hh.
|
private |
Definition at line 412 of file MediumMagboltz.hh.
|
private |
Ionisation potentials of each component.
Definition at line 419 of file MediumMagboltz.hh.
|
private |
Mean distance of Penning ionisation (by level)
Definition at line 371 of file MediumMagboltz.hh.
|
private |
Definition at line 244 of file MediumMagboltz.hh.
|
private |
Mass.
Definition at line 252 of file MediumMagboltz.hh.
|
private |
Minimum ionisation potential.
Definition at line 421 of file MediumMagboltz.hh.
|
private |
Mutex.
Definition at line 229 of file MediumMagboltz.hh.
|
private |
Definition at line 416 of file MediumMagboltz.hh.
|
private |
Number of Auger electrons produced in a collision.
Definition at line 310 of file MediumMagboltz.hh.
|
private |
Definition at line 311 of file MediumMagboltz.hh.
|
private |
Collision counters 0: elastic 1: ionisation 2: attachment 3: inelastic 4: excitation 5: super-elastic.
Definition at line 363 of file MediumMagboltz.hh.
|
private |
Number of collisions for each cross-section term.
Definition at line 365 of file MediumMagboltz.hh.
|
private |
Definition at line 315 of file MediumMagboltz.hh.
|
private |
Number of Penning ionisations.
Definition at line 373 of file MediumMagboltz.hh.
|
private |
Definition at line 440 of file MediumMagboltz.hh.
|
private |
Definition at line 429 of file MediumMagboltz.hh.
|
private |
Number of different cross-section types in the current gas mixture.
Definition at line 249 of file MediumMagboltz.hh.
|
private |
Green-Sawada splitting parameters [eV] (Γs, Γb, Ts, Ta, Tb).
Definition at line 260 of file MediumMagboltz.hh.
|
private |
Recoil energy parameter.
Definition at line 254 of file MediumMagboltz.hh.
|
private |
Penning transfer probability (by level)
Definition at line 369 of file MediumMagboltz.hh.
|
private |
Definition at line 255 of file MediumMagboltz.hh.
|
private |
Definition at line 424 of file MediumMagboltz.hh.
|
private |
Definition at line 320 of file MediumMagboltz.hh.
|
private |
Definition at line 322 of file MediumMagboltz.hh.
|
private |
Definition at line 323 of file MediumMagboltz.hh.
|
private |
Definition at line 319 of file MediumMagboltz.hh.
|
private |
Definition at line 321 of file MediumMagboltz.hh.
|
private |
Definition at line 336 of file MediumMagboltz.hh.
|
private |
Flag enabling/disabling output of cross-section table to file.
Definition at line 247 of file MediumMagboltz.hh.
|
private |
Flag enabling/disabling detailed simulation of de-excitation process.
Definition at line 377 of file MediumMagboltz.hh.
|
private |
Simulate thermal motion of the gas or not (when running Magboltz).
Definition at line 232 of file MediumMagboltz.hh.
|
private |
Sample secondary electron energies using Green-Sawada parameterisation.
Definition at line 266 of file MediumMagboltz.hh.
|
private |
Sample secondary electron energies using Opal-Beaty parameterisation.
Definition at line 264 of file MediumMagboltz.hh.
|
private |
Flag enabling/disable radiation trapping (absorption of photons discrete excitation lines)
Definition at line 380 of file MediumMagboltz.hh.
|
private |
Opal-Beaty-Peterson splitting parameter [eV].
Definition at line 257 of file MediumMagboltz.hh.
|
private |
Fluorescence yield.
Definition at line 308 of file MediumMagboltz.hh.
|
staticconstexprprivate |
Definition at line 219 of file MediumMagboltz.hh.
|
staticconstexprprivate |
Definition at line 222 of file MediumMagboltz.hh.
|
staticconstexprprivate |
Definition at line 218 of file MediumMagboltz.hh.
|
staticconstexprprivate |
Definition at line 217 of file MediumMagboltz.hh.