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

#include <MediumMagboltz.hh>

+ Inheritance diagram for Garfield::MediumMagboltz:

Public Member Functions

 MediumMagboltz ()
 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)
 
double GetMaxElectronEnergy () const
 Get the highest electron energy in the table of scattering rates.
 
bool SetMaxPhotonEnergy (const double e)
 
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.
 
bool EnablePenningTransfer () override
 
bool EnablePenningTransfer (const double r, const double lambda) override
 
bool EnablePenningTransfer (const double r, const double lambda, std::string gasname) override
 
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)
 
void PrintGas () override
 Print information about the present gas mixture and available data.
 
double GetElectronNullCollisionRate (const int band) override
 Get the overall null-collision rate [ns-1].
 
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< std::pair< Particle, double > > &secondaries, int &ndxc, int &band) override
 Sample the collision type.
 
void ComputeDeexcitation (int iLevel, int &fLevel)
 
unsigned int GetNumberOfDeexcitationProducts () const override
 
bool GetDeexcitationProduct (const unsigned int i, double &t, double &s, int &type, double &energy) const override
 
double GetPhotonCollisionRate (const double e) override
 
bool GetPhotonCollision (const double e, int &type, int &level, double &e1, double &ctheta, int &nsec, double &esec) 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)
 
void EnableAutoEnergyLimit (const bool on=true)
 
void RunMagboltz (const double e, const double b, const double btheta, const int ncoll, bool verbose, double &vx, double &vy, double &vz, double &dl, double &dt, double &alpha, double &eta, double &lor, double &vxerr, double &vyerr, double &vzerr, double &dlerr, double &dterr, double &alphaerr, double &etaerr, double &lorerr, double &alphatof, std::array< double, 6 > &difftens)
 
void GenerateGasTable (const int numCollisions=10, const bool verbose=true)
 
void PlotElectronCrossSections ()
 
- 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)
 
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 ()
 
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
 Is this medium a gas?
 
void GetComponent (const unsigned int i, std::string &label, double &f) override
 Get the name and fraction of a given component.
 
void SetAtomicNumber (const double z) override
 Set the effective atomic number.
 
double GetAtomicNumber () const override
 Get the effective atomic number.
 
void SetAtomicWeight (const double a) override
 Set the effective atomic weight.
 
double GetAtomicWeight () const override
 Get the effective atomic weight.
 
void SetNumberDensity (const double n) override
 Set the number density [cm-3].
 
double GetNumberDensity () const override
 Get the number density [cm-3].
 
void SetMassDensity (const double rho) override
 Set the mass density [g/cm3].
 
double GetMassDensity () const override
 Get the mass density [g/cm3].
 
void ResetTables () override
 Reset all tables of transport parameters.
 
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
 
- Public Member Functions inherited from Garfield::Medium
 Medium ()
 Constructor.
 
virtual ~Medium ()
 Destructor.
 
int GetId () const
 Return the id number of the class instance.
 
const std::string & GetName () const
 Get the medium name/identifier.
 
virtual bool IsSemiconductor () const
 Is this medium a semiconductor?
 
virtual bool IsConductor () const
 Is this medium a conductor?
 
void SetTemperature (const double t)
 Set the temperature [K].
 
double GetTemperature () const
 Get the temperature [K].
 
void SetPressure (const double p)
 
double GetPressure () const
 
void SetDielectricConstant (const double eps)
 Set the relative static dielectric constant.
 
double GetDielectricConstant () const
 Get the relative static dielectric constant.
 
unsigned int GetNumberOfComponents () const
 Get number of components of the medium.
 
virtual void EnableDrift (const bool on=true)
 Switch electron/ion/hole transport on/off.
 
virtual void EnablePrimaryIonisation (const bool on=true)
 Make the medium ionisable or non-ionisable.
 
bool IsDriftable () const
 Is charge carrier transport enabled in this medium?
 
bool IsMicroscopic () const
 Does the medium have electron scattering rates?
 
bool IsIonisable () const
 Is charge deposition by charged particles/photon enabled in this medium?
 
void SetW (const double w)
 Set the W value (average energy to produce an electron/ion or e/h pair).
 
double GetW () const
 Get the W value.
 
void SetFanoFactor (const double f)
 Set the Fano factor.
 
double GetFanoFactor () const
 Get the Fano factor.
 
void PlotVelocity (const std::string &carriers, TPad *pad)
 Plot the drift velocity as function of the electric field.
 
void PlotDiffusion (const std::string &carriers, TPad *pad)
 Plot the diffusion coefficients as function of the electric field.
 
void PlotTownsend (const std::string &carriers, TPad *pad)
 Plot the Townsend coefficient(s) as function of the electric field.
 
void PlotAttachment (const std::string &carriers, TPad *pad)
 Plot the attachment coefficient(s) as function of the electric field.
 
void PlotAlphaEta (const std::string &carriers, TPad *pad)
 Plot Townsend and attachment coefficients.
 
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)
 Drift velocity [cm / ns].
 
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)
 Longitudinal and transverse diffusion coefficients [cm1/2].
 
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)
 Ionisation coefficient [cm-1].
 
virtual bool ElectronAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 Attachment coefficient [cm-1].
 
virtual bool ElectronLorentzAngle (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
 Lorentz angle.
 
virtual double ElectronMobility ()
 Low-field mobility [cm2 V-1 ns-1].
 
virtual double GetElectronEnergy (const double px, const double py, const double pz, double &vx, double &vy, double &vz, const int band=0)
 Dispersion relation (energy vs. wave vector)
 
virtual void GetElectronMomentum (const double e, double &px, double &py, double &pz, int &band)
 
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)
 Drift velocity [cm / ns].
 
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)
 Longitudinal and transverse diffusion coefficients [cm1/2].
 
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])
 Diffusion tensor.
 
virtual bool HoleTownsend (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
 Ionisation coefficient [cm-1].
 
virtual bool HoleAttachment (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
 Attachment coefficient [cm-1].
 
virtual double HoleMobility ()
 Low-field mobility [cm2 V-1 ns-1].
 
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)
 Ion drift velocity [cm / ns].
 
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)
 Longitudinal and transverse diffusion coefficients [cm1/2].
 
virtual bool IonDissociation (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &diss)
 Dissociation coefficient.
 
virtual double IonMobility ()
 Low-field ion mobility [cm2 V-1 ns-1].
 
virtual bool NegativeIonVelocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
 Negative ion drift velocity [cm / ns].
 
virtual double NegativeIonMobility ()
 Low-field negative ion mobility [cm2 V-1 ns-1].
 
void SetFieldGrid (double emin, double emax, const size_t ne, bool logE, double bmin=0., double bmax=0., const size_t nb=1, double amin=HalfPi, double amax=HalfPi, const size_t na=1)
 Set the range of fields to be covered by the transport tables.
 
void SetFieldGrid (const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles)
 Set the fields and E-B angles to be used in the transport tables.
 
void GetFieldGrid (std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles)
 Get the fields and E-B angles used in the transport tables.
 
bool SetElectronVelocityE (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along E.
 
bool GetElectronVelocityE (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along E.
 
bool SetElectronVelocityExB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along ExB.
 
bool GetElectronVelocityExB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along ExB.
 
bool SetElectronVelocityB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along Btrans.
 
bool GetElectronVelocityB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along Btrans.
 
bool SetElectronLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetElectronLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetElectronTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetElectronTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetElectronTownsend (const size_t ie, const size_t ib, const size_t ia, const double alpha)
 Set an entry in the table of Townsend coefficients.
 
bool GetElectronTownsend (const size_t ie, const size_t ib, const size_t ia, double &alpha)
 Get an entry in the table of Townsend coefficients.
 
bool SetElectronAttachment (const size_t ie, const size_t ib, const size_t ia, const double eta)
 Set an entry in the table of attachment coefficients.
 
bool GetElectronAttachment (const size_t ie, const size_t ib, const size_t ia, double &eta)
 Get an entry in the table of attachment coefficients.
 
bool SetElectronLorentzAngle (const size_t ie, const size_t ib, const size_t ia, const double lor)
 Set an entry in the table of Lorentz angles.
 
bool GetElectronLorentzAngle (const size_t ie, const size_t ib, const size_t ia, double &lor)
 Get an entry in the table of Lorentz angles.
 
bool SetHoleVelocityE (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along E.
 
bool GetHoleVelocityE (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along E.
 
bool SetHoleVelocityExB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along ExB.
 
bool GetHoleVelocityExB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along ExB.
 
bool SetHoleVelocityB (const size_t ie, const size_t ib, const size_t ia, const double v)
 Set an entry in the table of drift speeds along Btrans.
 
bool GetHoleVelocityB (const size_t ie, const size_t ib, const size_t ia, double &v)
 Get an entry in the table of drift speeds along Btrans.
 
bool SetHoleLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetHoleLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetHoleTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetHoleTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetHoleTownsend (const size_t ie, const size_t ib, const size_t ia, const double alpha)
 Set an entry in the table of Townsend coefficients.
 
bool GetHoleTownsend (const size_t ie, const size_t ib, const size_t ia, double &alpha)
 Get an entry in the table of Townsend coefficients.
 
bool SetHoleAttachment (const size_t ie, const size_t ib, const size_t ia, const double eta)
 Set an entry in the table of attachment coefficients.
 
bool GetHoleAttachment (const size_t ie, const size_t ib, const size_t ia, double &eta)
 Get an entry in the table of attachment coefficients.
 
bool SetIonMobility (const std::vector< double > &fields, const std::vector< double > &mobilities, const bool negativeIons=false)
 
bool SetIonMobility (const size_t ie, const size_t ib, const size_t ia, const double mu)
 Set an entry in the table of ion mobilities.
 
bool GetIonMobility (const size_t ie, const size_t ib, const size_t ia, double &mu)
 Get an entry in the table of ion mobilities.
 
bool SetIonLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dl)
 Set an entry in the table of longitudinal diffusion coefficients.
 
bool GetIonLongitudinalDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dl)
 Get an entry in the table of longitudinal diffusion coefficients.
 
bool SetIonTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, const double dt)
 Set an entry in the table of transverse diffusion coefficients.
 
bool GetIonTransverseDiffusion (const size_t ie, const size_t ib, const size_t ia, double &dt)
 Get an entry in the table of transverse diffusion coefficients.
 
bool SetIonDissociation (const size_t ie, const size_t ib, const size_t ia, const double diss)
 Set an entry in the table of dissociation coefficients.
 
bool GetIonDissociation (const size_t ie, const size_t ib, const size_t ia, double &diss)
 Get an entry in the table of dissociation coefficients.
 
bool SetNegativeIonMobility (const size_t ie, const size_t ib, const size_t ia, const double mu)
 Set an entry in the table of negative ion mobilities.
 
bool GetNegativeIonMobility (const size_t ie, const size_t ib, const size_t ia, double &mu)
 Get an entry in the table of negative ion mobilities.
 
void ResetElectronVelocity ()
 
void ResetElectronDiffusion ()
 
void ResetElectronTownsend ()
 
void ResetElectronAttachment ()
 
void ResetElectronLorentzAngle ()
 
void ResetHoleVelocity ()
 
void ResetHoleDiffusion ()
 
void ResetHoleTownsend ()
 
void ResetHoleAttachment ()
 
void ResetIonMobility ()
 
void ResetIonDiffusion ()
 
void ResetIonDissociation ()
 
void ResetNegativeIonMobility ()
 
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)
 Set the degree of polynomial interpolation (usually 2).
 
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 ScaleVelocity (const double v) const
 
virtual double ScaleDissociation (const double diss) const
 
virtual bool GetOpticalDataRange (double &emin, double &emax, const unsigned int i=0)
 Get the energy range [eV] of the available optical data.
 
virtual bool GetDielectricFunction (const double e, double &eps1, double &eps2, const unsigned int i=0)
 Get the complex dielectric function at a given energy.
 
void EnableDebugging ()
 Switch on/off debugging messages.
 
void DisableDebugging ()
 

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.
 

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)
 
void ReadRecord1D (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)
 
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
 
- Protected Member Functions inherited from Garfield::Medium
bool Velocity (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &velE, const std::vector< std::vector< std::vector< double > > > &velB, const std::vector< std::vector< std::vector< double > > > &velX, const double q, double &vx, double &vy, double &vz) const
 
bool Diffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &difL, const std::vector< std::vector< std::vector< double > > > &difT, double &dl, double &dt) const
 
bool Diffusion (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< std::vector< double > > > > &diff, double cov[3][3]) const
 
bool Alpha (const double ex, const double ey, const double ez, const double bx, const double by, const double bz, const std::vector< std::vector< std::vector< double > > > &tab, unsigned int intp, const unsigned int thr, const std::pair< unsigned int, unsigned int > &extr, double &alpha) const
 
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
 
bool Interpolate (const double e, const double b, const double a, const std::vector< std::vector< std::vector< double > > > &table, double &y, const unsigned int intp, const std::pair< unsigned int, unsigned int > &extr) const
 
double Interpolate1D (const double e, const std::vector< double > &table, const std::vector< double > &fields, const unsigned int intpMeth, const std::pair< unsigned int, unsigned int > &extr) const
 
bool SetEntry (const size_t i, const size_t j, const size_t k, const std::string &fcn, std::vector< std::vector< std::vector< double > > > &tab, const double val)
 
bool GetEntry (const size_t i, const size_t j, const size_t k, const std::string &fcn, const std::vector< std::vector< std::vector< double > > > &tab, double &val) const
 
void SetExtrapolationMethod (const std::string &low, const std::string &high, std::pair< unsigned int, unsigned int > &extr, const std::string &fcn)
 
bool GetExtrapolationIndex (std::string str, unsigned int &nb) const
 
size_t SetThreshold (const std::vector< std::vector< std::vector< double > > > &tab) const
 
void Clone (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 std::pair< unsigned int, unsigned int > &extr, const double init, const std::string &label)
 
void Clone (std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const size_t n, const std::vector< double > &efields, const std::vector< double > &bfields, const std::vector< double > &angles, const unsigned int intp, const std::pair< unsigned int, unsigned int > &extr, const double init, const std::string &label)
 
void Init (const size_t nE, const size_t nB, const size_t nA, std::vector< std::vector< std::vector< double > > > &tab, const double val)
 
void Init (const size_t nE, const size_t nB, const size_t nA, const size_t nT, std::vector< std::vector< std::vector< std::vector< double > > > > &tab, const double val)
 
- 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)
 
- Static Protected Member Functions inherited from Garfield::Medium
static void Langevin (const double ex, const double ey, const double ez, double bx, double by, double bz, const double mu, double &vx, double &vy, double &vz)
 
static void Langevin (const double ex, const double ey, const double ez, double bx, double by, double bz, const double mu, const double muH, double &vx, double &vy, double &vz)
 
- Protected Attributes inherited from Garfield::MediumGas
std::array< std::string, m_nMaxGasesm_gas
 
std::array< double, m_nMaxGasesm_fraction
 
std::array< double, m_nMaxGasesm_atWeight
 
std::array< double, m_nMaxGasesm_atNum
 
bool m_usePenning = false
 
double m_rPenningGlobal = 0.
 
double m_lambdaPenningGlobal = 0.
 
std::array< double, m_nMaxGasesm_rPenningGas
 
std::array< double, m_nMaxGasesm_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< ExcLevelm_excLevels
 
std::vector< IonLevelm_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
 
- Protected Attributes inherited from Garfield::Medium
std::string m_className = "Medium"
 
int m_id
 
unsigned int m_nComponents = 1
 
std::string m_name = ""
 
double m_temperature = 293.15
 
double m_pressure = 760.
 
double m_epsilon = 1.
 
double m_z = 1.
 
double m_a = 0.
 
double m_density = 0.
 
double m_w = 0.
 
double m_fano = 0.
 
bool m_driftable = false
 
bool m_microscopic = false
 
bool m_ionisable = false
 
bool m_isChanged = true
 
bool m_debug = false
 
bool m_tab2d = false
 
std::vector< double > m_eFields
 
std::vector< double > m_bFields
 
std::vector< double > m_bAngles
 
std::vector< std::vector< std::vector< double > > > m_eVelE
 
std::vector< std::vector< std::vector< double > > > m_eVelX
 
std::vector< std::vector< std::vector< double > > > m_eVelB
 
std::vector< std::vector< std::vector< double > > > m_eDifL
 
std::vector< std::vector< std::vector< double > > > m_eDifT
 
std::vector< std::vector< std::vector< double > > > m_eAlp
 
std::vector< std::vector< std::vector< double > > > m_eAtt
 
std::vector< std::vector< std::vector< double > > > m_eLor
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
 
std::vector< std::vector< std::vector< double > > > m_hVelE
 
std::vector< std::vector< std::vector< double > > > m_hVelX
 
std::vector< std::vector< std::vector< double > > > m_hVelB
 
std::vector< std::vector< std::vector< double > > > m_hDifL
 
std::vector< std::vector< std::vector< double > > > m_hDifT
 
std::vector< std::vector< std::vector< double > > > m_hAlp
 
std::vector< std::vector< std::vector< double > > > m_hAtt
 
std::vector< std::vector< std::vector< std::vector< double > > > > m_hDifM
 
std::vector< std::vector< std::vector< double > > > m_iMob
 
std::vector< std::vector< std::vector< double > > > m_iDifL
 
std::vector< std::vector< std::vector< double > > > m_iDifT
 
std::vector< std::vector< std::vector< double > > > m_iDis
 
std::vector< std::vector< std::vector< double > > > m_nMob
 
unsigned int m_eThrAlp = 0
 
unsigned int m_eThrAtt = 0
 
unsigned int m_hThrAlp = 0
 
unsigned int m_hThrAtt = 0
 
unsigned int m_iThrDis = 0
 
std::pair< unsigned int, unsigned int > m_extrVel = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrDif = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrAlp = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrAtt = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrLor = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrMob = {0, 1}
 
std::pair< unsigned int, unsigned int > m_extrDis = {0, 1}
 
unsigned int m_intpVel = 2
 
unsigned int m_intpDif = 2
 
unsigned int m_intpAlp = 2
 
unsigned int m_intpAtt = 2
 
unsigned int m_intpLor = 2
 
unsigned int m_intpMob = 2
 
unsigned int m_intpDis = 2
 
- Static Protected Attributes inherited from Garfield::MediumGas
static constexpr unsigned int m_nMaxGases = 6
 
- Static Protected Attributes inherited from Garfield::Medium
static int m_idCounter = -1
 

Detailed Description

Interface to Magboltz (version 11).

Definition at line 15 of file MediumMagboltz.hh.

Constructor & Destructor Documentation

◆ MediumMagboltz() [1/2]

Garfield::MediumMagboltz::MediumMagboltz ( )

Default constructor.

Definition at line 79 of file MediumMagboltz.cc.

80 : MediumGas(),
81 m_eMax(40.),
82 m_eStep(m_eMax / Magboltz::nEnergySteps),
83 m_eStepInv(1. / m_eStep),
84 m_eHigh(400.),
85 m_eHighLog(log(m_eHigh)),
86 m_lnStep(1.),
87 m_eFinalGamma(20.),
88 m_eStepGamma(m_eFinalGamma / nEnergyStepsGamma) {
89 m_className = "MediumMagboltz";
90
91 // Set physical constants in Magboltz common blocks.
92 Magboltz::cnsts_.echarg = ElementaryCharge * 1.e-15;
93 Magboltz::cnsts_.emass = ElectronMassGramme;
94 Magboltz::cnsts_.amu = AtomicMassUnit;
95 Magboltz::cnsts_.pir2 = BohrRadius * BohrRadius * Pi;
96 Magboltz::inpt_.ary = RydbergEnergy;
97
98 // Set parameters in Magboltz common blocks.
101 // Select the scattering model.
102 Magboltz::inpt_.nAniso = 2;
103 // Max. energy [eV]
104 Magboltz::inpt_.efinal = m_eMax;
105 // Energy step size [eV]
106 Magboltz::inpt_.estep = m_eStep;
107 // Temperature and pressure
108 Magboltz::inpt_.akt = BoltzmannConstant * m_temperature;
109 Magboltz::inpt_.tempc = m_temperature - ZeroCelsius;
111 // Disable Penning transfer.
112 Magboltz::inpt_.ipen = 0;
113
114 m_description.assign(Magboltz::nMaxLevels,
115 std::string(Magboltz::nCharDescr, ' '));
116
117 m_cfTot.assign(Magboltz::nEnergySteps, 0.);
118 m_cfTotLog.assign(nEnergyStepsLog, 0.);
119 m_cf.assign(Magboltz::nEnergySteps,
120 std::vector<double>(Magboltz::nMaxLevels, 0.));
121 m_cfLog.assign(nEnergyStepsLog,
122 std::vector<double>(Magboltz::nMaxLevels, 0.));
123
124 m_isChanged = true;
125
126 m_driftable = true;
127 m_ionisable = true;
128 m_microscopic = true;
129
130 m_scaleExc.fill(1.);
131}
MediumGas()
Constructor.
Definition MediumGas.cc:115
double m_pressure
Definition Medium.hh:542
unsigned int m_nComponents
Definition Medium.hh:536
std::string m_className
Definition Medium.hh:529
double m_temperature
Definition Medium.hh:540
constexpr unsigned int nEnergySteps
struct Garfield::Magboltz::@272077330053140267135131315064054111215365003061 inpt_
constexpr unsigned int nCharDescr
struct Garfield::Magboltz::@057270266311332313123374024066154010246130331023 cnsts_
constexpr unsigned int nMaxLevels

Referenced by MediumMagboltz().

◆ MediumMagboltz() [2/2]

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.

Definition at line 132 of file MediumMagboltz.cc.

138 : MediumMagboltz() {
139 SetComposition(gas1, f1, gas2, f2, gas3, f3,
140 gas4, f4, gas5, f5, gas6, f6);
141}
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.
Definition MediumGas.cc:138
MediumMagboltz()
Default constructor.

◆ ~MediumMagboltz()

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

Destructor.

Definition at line 27 of file MediumMagboltz.hh.

27{}

Member Function Documentation

◆ ComputeDeexcitation()

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

Definition at line 2403 of file MediumMagboltz.cc.

2403 {
2404 if (!m_useDeexcitation) {
2405 std::cerr << m_className << "::ComputeDeexcitation: Not enabled.\n";
2406 return;
2407 }
2408
2409 // Make sure that the tables are updated.
2410 if (!Update()) return;
2411
2412 if (iLevel < 0 || iLevel >= (int)m_nTerms) {
2413 std::cerr << m_className << "::ComputeDeexcitation: Index out of range.\n";
2414 return;
2415 }
2416
2417 iLevel = m_iDeexcitation[iLevel];
2418 if (iLevel < 0 || iLevel >= (int)m_deexcitations.size()) {
2419 std::cerr << m_className << "::ComputeDeexcitation:\n"
2420 << " Level is not deexcitable.\n";
2421 return;
2422 }
2423
2424 ComputeDeexcitationInternal(iLevel, fLevel);
2425 if (fLevel >= 0 && fLevel < (int)m_deexcitations.size()) {
2426 fLevel = m_deexcitations[fLevel].level;
2427 }
2428}

◆ DisableDeexcitation()

void Garfield::MediumMagboltz::DisableDeexcitation ( )
inline

Switch off (microscopic) de-excitation handling.

Definition at line 57 of file MediumMagboltz.hh.

57{ m_useDeexcitation = false; }

◆ DisablePenningTransfer() [1/2]

void Garfield::MediumMagboltz::DisablePenningTransfer ( )
overridevirtual

Switch the simulation of Penning transfers off globally.

Reimplemented from Garfield::MediumGas.

Definition at line 335 of file MediumMagboltz.cc.

335 {
336
338 m_rPenning.fill(0.);
339 m_lambdaPenning.fill(0.);
340
341 m_usePenning = false;
342}
virtual void DisablePenningTransfer()
Switch the simulation of Penning transfers off globally.

◆ DisablePenningTransfer() [2/2]

bool Garfield::MediumMagboltz::DisablePenningTransfer ( std::string gasname)
overridevirtual

Switch the simulation of Penning transfers off for a given component.

Reimplemented from Garfield::MediumGas.

Definition at line 344 of file MediumMagboltz.cc.

344 {
345
346 if (!MediumGas::DisablePenningTransfer(gasname)) return false;
347 // Get the "standard" name of this gas.
348 gasname = GetGasName(gasname);
349 if (gasname.empty()) return false;
350
351 // Look (again) for this gas in the present gas mixture.
352 int iGas = -1;
353 for (unsigned int i = 0; i < m_nComponents; ++i) {
354 if (m_gas[i] == gasname) {
355 iGas = i;
356 break;
357 }
358 }
359
360 if (iGas < 0) return false;
361
362 unsigned int nLevelsFound = 0;
363 for (unsigned int i = 0; i < m_nTerms; ++i) {
364 if (int(m_csType[i] / nCsTypes) == iGas) {
365 m_rPenning[i] = 0.;
366 m_lambdaPenning[i] = 0.;
367 } else {
368 if (m_csType[i] % nCsTypes == ElectronCollisionTypeExcitation &&
369 m_rPenning[i] > Small) {
370 ++nLevelsFound;
371 }
372 }
373 }
374
375 if (nLevelsFound == 0) {
376 // There are no more excitation levels with r > 0.
377 std::cout << m_className << "::DisablePenningTransfer:\n"
378 << " Penning transfer switched off for all excitations.\n";
379 m_usePenning = false;
380 }
381 return true;
382}
static std::string GetGasName(const int gasnumber, const int version)
std::array< std::string, m_nMaxGases > m_gas
Definition MediumGas.hh:164

◆ DisableRadiationTrapping()

void Garfield::MediumMagboltz::DisableRadiationTrapping ( )
inline

Switch off discrete photoabsorption levels.

Definition at line 61 of file MediumMagboltz.hh.

61{ m_useRadTrap = false; }

◆ ElectronCollision()

bool Garfield::MediumMagboltz::ElectronCollision ( const double e,
int & type,
int & level,
double & e1,
double & dx,
double & dy,
double & dz,
std::vector< std::pair< Particle, double > > & secondaries,
int & ndxc,
int & band )
overridevirtual

Sample the collision type.

Reimplemented from Garfield::Medium.

Definition at line 589 of file MediumMagboltz.cc.

592 {
593 band = 0;
594 ndxc = 0;
595 if (e <= 0.) {
596 std::cerr << m_className << "::ElectronCollision: Invalid energy.\n";
597 return false;
598 }
599 // Check if the electron energy is within the currently set range.
600 if (e > m_eMax) {
601 std::cerr << m_className << "::ElectronCollision:\n Requested energy ("
602 << e << " eV) exceeds current energy range.\n"
603 << " Increasing energy range to " << 1.05 * e << " eV.\n";
604 SetMaxElectronEnergy(1.05 * e);
605 }
606
607 // If necessary, update the collision rates table.
608 if (!Update()) return false;
609
610 double angCut = 1.;
611 double angPar = 0.5;
612
613 if (e < m_eHigh) {
614 // Linear binning
615 // Get the energy interval.
616 const int iE = int(e * m_eStepInv);
617
618 // Sample the scattering process.
619 const double r = RndmUniform();
620 level = 0;
621 if (r > m_cf[iE][0]) {
622 int iLow = 0;
623 int iUp = m_nTerms - 1;
624 while (iUp - iLow > 1) {
625 int iMid = (iLow + iUp) >> 1;
626 if (r < m_cf[iE][iMid]) {
627 iUp = iMid;
628 } else {
629 iLow = iMid;
630 }
631 }
632 level = iUp;
633 }
634 // Get the angular distribution parameters.
635 angCut = m_scatCut[iE][level];
636 angPar = m_scatPar[iE][level];
637 } else {
638 // Logarithmic binning
639 // Get the energy interval.
640 const int iE = std::min(std::max(int(log(e / m_eHigh) / m_lnStep), 0),
641 nEnergyStepsLog - 1);
642 // Sample the scattering process.
643 const double r = RndmUniform();
644 level = 0;
645 if (r > m_cfLog[iE][0]) {
646 int iLow = 0;
647 int iUp = m_nTerms - 1;
648 while (iUp - iLow > 1) {
649 int iMid = (iLow + iUp) >> 1;
650 if (r < m_cfLog[iE][iMid]) {
651 iUp = iMid;
652 } else {
653 iLow = iMid;
654 }
655 }
656 level = iUp;
657 }
658 // Get the angular distribution parameters.
659 angCut = m_scatCutLog[iE][level];
660 angPar = m_scatParLog[iE][level];
661 }
662
663 // Extract the collision type.
664 type = m_csType[level] % nCsTypes;
665 const int igas = int(m_csType[level] / nCsTypes);
666 // Increase the collision counters.
667 ++m_nCollisions[type];
668 ++m_nCollisionsDetailed[level];
669
670 // Get the energy loss for this process.
671 double loss = m_energyLoss[level];
672
673 if (type == ElectronCollisionTypeVirtual) return true;
674
675 if (type == ElectronCollisionTypeIonisation) {
676 // Sample the secondary electron energy according to
677 // the Opal-Beaty-Peterson parameterisation.
678 double esec = 0.;
679 if (e < loss) loss = e - 0.0001;
680 if (m_useOpalBeaty) {
681 // Get the splitting parameter.
682 const double w = m_wOpalBeaty[level];
683 esec = w * tan(RndmUniform() * atan(0.5 * (e - loss) / w));
684 // Rescaling (SST)
685 // esec = w * pow(esec / w, 0.9524);
686 } else if (m_useGreenSawada) {
687 const double gs = m_parGreenSawada[igas][0];
688 const double gb = m_parGreenSawada[igas][1];
689 const double w = gs * e / (e + gb);
690 const double ts = m_parGreenSawada[igas][2];
691 const double ta = m_parGreenSawada[igas][3];
692 const double tb = m_parGreenSawada[igas][4];
693 const double esec0 = ts - ta / (e + tb);
694 const double r = RndmUniform();
695 esec = esec0 +
696 w * tan((r - 1.) * atan(esec0 / w) +
697 r * atan((0.5 * (e - loss) - esec0) / w));
698 } else {
699 esec = RndmUniform() * (e - loss);
700 }
701 if (esec <= 0) esec = Small;
702 loss += esec;
703 // Add the secondary electron.
704 secondaries.emplace_back(std::make_pair(Particle::Electron, esec));
705 // Add the ion.
706 secondaries.emplace_back(std::make_pair(Particle::Ion, 0.));
707 bool fluorescence = false;
708 if (m_yFluorescence[level] > Small) {
709 if (RndmUniform() < m_yFluorescence[level]) fluorescence = true;
710 }
711 // Add Auger and photo electrons (if any).
712 if (fluorescence) {
713 if (m_nAuger2[level] > 0) {
714 const double eav = m_eAuger2[level] / m_nAuger2[level];
715 for (unsigned int i = 0; i < m_nAuger2[level]; ++i) {
716 secondaries.emplace_back(std::make_pair(Particle::Electron, eav));
717 }
718 }
719 if (m_nFluorescence[level] > 0) {
720 const double eav = m_eFluorescence[level] / m_nFluorescence[level];
721 for (unsigned int i = 0; i < m_nFluorescence[level]; ++i) {
722 secondaries.emplace_back(std::make_pair(Particle::Electron, eav));
723 }
724 }
725 } else if (m_nAuger1[level] > 0) {
726 const double eav = m_eAuger1[level] / m_nAuger1[level];
727 for (unsigned int i = 0; i < m_nAuger1[level]; ++i) {
728 secondaries.emplace_back(std::make_pair(Particle::Electron, eav));
729 }
730 }
731 } else if (type == ElectronCollisionTypeExcitation) {
732 // Follow the de-excitation cascade (if switched on).
733 if (m_useDeexcitation && m_iDeexcitation[level] >= 0) {
734 int fLevel = 0;
735 ComputeDeexcitationInternal(m_iDeexcitation[level], fLevel);
736 ndxc = m_dxcProducts.size();
737 } else if (m_usePenning) {
738 m_dxcProducts.clear();
739 // Simplified treatment of Penning ionisation.
740 // If the energy threshold of this level exceeds the
741 // ionisation potential of one of the gases,
742 // create a new electron (with probability rPenning).
743 if (m_debug) {
744 std::cout << m_className << "::ElectronCollision:\n"
745 << " Level: " << level << "\n"
746 << " Ionization potential: " << m_minIonPot << "\n"
747 << " Excitation energy: " << loss * m_rgas[igas] << "\n"
748 << " Penning probability: " << m_rPenning[level] << "\n";
749 }
750 if (loss * m_rgas[igas] > m_minIonPot &&
751 RndmUniform() < m_rPenning[level]) {
752 // The energy of the secondary electron is assumed to be given by
753 // the difference of excitation and ionisation threshold.
754 double esec = loss * m_rgas[igas] - m_minIonPot;
755 if (esec <= 0) esec = Small;
756 // Add the secondary electron to the list.
757 dxcProd newDxcProd;
758 newDxcProd.t = 0.;
759 newDxcProd.s = 0.;
760 if (m_lambdaPenning[level] > Small) {
761 // Uniform distribution within a sphere of radius lambda
762 newDxcProd.s = m_lambdaPenning[level] * std::cbrt(RndmUniformPos());
763 }
764 newDxcProd.energy = esec;
765 newDxcProd.type = DxcProdTypeElectron;
766 m_dxcProducts.push_back(std::move(newDxcProd));
767 ndxc = 1;
768 ++m_nPenning;
769 }
770 }
771 }
772
773 if (e < loss) loss = e - 0.0001;
774
775 // Determine the scattering angle.
776 double ctheta0 = 1. - 2. * RndmUniform();
777 if (m_useAnisotropic) {
778 switch (m_scatModel[level]) {
779 case 0:
780 break;
781 case 1:
782 ctheta0 = 1. - RndmUniform() * angCut;
783 if (RndmUniform() > angPar) ctheta0 = -ctheta0;
784 break;
785 case 2:
786 ctheta0 = (ctheta0 + angPar) / (1. + angPar * ctheta0);
787 break;
788 default:
789 std::cerr << m_className << "::ElectronCollision:\n"
790 << " Unknown scattering model.\n"
791 << " Using isotropic distribution.\n";
792 break;
793 }
794 }
795
796 const double s1 = m_rgas[igas];
797 const double theta0 = acos(ctheta0);
798 const double arg = std::max(1. - s1 * loss / e, Small);
799 const double d = 1. - ctheta0 * sqrt(arg);
800
801 // Update the energy.
802 e1 = std::max(e * (1. - loss / (s1 * e) - 2. * d * m_s2[igas]), Small);
803 double q = std::min(sqrt((e / e1) * arg) / s1, 1.);
804 const double theta = asin(q * sin(theta0));
805 double ctheta = cos(theta);
806 if (ctheta0 < 0.) {
807 const double u = (s1 - 1.) * (s1 - 1.) / arg;
808 if (ctheta0 * ctheta0 > u) ctheta = -ctheta;
809 }
810 const double stheta = sin(theta);
811 // Calculate the direction after the collision.
812 dz = std::min(dz, 1.);
813 const double argZ = sqrt(dx * dx + dy * dy);
814
815 // Azimuth is chosen at random.
816 const double phi = TwoPi * RndmUniform();
817 const double cphi = cos(phi);
818 const double sphi = sin(phi);
819 if (argZ > 0.) {
820 const double a = stheta / argZ;
821 const double dz1 = dz * ctheta + argZ * stheta * sphi;
822 const double dy1 = dy * ctheta + a * (dx * cphi - dy * dz * sphi);
823 const double dx1 = dx * ctheta - a * (dy * cphi + dx * dz * sphi);
824 dz = dz1;
825 dy = dy1;
826 dx = dx1;
827 } else {
828 dz = ctheta;
829 dx = cphi * stheta;
830 dy = sphi * stheta;
831 }
832 return true;
833}
bool SetMaxElectronEnergy(const double e)
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition Random.hh:14
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition Random.hh:17
DoubleAc cos(const DoubleAc &f)
Definition DoubleAc.cpp:432
DoubleAc acos(const DoubleAc &f)
Definition DoubleAc.cpp:490
DoubleAc sin(const DoubleAc &f)
Definition DoubleAc.cpp:384
DoubleAc asin(const DoubleAc &f)
Definition DoubleAc.cpp:470
DoubleAc sqrt(const DoubleAc &f)
Definition DoubleAc.cpp:314

◆ EnableAnisotropicScattering()

void Garfield::MediumMagboltz::EnableAnisotropicScattering ( const bool on = true)
inline

Switch on/off anisotropic scattering (enabled by default)

Definition at line 42 of file MediumMagboltz.hh.

42 {
43 m_useAnisotropic = on;
44 m_isChanged = true;
45 }

◆ EnableAutoEnergyLimit()

void Garfield::MediumMagboltz::EnableAutoEnergyLimit ( const bool on = true)
inline

Let Magboltz determine the upper energy limit (this is the default) or use the energy limit specified using SetMaxElectronEnergy).

Definition at line 144 of file MediumMagboltz.hh.

144{ m_autoEnergyLimit = on; }

◆ EnableCrossSectionOutput()

void Garfield::MediumMagboltz::EnableCrossSectionOutput ( const bool on = true)
inline

Write the gas cross-section table to a file during the initialisation.

Definition at line 71 of file MediumMagboltz.hh.

71{ m_useCsOutput = on; }

◆ EnableDeexcitation()

void Garfield::MediumMagboltz::EnableDeexcitation ( )

Switch on (microscopic) de-excitation handling.

Definition at line 205 of file MediumMagboltz.cc.

205 {
206 if (m_usePenning) {
207 std::cout << m_className << "::EnableDeexcitation:\n"
208 << " Penning transfer will be switched off.\n";
209 }
210 // if (m_useRadTrap) {
211 // std::cout << " Radiation trapping is switched on.\n";
212 // } else {
213 // std::cout << " Radiation trapping is switched off.\n";
214 // }
215 m_usePenning = false;
216 m_useDeexcitation = true;
217 m_isChanged = true;
218 m_dxcProducts.clear();
219}

◆ EnablePenningTransfer() [1/3]

bool Garfield::MediumMagboltz::EnablePenningTransfer ( )
overridevirtual

Switch on simulation of Penning transfers, using pre-implemented parameterisations of the transfer probability (if available).

Reimplemented from Garfield::MediumGas.

Definition at line 231 of file MediumMagboltz.cc.

231 {
232 m_rPenning.fill(0.);
233 m_lambdaPenning.fill(0.);
234 if (!MediumGas::EnablePenningTransfer()) return false;
235
236 m_usePenning = true;
237 return true;
238}
virtual bool EnablePenningTransfer()

◆ EnablePenningTransfer() [2/3]

bool Garfield::MediumMagboltz::EnablePenningTransfer ( const double r,
const double lambda )
overridevirtual

Switch on simulation of Penning transfers by means of transfer probabilities, for all excitation levels in the mixture.

Parameters
rtransfer probability [0, 1]
lambdaparameter for sampling the distance of the Penning electron with respect to the excitation.

Reimplemented from Garfield::MediumGas.

Definition at line 240 of file MediumMagboltz.cc.

241 {
242
243 if (!MediumGas::EnablePenningTransfer(r, lambda)) return false;
244
245 m_rPenning.fill(0.);
246 m_lambdaPenning.fill(0.);
247
248 // Make sure that the collision rate table is up to date.
249 if (!Update()) return false;
250 unsigned int nLevelsFound = 0;
251 for (unsigned int i = 0; i < m_nTerms; ++i) {
252 if (m_csType[i] % nCsTypes == ElectronCollisionTypeExcitation) {
253 ++nLevelsFound;
254 }
255 m_rPenning[i] = m_rPenningGlobal;
256 m_lambdaPenning[i] = m_lambdaPenningGlobal;
257 }
258
259 if (nLevelsFound > 0) {
260 std::cout << m_className << "::EnablePenningTransfer:\n "
261 << "Updated Penning transfer parameters for " << nLevelsFound
262 << " excitation cross-sections.\n";
263 if (nLevelsFound != m_excLevels.size() && !m_excLevels.empty()) {
264 std::cerr << m_className << "::EnablePenningTransfer:\n Warning: "
265 << "mismatch between number of excitation cross-sections ("
266 << nLevelsFound << ")\n and number of excitation rates in "
267 << "the gas table (" << m_excLevels.size() << ").\n "
268 << "The gas table was probably calculated using a different "
269 << "version of Magboltz.\n";
270 }
271 } else {
272 std::cerr << m_className << "::EnablePenningTransfer:\n "
273 << "No excitation cross-sections in the present energy range.\n";
274 }
275
276 if (m_useDeexcitation) {
277 std::cout << m_className << "::EnablePenningTransfer:\n "
278 << "Deexcitation handling will be switched off.\n";
279 }
280 m_usePenning = true;
281 return true;
282}
double m_lambdaPenningGlobal
Definition MediumGas.hh:175
std::vector< ExcLevel > m_excLevels
Definition MediumGas.hh:201

◆ EnablePenningTransfer() [3/3]

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

Switch on simulation of Penning transfers by means of transfer probabilities, for all excitations of a given component.

Reimplemented from Garfield::MediumGas.

Definition at line 284 of file MediumMagboltz.cc.

285 {
286
287 if (!MediumGas::EnablePenningTransfer(r, lambda, gasname)) return false;
288
289 // Get (again) the "standard" name of this gas.
290 gasname = GetGasName(gasname);
291 if (gasname.empty()) return false;
292
293 // Look (again) for this gas in the present mixture.
294 int iGas = -1;
295 for (unsigned int i = 0; i < m_nComponents; ++i) {
296 if (m_gas[i] == gasname) {
297 iGas = i;
298 break;
299 }
300 }
301
302 // Make sure that the collision rate table is up to date.
303 if (!Update()) return false;
304 unsigned int nLevelsFound = 0;
305 for (unsigned int i = 0; i < m_nTerms; ++i) {
306 if (int(m_csType[i] / nCsTypes) != iGas) continue;
307 if (m_csType[i] % nCsTypes == ElectronCollisionTypeExcitation) {
308 ++nLevelsFound;
309 }
310 m_rPenning[i] = m_rPenningGas[iGas];
311 m_lambdaPenning[i] = m_lambdaPenningGas[iGas];
312 }
313
314 if (nLevelsFound > 0) {
315 std::cout << m_className << "::EnablePenningTransfer:\n";
316 if (m_lambdaPenningGas[iGas] > 0.) {
317 std::cout << " Penning transfer parameters for " << nLevelsFound
318 << " " << gasname << " excitation levels set to:\n"
319 << " r = " << m_rPenningGas[iGas] << ", lambda = "
320 << m_lambdaPenningGas[iGas] << " cm\n";
321 } else {
322 std::cout << " Penning transfer probability for " << nLevelsFound
323 << " " << gasname << " excitation levels set to r = "
324 << m_rPenningGas[iGas] << "\n";
325 }
326 } else {
327 std::cerr << m_className << "::EnablePenningTransfer:\n " << gasname
328 << " has no excitation levels in the present energy range.\n";
329 }
330
331 m_usePenning = true;
332 return true;
333}
std::array< double, m_nMaxGases > m_rPenningGas
Definition MediumGas.hh:177
std::array< double, m_nMaxGases > m_lambdaPenningGas
Definition MediumGas.hh:179

◆ EnableRadiationTrapping()

void Garfield::MediumMagboltz::EnableRadiationTrapping ( )

Switch on discrete photoabsorption levels.

Definition at line 221 of file MediumMagboltz.cc.

221 {
222 m_useRadTrap = true;
223 if (!m_useDeexcitation) {
224 std::cout << m_className << "::EnableRadiationTrapping:\n "
225 << "Radiation trapping is enabled but de-excitation is not.\n";
226 } else {
227 m_isChanged = true;
228 }
229}

◆ EnableThermalMotion()

void Garfield::MediumMagboltz::EnableThermalMotion ( const bool on = true)
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 141 of file MediumMagboltz.hh.

141{ m_useGasMotion = on; }

◆ GenerateGasTable()

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.

Definition at line 2787 of file MediumMagboltz.cc.

2787 {
2788 // Set the reference pressure and temperature.
2791
2792 // Initialize the parameter arrays.
2793 const unsigned int nEfields = m_eFields.size();
2794 const unsigned int nBfields = m_bFields.size();
2795 const unsigned int nAngles = m_bAngles.size();
2796 Init(nEfields, nBfields, nAngles, m_eVelE, 0.);
2797 Init(nEfields, nBfields, nAngles, m_eVelB, 0.);
2798 Init(nEfields, nBfields, nAngles, m_eVelX, 0.);
2799 Init(nEfields, nBfields, nAngles, m_eDifL, 0.);
2800 Init(nEfields, nBfields, nAngles, m_eDifT, 0.);
2801 Init(nEfields, nBfields, nAngles, m_eLor, 0.);
2802 Init(nEfields, nBfields, nAngles, m_eAlp, -30.);
2803 Init(nEfields, nBfields, nAngles, m_eAlp0, -30.);
2804 Init(nEfields, nBfields, nAngles, m_eAtt, -30.);
2805 Init(nEfields, nBfields, nAngles, 6, m_eDifM, 0.);
2806
2807 m_excRates.clear();
2808 m_ionRates.clear();
2809 // Retrieve the excitation and ionisation cross-sections in the gas mixture.
2810 GetExcitationIonisationLevels();
2811 std::cout << m_className << "::GenerateGasTable: Found "
2812 << m_excLevels.size() << " excitations and "
2813 << m_ionLevels.size() << " ionisations.\n";
2814 for (const auto& exc : m_excLevels) {
2815 std::cout << " " << exc.label << ", energy = " << exc.energy << " eV.\n";
2816 }
2817 for (const auto& ion : m_ionLevels) {
2818 std::cout << " " << ion.label << ", energy = " << ion.energy << " eV.\n";
2819 }
2820 if (!m_excLevels.empty()) {
2821 Init(nEfields, nBfields, nAngles, m_excLevels.size(), m_excRates, 0.);
2822 }
2823 if (!m_ionLevels.empty()) {
2824 Init(nEfields, nBfields, nAngles, m_ionLevels.size(), m_ionRates, 0.);
2825 }
2826 double vx = 0., vy = 0., vz = 0.;
2827 double difl = 0., dift = 0.;
2828 double alpha = 0., eta = 0.;
2829 double lor = 0.;
2830 double vxerr = 0., vyerr = 0., vzerr = 0.;
2831 double diflerr = 0., difterr = 0.;
2832 double alphaerr = 0., etaerr = 0.;
2833 double alphatof = 0.;
2834 double lorerr = 0.;
2835 std::array<double, 6> difftens;
2836
2837 // Run through the grid of E- and B-fields and angles.
2838 for (unsigned int i = 0; i < nEfields; ++i) {
2839 const double e = m_eFields[i];
2840 for (unsigned int j = 0; j < nAngles; ++j) {
2841 const double a = m_bAngles[j];
2842 for (unsigned int k = 0; k < nBfields; ++k) {
2843 const double b = m_bFields[k];
2844 std::cout << m_className << "::GenerateGasTable: E = " << e
2845 << " V/cm, B = " << b << " T, angle: " << a << " rad\n";
2846 RunMagboltz(e, b, a, numColl, verbose, vx, vy, vz, difl, dift, alpha,
2847 eta, lor, vxerr, vyerr, vzerr, diflerr, difterr, alphaerr,
2848 etaerr, lorerr, alphatof, difftens);
2849 m_eVelE[j][k][i] = vz;
2850 m_eVelX[j][k][i] = vy;
2851 m_eVelB[j][k][i] = vx;
2852 m_eDifL[j][k][i] = difl;
2853 m_eDifT[j][k][i] = dift;
2854 m_eLor[j][k][i] = lor;
2855 m_eAlp[j][k][i] = alpha > 0. ? log(alpha) : -30.;
2856 m_eAlp0[j][k][i] = m_eAlp[j][k][i];
2857 m_eAtt[j][k][i] = eta > 0. ? log(eta) : -30.;
2858 for (unsigned int l = 0; l < 6; ++l) {
2859 m_eDifM[l][j][k][i] = difftens[l];
2860 }
2861 // Retrieve the excitation and ionisation rates.
2862 unsigned int nNonZero = 0;
2863 if (m_useGasMotion) {
2864 // Retrieve the total collision frequency and number of collisions.
2865 double ftot = 0., fel = 0., fion = 0., fatt = 0., fin = 0.;
2866 std::int64_t ntotal = 0;
2867 Magboltz::colft_(&ftot, &fel, &fion, &fatt, &fin, &ntotal);
2868 if (ntotal == 0) continue;
2869 // Convert from ps-1 to ns-1.
2870 const double scale = 1.e3 * ftot / ntotal;
2871 for (unsigned int ig = 0; ig < m_nComponents; ++ig) {
2872 const auto nL = Magboltz::larget_.last[ig];
2873 for (std::int64_t il = 0; il < nL; ++il) {
2874 if (Magboltz::larget_.iarry[il][ig] <= 0) break;
2875 // Skip levels that are not ionisations or inelastic collisions.
2876 const int cstype = (Magboltz::larget_.iarry[il][ig] - 1) % 5;
2877 if (cstype != 1 && cstype != 3) continue;
2878 // const int igas = int((Magboltz::larget_.iarry[il][ig] - 1) / 5);
2879 auto descr = GetDescription(il, ig, Magboltz::script_.dscrpt);
2880 descr = m_gas[ig] + descr;
2881 if (cstype == 3) {
2882 const unsigned int nExc = m_excLevels.size();
2883 for (unsigned int ie = 0; ie < nExc; ++ie) {
2884 if (descr != m_excLevels[ie].label) continue;
2885 const auto ncoll = Magboltz::outptt_.icoln[il][ig];
2886 m_excRates[ie][j][k][i] = scale * ncoll;
2887 if (ncoll > 0) ++nNonZero;
2888 break;
2889 }
2890 } else if (cstype == 1) {
2891 const unsigned int nIon = m_ionLevels.size();
2892 for (unsigned int ii = 0; ii < nIon; ++ii) {
2893 if (descr != m_ionLevels[ii].label) continue;
2894 const auto ncoll = Magboltz::outptt_.icoln[il][ig];
2895 m_ionRates[ii][j][k][i] = scale * ncoll;
2896 if (ncoll > 0) ++nNonZero;
2897 break;
2898 }
2899 }
2900 }
2901 }
2902 } else {
2903 // Retrieve the total collision frequency and number of collisions.
2904 double ftot = 0., fel = 0., fion = 0., fatt = 0., fin = 0.;
2905 std::int64_t ntotal = 0;
2906 Magboltz::colf_(&ftot, &fel, &fion, &fatt, &fin, &ntotal);
2907 if (ntotal == 0) continue;
2908 // Convert from ps-1 to ns-1.
2909 const double scale = 1.e3 * ftot / ntotal;
2910 for (std::int64_t il = 0; il < Magboltz::nMaxLevels; ++il) {
2911 if (Magboltz::large_.iarry[il] <= 0) break;
2912 // Skip levels that are not ionisations or inelastic collisions.
2913 const int cstype = (Magboltz::large_.iarry[il] - 1) % 5;
2914 if (cstype != 1 && cstype != 3) continue;
2915 const int igas = int((Magboltz::large_.iarry[il] - 1) / 5);
2916 std::string descr = GetDescription(il, Magboltz::scrip_.dscrpt);
2917 descr = m_gas[igas] + descr;
2918 if (cstype == 3) {
2919 const unsigned int nExc = m_excLevels.size();
2920 for (unsigned int ie = 0; ie < nExc; ++ie) {
2921 if (descr != m_excLevels[ie].label) continue;
2922 m_excRates[ie][j][k][i] = scale * Magboltz::outpt_.icoln[il];
2923 if (Magboltz::outpt_.icoln[il] > 0) ++nNonZero;
2924 break;
2925 }
2926 } else if (cstype == 1) {
2927 const unsigned int nIon = m_ionLevels.size();
2928 for (unsigned int ii = 0; ii < nIon; ++ii) {
2929 if (descr != m_ionLevels[ii].label) continue;
2930 m_ionRates[ii][j][k][i] = scale * Magboltz::outpt_.icoln[il];
2931 if (Magboltz::outpt_.icoln[il] > 0) ++nNonZero;
2932 break;
2933 }
2934 }
2935 }
2936 }
2937 if (nNonZero > 0) {
2938 std::cout << " Excitation and ionisation rates:\n";
2939 std::cout << " Level "
2940 << " Rate [ns-1]\n";
2941 const unsigned int nExc = m_excLevels.size();
2942 for (unsigned int ie = 0; ie < nExc; ++ie) {
2943 if (m_excRates[ie][j][k][i] <= 0) continue;
2944 std::cout << std::setw(60) << m_excLevels[ie].label;
2945 std::printf(" %15.8f\n", m_excRates[ie][j][k][i]);
2946 }
2947 const unsigned int nIon = m_ionLevels.size();
2948 for (unsigned int ii = 0; ii < nIon; ++ii) {
2949 if (m_ionRates[ii][j][k][i] <= 0) continue;
2950 std::cout << std::setw(60) << m_ionLevels[ii].label;
2951 std::printf(" %15.8f\n", m_ionRates[ii][j][k][i]);
2952 }
2953 }
2954 }
2955 }
2956 }
2957 // Set the threshold indices.
2960}
std::vector< std::vector< std::vector< std::vector< double > > > > m_excRates
Definition MediumGas.hh:190
std::vector< IonLevel > m_ionLevels
Definition MediumGas.hh:207
std::vector< std::vector< std::vector< std::vector< double > > > > m_ionRates
Definition MediumGas.hh:191
std::vector< std::vector< std::vector< double > > > m_eAlp0
Definition MediumGas.hh:187
void RunMagboltz(const double e, const double b, const double btheta, const int ncoll, bool verbose, double &vx, double &vy, double &vz, double &dl, double &dt, double &alpha, double &eta, double &lor, double &vxerr, double &vyerr, double &vzerr, double &dlerr, double &dterr, double &alphaerr, double &etaerr, double &lorerr, double &alphatof, std::array< double, 6 > &difftens)
std::vector< double > m_bFields
Definition Medium.hh:573
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition Medium.hh:582
size_t SetThreshold(const std::vector< std::vector< std::vector< double > > > &tab) const
Definition Medium.cc:1252
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition Medium.hh:577
std::vector< std::vector< std::vector< double > > > m_eVelX
Definition Medium.hh:578
std::vector< std::vector< std::vector< double > > > m_eDifL
Definition Medium.hh:580
void Init(const size_t nE, const size_t nB, const size_t nA, std::vector< std::vector< std::vector< double > > > &tab, const double val)
Definition Medium.cc:1408
std::vector< double > m_eFields
Definition Medium.hh:572
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition Medium.hh:583
std::vector< double > m_bAngles
Definition Medium.hh:574
std::vector< std::vector< std::vector< double > > > m_eLor
Definition Medium.hh:584
std::vector< std::vector< std::vector< double > > > m_eDifT
Definition Medium.hh:581
std::vector< std::vector< std::vector< double > > > m_eVelB
Definition Medium.hh:579
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
Definition Medium.hh:586
struct Garfield::Magboltz::@071266015331337350342074140245340115315173206173 script_
struct Garfield::Magboltz::@236122006074061251067327101212023331117065272002 larget_
struct Garfield::Magboltz::@016351066267021123171150005363312076230027317210 outpt_
struct Garfield::Magboltz::@304340026334031267143114007126147372164074306062 outptt_
void colf_(double *freq, double *freel, double *freion, double *freatt, double *frein, std::int64_t *ntotal)
struct Garfield::Magboltz::@363161251140110014365135052375245247262302340235 large_
void colft_(double *freq, double *freel, double *freion, double *freatt, double *frein, std::int64_t *ntotal)
struct Garfield::Magboltz::@024075325050032167066064232170013161226346016145 scrip_

◆ GetDeexcitationProduct()

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

Reimplemented from Garfield::Medium.

Definition at line 835 of file MediumMagboltz.cc.

837 {
838 if (i >= m_dxcProducts.size() || !(m_useDeexcitation || m_usePenning)) {
839 return false;
840 }
841 t = m_dxcProducts[i].t;
842 s = m_dxcProducts[i].s;
843 type = m_dxcProducts[i].type;
844 energy = m_dxcProducts[i].energy;
845 return true;
846}

◆ GetElectronCollisionRate() [1/2]

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

Get the (real) collision rate [ns-1] at a given electron energy e [eV].

Reimplemented from Garfield::Medium.

Definition at line 517 of file MediumMagboltz.cc.

518 {
519 // Check if the electron energy is within the currently set range.
520 if (e <= 0.) {
521 std::cerr << m_className << "::GetElectronCollisionRate: Invalid energy.\n";
522 return m_cfTot[0];
523 }
524 if (e > m_eMax) {
525 std::cerr << m_className << "::GetElectronCollisionRate:\n Rate at " << e
526 << " eV is not included in the current table.\n "
527 << "Increasing energy range to " << 1.05 * e << " eV.\n";
528 SetMaxElectronEnergy(1.05 * e);
529 }
530
531 // If necessary, update the collision rates table.
532 if (!Update()) return 0.;
533
534 // Get the energy interval.
535 if (e < m_eHigh) {
536 // Linear binning
537 return m_cfTot[int(e * m_eStepInv)];
538 }
539
540 // Logarithmic binning
541 const double eLog = log(e);
542 int iE = int((eLog - m_eHighLog) / m_lnStep);
543 // Calculate the collision rate by log-log interpolation.
544 const double fmax = m_cfTotLog[iE];
545 const double fmin = iE == 0 ? log(m_cfTot.back()) : m_cfTotLog[iE - 1];
546 const double emin = m_eHighLog + iE * m_lnStep;
547 const double f = fmin + (eLog - emin) * (fmax - fmin) / m_lnStep;
548 return exp(f);
549}
DoubleAc exp(const DoubleAc &f)
Definition DoubleAc.cpp:377

Referenced by GetElectronCollisionRate().

◆ GetElectronCollisionRate() [2/2]

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

Get the collision rate [ns-1] for a specific level.

Definition at line 551 of file MediumMagboltz.cc.

553 {
554 // Check if the electron energy is within the currently set range.
555 if (e <= 0.) {
556 std::cerr << m_className << "::GetElectronCollisionRate: Invalid energy.\n";
557 return 0.;
558 }
559
560 // Check if the level exists.
561 if (level >= m_nTerms) {
562 std::cerr << m_className << "::GetElectronCollisionRate: Invalid level.\n";
563 return 0.;
564 }
565
566 // Get the total scattering rate.
567 double rate = GetElectronCollisionRate(e, band);
568 // Get the energy interval.
569 if (e < m_eHigh) {
570 // Linear binning
571 const int iE = int(e * m_eStepInv);
572 if (level == 0) {
573 rate *= m_cf[iE][0];
574 } else {
575 rate *= m_cf[iE][level] - m_cf[iE][level - 1];
576 }
577 } else {
578 // Logarithmic binning
579 const int iE = int((log(e) - m_eHighLog) / m_lnStep);
580 if (level == 0) {
581 rate *= m_cfLog[iE][0];
582 } else {
583 rate *= m_cfLog[iE][level] - m_cfLog[iE][level - 1];
584 }
585 }
586 return rate;
587}
double GetElectronCollisionRate(const double e, const int band) override
Get the (real) collision rate [ns-1] at a given electron energy e [eV].

◆ GetElectronNullCollisionRate()

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

Get the overall null-collision rate [ns-1].

Reimplemented from Garfield::Medium.

Definition at line 511 of file MediumMagboltz.cc.

511 {
512 // If necessary, update the collision rates table.
513 if (!Update()) return 0.;
514 return m_cfNull;
515}

◆ GetGasNumberMagboltz()

int Garfield::MediumMagboltz::GetGasNumberMagboltz ( const std::string & input)
static

Definition at line 1099 of file MediumMagboltz.cc.

1099 {
1100
1101 if (input.empty()) return 0;
1102
1103 if (input == "CF4") {
1104 return 1;
1105 } else if (input == "Ar") {
1106 return 2;
1107 } else if (input == "He" || input == "He-4") {
1108 // Helium 4
1109 return 3;
1110 } else if (input == "He-3") {
1111 // Helium 3
1112 return 4;
1113 } else if (input == "Ne") {
1114 return 5;
1115 } else if (input == "Kr") {
1116 return 6;
1117 } else if (input == "Xe") {
1118 return 7;
1119 } else if (input == "CH4") {
1120 // Methane
1121 return 8;
1122 } else if (input == "C2H6") {
1123 // Ethane
1124 return 9;
1125 } else if (input == "C3H8") {
1126 // Propane
1127 return 10;
1128 } else if (input == "iC4H10") {
1129 // Isobutane
1130 return 11;
1131 } else if (input == "CO2") {
1132 return 12;
1133 } else if (input == "neoC5H12") {
1134 // Neopentane
1135 return 13;
1136 } else if (input == "H2O") {
1137 return 14;
1138 } else if (input == "O2") {
1139 return 15;
1140 } else if (input == "N2") {
1141 return 16;
1142 } else if (input == "NO") {
1143 // Nitric oxide (NO)
1144 return 17;
1145 } else if (input == "N2O") {
1146 // Nitrous oxide (N2O)
1147 return 18;
1148 } else if (input == "C2H4") {
1149 // Ethene (C2H4)
1150 return 19;
1151 } else if (input == "C2H2") {
1152 // Acetylene (C2H2)
1153 return 20;
1154 } else if (input == "H2") {
1155 // Hydrogen
1156 return 21;
1157 } else if (input == "D2") {
1158 // Deuterium
1159 return 22;
1160 } else if (input == "CO") {
1161 // Carbon monoxide (CO)
1162 return 23;
1163 } else if (input == "Methylal") {
1164 // Methylal (dimethoxymethane, CH3-O-CH2-O-CH3, "hot" version)
1165 return 24;
1166 } else if (input == "DME") {
1167 return 25;
1168 } else if (input == "Reid-Step") {
1169 return 26;
1170 } else if (input == "Maxwell-Model") {
1171 return 27;
1172 } else if (input == "Reid-Ramp") {
1173 return 28;
1174 } else if (input == "C2F6") {
1175 return 29;
1176 } else if (input == "SF6") {
1177 return 30;
1178 } else if (input == "NH3") {
1179 return 31;
1180 } else if (input == "C3H6") {
1181 // Propene
1182 return 32;
1183 } else if (input == "cC3H6") {
1184 // Cyclopropane
1185 return 33;
1186 } else if (input == "CH3OH") {
1187 // Methanol
1188 return 34;
1189 } else if (input == "C2H5OH") {
1190 // Ethanol
1191 return 35;
1192 } else if (input == "C3H7OH") {
1193 // Propanol
1194 return 36;
1195 } else if (input == "Cs") {
1196 return 37;
1197 } else if (input == "F2") {
1198 // Fluorine
1199 return 38;
1200 } else if (input == "CS2") {
1201 return 39;
1202 } else if (input == "COS") {
1203 return 40;
1204 } else if (input == "CD4") {
1205 // Deuterated methane
1206 return 41;
1207 } else if (input == "BF3") {
1208 return 42;
1209 } else if (input == "C2HF5" || input == "C2H2F4") {
1210 return 43;
1211 } else if (input == "TMA") {
1212 return 44;
1213 } else if (input == "nC3H7OH") {
1214 // n-propanol
1215 return 46;
1216 } else if (input == "paraH2") {
1217 // Para hydrogen
1218 return 47;
1219 } else if (input == "orthoD2") {
1220 // Ortho deuterium
1221 return 48;
1222 } else if (input == "CHF3") {
1223 return 50;
1224 } else if (input == "CF3Br") {
1225 return 51;
1226 } else if (input == "C3F8") {
1227 return 52;
1228 } else if (input == "O3") {
1229 // Ozone
1230 return 53;
1231 } else if (input == "Hg") {
1232 // Mercury
1233 return 54;
1234 } else if (input == "H2S") {
1235 return 55;
1236 } else if (input == "nC4H10") {
1237 // n-Butane
1238 return 56;
1239 } else if (input == "nC5H12") {
1240 // n-Pentane
1241 return 57;
1242 } else if (input == "N2 (Phelps)") {
1243 return 58;
1244 } else if (input == "GeH4") {
1245 // Germane, GeH4
1246 return 59;
1247 } else if (input == "SiH4") {
1248 // Silane, SiH4
1249 return 60;
1250 } else if (input == "CCl4") {
1251 return 61;
1252 }
1253
1254 std::cerr << "MediumMagboltz::GetGasNumberMagboltz:\n"
1255 << " Gas " << input << " is not defined.\n";
1256 return 0;
1257}

Referenced by Garfield::TrackDegrade::Initialise(), and RunMagboltz().

◆ GetLevel()

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.

Definition at line 996 of file MediumMagboltz.cc.

997 {
998 if (!Update()) return false;
999
1000 if (i >= m_nTerms) {
1001 std::cerr << m_className << "::GetLevel: Index out of range.\n";
1002 return false;
1003 }
1004
1005 // Collision type
1006 type = m_csType[i] % nCsTypes;
1007 ngas = int(m_csType[i] / nCsTypes);
1008 // Description (from Magboltz)
1009 descr = m_description[i];
1010 // Threshold energy
1011 e = m_rgas[ngas] * m_energyLoss[i];
1012 if (m_debug) {
1013 std::cout << m_className << "::GetLevel:\n"
1014 << " Level " << i << ": " << descr << "\n"
1015 << " Type " << type << "\n"
1016 << " Threshold energy: " << e << " eV\n";
1017 if (type == ElectronCollisionTypeExcitation && m_usePenning &&
1018 e > m_minIonPot) {
1019 std::cout << " Penning transfer coefficient: " << m_rPenning[i]
1020 << "\n";
1021 } else if (type == ElectronCollisionTypeExcitation && m_useDeexcitation) {
1022 const int idxc = m_iDeexcitation[i];
1023 if (idxc < 0 || idxc >= (int)m_deexcitations.size()) {
1024 std::cout << " Deexcitation cascade not implemented.\n";
1025 return true;
1026 }
1027 const auto& dxc = m_deexcitations[idxc];
1028 if (dxc.osc > 0.) {
1029 std::cout << " Oscillator strength: " << dxc.osc << "\n";
1030 }
1031 std::cout << " Decay channels:\n";
1032 const int nChannels = dxc.type.size();
1033 for (int j = 0; j < nChannels; ++j) {
1034 if (dxc.type[j] == DxcTypeRad) {
1035 std::cout << " Radiative decay to ";
1036 if (dxc.final[j] < 0) {
1037 std::cout << "ground state: ";
1038 } else {
1039 std::cout << m_deexcitations[dxc.final[j]].label << ": ";
1040 }
1041 } else if (dxc.type[j] == DxcTypeCollIon) {
1042 if (dxc.final[j] < 0) {
1043 std::cout << " Penning ionisation: ";
1044 } else {
1045 std::cout << " Associative ionisation: ";
1046 }
1047 } else if (dxc.type[j] == DxcTypeCollNonIon) {
1048 if (dxc.final[j] >= 0) {
1049 std::cout << " Collision-induced transition to "
1050 << m_deexcitations[dxc.final[j]].label << ": ";
1051 } else {
1052 std::cout << " Loss: ";
1053 }
1054 }
1055 const double br = j == 0 ? dxc.p[j] : dxc.p[j] - dxc.p[j - 1];
1056 std::cout << std::setprecision(5) << br * 100. << "%\n";
1057 }
1058 }
1059 }
1060
1061 return true;
1062}

◆ GetMaxElectronEnergy()

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

Get the highest electron energy in the table of scattering rates.

Definition at line 33 of file MediumMagboltz.hh.

33{ return m_eMax; }

◆ GetMaxPhotonEnergy()

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

Get the highest photon energy in the table of scattering rates.

Definition at line 39 of file MediumMagboltz.hh.

39{ return m_eFinalGamma; }

◆ GetNumberOfDeexcitationProducts()

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

Reimplemented from Garfield::Medium.

Definition at line 96 of file MediumMagboltz.hh.

96 {
97 return m_dxcProducts.size();
98 }

◆ GetNumberOfElectronCollisions() [1/3]

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

Get the total number of electron collisions.

Definition at line 973 of file MediumMagboltz.cc.

973 {
974 return std::accumulate(std::begin(m_nCollisions), std::end(m_nCollisions), 0);
975}

◆ GetNumberOfElectronCollisions() [2/3]

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

Get the number of collisions for a specific cross-section term.

Definition at line 1075 of file MediumMagboltz.cc.

1076 {
1077 if (level >= m_nTerms) {
1078 std::cerr << m_className << "::GetNumberOfElectronCollisions: "
1079 << "Level " << level << " does not exist.\n";
1080 return 0;
1081 }
1082 return m_nCollisionsDetailed[level];
1083}

◆ GetNumberOfElectronCollisions() [3/3]

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.

Definition at line 977 of file MediumMagboltz.cc.

980 {
981 nElastic = m_nCollisions[ElectronCollisionTypeElastic];
982 nIonisation = m_nCollisions[ElectronCollisionTypeIonisation];
983 nAttachment = m_nCollisions[ElectronCollisionTypeAttachment];
984 nInelastic = m_nCollisions[ElectronCollisionTypeInelastic];
985 nExcitation = m_nCollisions[ElectronCollisionTypeExcitation];
986 nSuperelastic = m_nCollisions[ElectronCollisionTypeSuperelastic];
987 return nElastic + nIonisation + nAttachment + nInelastic + nExcitation +
988 nSuperelastic;
989}

◆ GetNumberOfLevels()

unsigned int Garfield::MediumMagboltz::GetNumberOfLevels ( )

Get the number of cross-section terms.

Definition at line 991 of file MediumMagboltz.cc.

991 {
992 if (!Update()) return 0;
993 return m_nTerms;
994}

◆ GetNumberOfPenningTransfers()

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

Get the number of Penning transfers that occured since the last reset.

Definition at line 129 of file MediumMagboltz.hh.

129{ return m_nPenning; }

◆ GetNumberOfPhotonCollisions() [1/2]

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

Get the total number of photon collisions.

Definition at line 1085 of file MediumMagboltz.cc.

1085 {
1086 return std::accumulate(std::begin(m_nPhotonCollisions),
1087 std::end(m_nPhotonCollisions), 0);
1088}

◆ GetNumberOfPhotonCollisions() [2/2]

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

Get number of photon collisions by collision type.

Definition at line 1090 of file MediumMagboltz.cc.

1092 {
1093 nElastic = m_nPhotonCollisions[0];
1094 nIonising = m_nPhotonCollisions[1];
1095 nInelastic = m_nPhotonCollisions[2];
1096 return nElastic + nIonising + nInelastic;
1097}

◆ GetPenningTransfer()

bool Garfield::MediumMagboltz::GetPenningTransfer ( const unsigned int i,
double & r,
double & lambda )

Get the Penning transfer probability and distance of a specific level.

Definition at line 1064 of file MediumMagboltz.cc.

1065 {
1066 r = 0.;
1067 lambda = 0.;
1068 if (!Update()) return false;
1069 if (i >= m_nTerms) return false;
1070 r = m_rPenning[i];
1071 lambda = m_lambdaPenning[i];
1072 return true;
1073}

◆ GetPhotonCollision()

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

Reimplemented from Garfield::Medium.

Definition at line 879 of file MediumMagboltz.cc.

881 {
882 if (e <= 0.) {
883 std::cerr << m_className << "::GetPhotonCollision: Invalid energy.\n";
884 return false;
885 }
886 if (e > m_eFinalGamma) {
887 std::cerr << m_className << "::GetPhotonCollision:\n Provided energy ("
888 << e << " eV) exceeds current energy range.\n"
889 << " Increasing energy range to " << 1.05 * e << " eV.\n";
890 SetMaxPhotonEnergy(1.05 * e);
891 }
892
893 if (!Update()) return false;
894
895 // Energy interval
896 const int iE =
897 std::min(std::max(int(e / m_eStepGamma), 0), nEnergyStepsGamma - 1);
898
899 double r = m_cfTotGamma[iE];
900 if (m_useDeexcitation && m_useRadTrap && !m_deexcitations.empty()) {
901 int nLines = 0;
902 std::vector<double> pLine(0);
903 std::vector<int> iLine(0);
904 // Loop over the excitations.
905 const unsigned int nDeexcitations = m_deexcitations.size();
906 for (unsigned int i = 0; i < nDeexcitations; ++i) {
907 const auto& dxc = m_deexcitations[i];
908 if (dxc.cf > 0. && fabs(e - dxc.energy) <= dxc.width) {
909 r += dxc.cf *
910 TMath::Voigt(e - dxc.energy, dxc.sDoppler, 2 * dxc.gPressure);
911 pLine.push_back(r);
912 iLine.push_back(i);
913 ++nLines;
914 }
915 }
916 r *= RndmUniform();
917 if (nLines > 0 && r >= m_cfTotGamma[iE]) {
918 // Photon is absorbed by a discrete line.
919 for (int i = 0; i < nLines; ++i) {
920 if (r <= pLine[i]) {
921 ++m_nPhotonCollisions[PhotonCollisionTypeExcitation];
922 int fLevel = 0;
923 ComputeDeexcitationInternal(iLine[i], fLevel);
924 type = PhotonCollisionTypeExcitation;
925 nsec = m_dxcProducts.size();
926 return true;
927 }
928 }
929 std::cerr << m_className << "::GetPhotonCollision:\n";
930 std::cerr << " Random sampling of deexcitation line failed.\n";
931 std::cerr << " Program bug!\n";
932 return false;
933 }
934 } else {
935 r *= RndmUniform();
936 }
937
938 if (r <= m_cfGamma[iE][0]) {
939 level = 0;
940 } else if (r >= m_cfGamma[iE][m_nPhotonTerms - 1]) {
941 level = m_nPhotonTerms - 1;
942 } else {
943 const auto begin = m_cfGamma[iE].cbegin();
944 level = std::lower_bound(begin, begin + m_nPhotonTerms, r) - begin;
945 }
946
947 nsec = 0;
948 esec = e1 = 0.;
949 type = csTypeGamma[level];
950 // Collision type
951 type = type % nCsTypesGamma;
952 int ngas = int(csTypeGamma[level] / nCsTypesGamma);
953 ++m_nPhotonCollisions[type];
954 // Ionising collision
955 if (type == 1) {
956 esec = std::max(e - m_ionPot[ngas], Small);
957 nsec = 1;
958 }
959
960 // Determine the scattering angle
961 ctheta = 2 * RndmUniform() - 1.;
962
963 return true;
964}
bool SetMaxPhotonEnergy(const double e)
DoubleAc fabs(const DoubleAc &f)
Definition DoubleAc.h:615

◆ GetPhotonCollisionRate()

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

Reimplemented from Garfield::Medium.

Definition at line 848 of file MediumMagboltz.cc.

848 {
849 if (e <= 0.) {
850 std::cerr << m_className << "::GetPhotonCollisionRate: Invalid energy.\n";
851 return m_cfTotGamma[0];
852 }
853 if (e > m_eFinalGamma) {
854 std::cerr << m_className << "::GetPhotonCollisionRate:\n Rate at " << e
855 << " eV is not included in the current table.\n"
856 << " Increasing energy range to " << 1.05 * e << " eV.\n";
857 SetMaxPhotonEnergy(1.05 * e);
858 }
859
860 if (!Update()) return 0.;
861
862 const int iE =
863 std::min(std::max(int(e / m_eStepGamma), 0), nEnergyStepsGamma - 1);
864
865 double cfSum = m_cfTotGamma[iE];
866 if (m_useDeexcitation && m_useRadTrap && !m_deexcitations.empty()) {
867 // Loop over the excitations.
868 for (const auto& dxc : m_deexcitations) {
869 if (dxc.cf > 0. && fabs(e - dxc.energy) <= dxc.width) {
870 cfSum += dxc.cf *
871 TMath::Voigt(e - dxc.energy, dxc.sDoppler, 2 * dxc.gPressure);
872 }
873 }
874 }
875
876 return cfSum;
877}

◆ Initialise()

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

Definition at line 418 of file MediumMagboltz.cc.

418 {
419 if (!m_isChanged) {
420 if (m_debug) {
421 std::cerr << m_className << "::Initialise: Nothing changed.\n";
422 }
423 return true;
424 }
425 return Update(verbose);
426}

Referenced by PrintGas().

◆ PlotElectronCrossSections()

void Garfield::MediumMagboltz::PlotElectronCrossSections ( )

Definition at line 1871 of file MediumMagboltz.cc.

1871 {
1872
1873 if (!Update()) return;
1874
1875 std::array<float, Magboltz::nEnergySteps> en;
1876 for (unsigned int k = 0; k < Magboltz::nEnergySteps; ++k) {
1877 en[k] = (k + 0.5) * m_eStep;
1878 }
1879 std::array<std::array<float, Magboltz::nEnergySteps>, 5> cs;
1880 for (unsigned int i = 0; i < m_nComponents; ++i) {
1881 for (size_t j = 0; j < 5; ++j) cs[j].fill(0.);
1882 const double density = GetNumberDensity() * m_fraction[i];
1883 const double scale = sqrt(0.5 * ElectronMass) / (density * SpeedOfLight);
1884 for (unsigned int j = 0; j < m_nTerms; ++j) {
1885 if (int(m_csType[j] / nCsTypes) != int(i)) continue;
1886 int cstype = m_csType[j] % nCsTypes;
1887 if (cstype >= ElectronCollisionTypeVirtual) continue;
1888 // Group inelastic and superelastic collisions.
1889 if (cstype == 5) cstype = 3;
1890 for (unsigned int k = 0; k < Magboltz::nEnergySteps; ++k) {
1891 double cf = m_cf[k][j];
1892 if (j > 0) cf -= m_cf[k][j - 1];
1893 cs[cstype][k] += cf * 1.e18 * scale;
1894 }
1895 }
1896 const std::string name = ViewBase::FindUnusedCanvasName("cCs");
1897 TCanvas* canvas = new TCanvas(name.c_str(), m_gas[i].c_str(), 800, 600);
1898 canvas->cd();
1899 canvas->SetLogx();
1900 canvas->SetLogy();
1901 canvas->SetGridx();
1902 canvas->SetGridy();
1903 canvas->DrawFrame(en[0], 0.01, en.back(), 100.,
1904 ";energy [eV];#sigma [Mbarn]");
1905 TGraph gr(Magboltz::nEnergySteps);
1906 gr.SetLineWidth(3);
1907 const std::array<short, 5> cols = {kBlack, kCyan - 2, kRed + 2,
1908 kGreen + 3, kMagenta + 3};
1909 for (size_t j = 0; j < 5; ++j) {
1910 if (*std::max_element(cs[j].begin(), cs[j].end()) < 1.e-10) continue;
1911 gr.SetLineColor(cols[j]);
1912 gr.DrawGraph(Magboltz::nEnergySteps, en.data(), cs[j].data(), "lsame");
1913 }
1914 canvas->Update();
1915 }
1916
1917}
double GetNumberDensity() const override
Get the number density [cm-3].
Definition MediumGas.cc:304
std::array< double, m_nMaxGases > m_fraction
Definition MediumGas.hh:165
static std::string FindUnusedCanvasName(const std::string &s)
Find an unused canvas name.
Definition ViewBase.cc:337
double cf[nMaxLevels][nEnergySteps]

◆ PrintGas()

void Garfield::MediumMagboltz::PrintGas ( )
overridevirtual

Print information about the present gas mixture and available data.

Reimplemented from Garfield::MediumGas.

Definition at line 428 of file MediumMagboltz.cc.

428 {
430
431 if (m_isChanged) {
432 if (!Initialise()) return;
433 }
434
435 std::cout << " Electron cross-sections:\n";
436 int igas = -1;
437 for (unsigned int i = 0; i < m_nTerms; ++i) {
438 // Collision type
439 int type = m_csType[i] % nCsTypes;
440 if (igas != int(m_csType[i] / nCsTypes)) {
441 igas = int(m_csType[i] / nCsTypes);
442 std::cout << " " << m_gas[igas] << "\n";
443 }
444 // Description (from Magboltz)
445 // Threshold energy
446 double e = m_rgas[igas] * m_energyLoss[i];
447 std::cout << " Level " << i << ": " << m_description[i] << "\n";
448 std::cout << " Type " << type;
449 if (type == ElectronCollisionTypeElastic) {
450 std::cout << " (elastic)\n";
451 } else if (type == ElectronCollisionTypeIonisation) {
452 std::cout << " (ionisation). Ionisation threshold: " << e << " eV.\n";
453 } else if (type == ElectronCollisionTypeAttachment) {
454 std::cout << " (attachment)\n";
455 } else if (type == ElectronCollisionTypeInelastic) {
456 std::cout << " (inelastic). Energy loss: " << e << " eV.\n";
457 } else if (type == ElectronCollisionTypeExcitation) {
458 std::cout << " (excitation). Excitation energy: " << e << " eV.\n";
459 } else if (type == ElectronCollisionTypeSuperelastic) {
460 std::cout << " (super-elastic). Energy gain: " << -e << " eV.\n";
461 } else if (type == ElectronCollisionTypeVirtual) {
462 std::cout << " (virtual)\n";
463 } else {
464 std::cout << " (unknown)\n";
465 }
466 if (type == ElectronCollisionTypeExcitation && m_usePenning &&
467 e > m_minIonPot) {
468 std::cout << " Penning transfer coefficient: "
469 << m_rPenning[i] << "\n";
470 } else if (type == ElectronCollisionTypeExcitation && m_useDeexcitation) {
471 const int idxc = m_iDeexcitation[i];
472 if (idxc < 0 || idxc >= (int)m_deexcitations.size()) {
473 std::cout << " Deexcitation cascade not implemented.\n";
474 continue;
475 }
476 const auto& dxc = m_deexcitations[idxc];
477 if (dxc.osc > 0.) {
478 std::cout << " Oscillator strength: " << dxc.osc << "\n";
479 }
480 std::cout << " Decay channels:\n";
481 const int nChannels = dxc.type.size();
482 for (int j = 0; j < nChannels; ++j) {
483 if (dxc.type[j] == DxcTypeRad) {
484 std::cout << " Radiative decay to ";
485 if (dxc.final[j] < 0) {
486 std::cout << "ground state: ";
487 } else {
488 std::cout << m_deexcitations[dxc.final[j]].label << ": ";
489 }
490 } else if (dxc.type[j] == DxcTypeCollIon) {
491 if (dxc.final[j] < 0) {
492 std::cout << " Penning ionisation: ";
493 } else {
494 std::cout << " Associative ionisation: ";
495 }
496 } else if (dxc.type[j] == DxcTypeCollNonIon) {
497 if (dxc.final[j] >= 0) {
498 std::cout << " Collision-induced transition to "
499 << m_deexcitations[dxc.final[j]].label << ": ";
500 } else {
501 std::cout << " Loss: ";
502 }
503 }
504 const double br = j == 0 ? dxc.p[j] : dxc.p[j] - dxc.p[j - 1];
505 std::cout << std::setprecision(5) << br * 100. << "%\n";
506 }
507 }
508 }
509}
virtual void PrintGas()
Print information about the present gas mixture and available data.
bool Initialise(const bool verbose=false)

◆ ResetCollisionCounters()

void Garfield::MediumMagboltz::ResetCollisionCounters ( )

Reset the collision counters.

Definition at line 966 of file MediumMagboltz.cc.

966 {
967 m_nCollisions.fill(0);
968 m_nCollisionsDetailed.assign(m_nTerms, 0);
969 m_nPenning = 0;
970 m_nPhotonCollisions.fill(0);
971}

◆ RunMagboltz()

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

Run Magboltz for a given electric field, magnetic field and angle.

Parameters
[in]eelectric field
[in]bmagnetic field
[in]bthetaangle between electric and magnetic field
[in]ncollnumber of collisions (in multiples of 107) to be simulated
[in]verboseverbosity flag
[out]vx,vy,vzdrift velocity vector
[out]dl,dtdiffusion cofficients
[out]alphaTownsend cofficient
[out]etaattachment cofficient
[out]lorLorentz angle
[out]vxerr,vyerr,vzerrerrors on drift velocity
[out]dlerr,dterrerrors on diffusion coefficients
[out]alphaerr,etaerrerrors on Townsend/attachment coefficients
[out]lorerrerror on Lorentz angle
[out]alphatofeffective Townsend coefficient $(\alpha - \eta)$ calculated using time-of-flight method
[out]difftenscomponents of the diffusion tensor (zz, xx, yy, xz, yz, xy)

Definition at line 2661 of file MediumMagboltz.cc.

2667 {
2668 // Initialize the values.
2669 vx = vy = vz = 0.;
2670 dl = dt = 0.;
2671 alpha = eta = alphatof = 0.;
2672 lor = 0.;
2673 vxerr = vyerr = vzerr = 0.;
2674 dlerr = dterr = 0.;
2675 alphaerr = etaerr = 0.;
2676 lorerr = 0.;
2677
2678 // Set the input parameters in the Magboltz common blocks.
2680 Magboltz::inpt_.nAniso = 2;
2681 if (m_autoEnergyLimit) {
2682 Magboltz::inpt_.efinal = 0.;
2683 } else {
2684 Magboltz::inpt_.efinal = std::min(m_eMax, m_eHigh);
2685 }
2686 Magboltz::inpt_.tempc = m_temperature - ZeroCelsius;
2688 Magboltz::inpt_.ipen = 0;
2689 Magboltz::setp_.nmax = ncoll;
2690
2691 Magboltz::thrm_.ithrm = m_useGasMotion ? 1 : 0;
2692
2693 Magboltz::setp_.efield = emag;
2694 // Convert from Tesla to kGauss.
2695 Magboltz::bfld_.bmag = bmag * 10.;
2696 // Convert from radians to degree.
2697 Magboltz::bfld_.btheta = btheta * RadToDegree;
2698
2699 // Set the gas composition in Magboltz.
2700 for (unsigned int i = 0; i < m_nComponents; ++i) {
2701 const int ng = GetGasNumberMagboltz(m_gas[i]);
2702 if (ng <= 0) {
2703 std::cerr << m_className << "::RunMagboltz:\n Gas " << m_gas[i]
2704 << " does not have a gas number in Magboltz.\n";
2705 return;
2706 }
2707 Magboltz::gasn_.ngasn[i] = ng;
2708 Magboltz::ratio_.frac[i] = 100 * m_fraction[i];
2709 }
2710
2711 // Run Magboltz.
2713
2714 // Velocities. Convert to cm / ns.
2715 vx = Magboltz::vel_.wx * 1.e-9;
2716 vxerr = Magboltz::velerr_.dwx;
2717 vy = Magboltz::vel_.wy * 1.e-9;
2718 vyerr = Magboltz::velerr_.dwy;
2719 vz = Magboltz::vel_.wz * 1.e-9;
2720 vzerr = Magboltz::velerr_.dwz;
2721
2722 // Calculate the Lorentz angle.
2723 const double vt = sqrt(vx * vx + vy * vy);
2724 const double v2 = (vx * vx + vy * vy + vz * vz);
2725 lor = atan2(vt, vz);
2726 if (vt > 0. && v2 > 0. && fabs(lor) > 0.) {
2727 const double dvx = vx * vxerr;
2728 const double dvy = vy * vyerr;
2729 const double dvz = vz * vzerr;
2730 const double a = vz / vt;
2731 lorerr = sqrt(a * a * (vx * vx * dvx * dvx + vy * vy * dvy * dvy) +
2732 vt * vt * dvz * dvz) /
2733 v2;
2734 lorerr /= lor;
2735 }
2736
2737 // Diffusion coefficients.
2738 dt = sqrt(0.2 * 0.5 * (Magboltz::diflab_.difxx + Magboltz::diflab_.difyy) /
2739 vz) *
2740 1.e-4;
2741 dterr = 0.5 * sqrt(Magboltz::diferl_.dfter * Magboltz::diferl_.dfter +
2742 vzerr * vzerr);
2743 dl = sqrt(0.2 * Magboltz::diflab_.difzz / vz) * 1.e-4;
2744 dlerr = 0.5 * sqrt(Magboltz::diferl_.dfler * Magboltz::diferl_.dfler +
2745 vzerr * vzerr);
2746 // Diffusion tensor.
2747 difftens[0] = 0.2e-4 * Magboltz::diflab_.difzz / vz;
2748 difftens[1] = 0.2e-4 * Magboltz::diflab_.difxx / vz;
2749 difftens[2] = 0.2e-4 * Magboltz::diflab_.difyy / vz;
2750 difftens[3] = 0.2e-4 * Magboltz::diflab_.difxz / vz;
2751 difftens[4] = 0.2e-4 * Magboltz::diflab_.difyz / vz;
2752 difftens[5] = 0.2e-4 * Magboltz::diflab_.difxy / vz;
2753 // Townsend and attachment coefficients.
2754 alpha = Magboltz::ctowns_.alpha;
2755 alphaerr = Magboltz::ctwner_.alper;
2756 eta = Magboltz::ctowns_.att;
2757 etaerr = Magboltz::ctwner_.atter;
2758
2759 // Calculate effective Townsend SST coefficient from TOF results.
2760 if (fabs(Magboltz::tofout_.tofdl) > 0.) {
2761 const double wrzn = 1.e5 * Magboltz::tofout_.tofwr;
2762 const double fc1 = 0.5 * wrzn / Magboltz::tofout_.tofdl;
2763 const double fc2 = (Magboltz::tofout_.ralpha -
2764 Magboltz::tofout_.rattof) * 1.e12 /
2765 Magboltz::tofout_.tofdl;
2766 alphatof = fc1 - sqrt(fc1 * fc1 - fc2);
2767 }
2768 // Print the results.
2769 if (!(m_debug || verbose)) return;
2770 std::cout << m_className << "::RunMagboltz: Results:\n";
2771 printf(" Drift velocity along E: %12.8f cm/ns +/- %5.2f%%\n", vz, vzerr);
2772 printf(" Drift velocity along Bt: %12.8f cm/ns +/- %5.2f%%\n", vx, vxerr);
2773 printf(" Drift velocity along ExB: %12.8f cm/ns +/- %5.2f%%\n", vy, vyerr);
2774 printf(" Lorentz angle: %12.3f degree\n", lor * RadToDegree);
2775 printf(" Longitudinal diffusion: %12.8f cm1/2 +/- %5.2f%%\n", dl, dlerr);
2776 printf(" Transverse diffusion: %12.8f cm1/2 +/- %5.2f%%\n", dt, dterr);
2777 printf(" Townsend coefficient: %12.4f cm-1 +/- %5.2f%%\n", alpha,
2778 alphaerr);
2779 printf(" Attachment coefficient: %12.4f cm-1 +/- %5.2f%%\n", eta,
2780 etaerr);
2781 if (alphatof > 0.) {
2782 printf(" TOF effective Townsend: %12.4f cm-1 (alpha - eta)\n",
2783 alphatof);
2784 }
2785}
static int GetGasNumberMagboltz(const std::string &input)
struct Garfield::Magboltz::@060222016335007357276230325267072324307245075335 velerr_
struct Garfield::Magboltz::@367143224022215036160202274360277173043324053167 diferl_
struct Garfield::Magboltz::@157012065177270115006254216365004043341121046005 vel_
struct Garfield::Magboltz::@245266170301200003106052333313347204035061117342 ratio_
struct Garfield::Magboltz::@010063252323105360345105253327133025274175326133 bfld_
struct Garfield::Magboltz::@367172176335371151066116323064254124322366107157 diflab_
struct Garfield::Magboltz::@043341123353356340117075350116005072335377127162 gasn_
struct Garfield::Magboltz::@152023163236327126074317062211255304316001105103 tofout_
struct Garfield::Magboltz::@261142307054266253060114162042322335114302253152 ctowns_
struct Garfield::Magboltz::@343242004052315221243247367236006336135121373020 setp_
struct Garfield::Magboltz::@243337270271362237030277156013024035163331176157 thrm_
struct Garfield::Magboltz::@301116172155326317070225220046057374317172177377 ctwner_

Referenced by GenerateGasTable().

◆ SetExcitationScaling()

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

Multiply all excitation cross-sections by a uniform scaling factor.

Definition at line 384 of file MediumMagboltz.cc.

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

◆ SetMaxElectronEnergy()

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

Set the highest electron energy to be included in the table of scattering rates.

Definition at line 143 of file MediumMagboltz.cc.

143 {
144 if (e <= Small) {
145 std::cerr << m_className << "::SetMaxElectronEnergy: Invalid energy.\n";
146 return false;
147 }
148 m_eMax = e;
149
150 std::lock_guard<std::mutex> guard(m_mutex);
151 // Determine the energy interval size.
152 m_eStep = std::min(m_eMax, m_eHigh) / Magboltz::nEnergySteps;
153 m_eStepInv = 1. / m_eStep;
154
155 // Force recalculation of the scattering rates table.
156 m_isChanged = true;
157
158 return true;
159}

Referenced by ElectronCollision(), and GetElectronCollisionRate().

◆ SetMaxPhotonEnergy()

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

Set the highest photon energy to be included in the table of scattering rates.

Definition at line 161 of file MediumMagboltz.cc.

161 {
162 if (e <= Small) {
163 std::cerr << m_className << "::SetMaxPhotonEnergy: Invalid energy.\n";
164 return false;
165 }
166 m_eFinalGamma = e;
167
168 // Determine the energy interval size.
169 m_eStepGamma = m_eFinalGamma / nEnergyStepsGamma;
170
171 // Force recalculation of the scattering rates table.
172 m_isChanged = true;
173
174 return true;
175}

Referenced by GetPhotonCollision(), and GetPhotonCollisionRate().

◆ SetSplittingFunctionFlat()

void Garfield::MediumMagboltz::SetSplittingFunctionFlat ( )

Sample the secondary electron energy from a flat distribution.

Definition at line 200 of file MediumMagboltz.cc.

200 {
201 m_useOpalBeaty = false;
202 m_useGreenSawada = false;
203}

◆ SetSplittingFunctionGreenSawada()

void Garfield::MediumMagboltz::SetSplittingFunctionGreenSawada ( )

Sample the secondary electron energy according to the Green-Sawada model.

Definition at line 182 of file MediumMagboltz.cc.

182 {
183 m_useOpalBeaty = false;
184 m_useGreenSawada = true;
185 if (m_isChanged) return;
186
187 bool allset = true;
188 for (unsigned int i = 0; i < m_nComponents; ++i) {
189 if (!m_hasGreenSawada[i]) {
190 if (allset) {
191 std::cout << m_className << "::SetSplittingFunctionGreenSawada:\n";
192 allset = false;
193 }
194 std::cout << " Fit parameters for " << m_gas[i] << " not available.\n"
195 << " Using Opal-Beaty formula instead.\n";
196 }
197 }
198}

◆ SetSplittingFunctionOpalBeaty()

void Garfield::MediumMagboltz::SetSplittingFunctionOpalBeaty ( )

Sample the secondary electron energy according to the Opal-Beaty model.

Definition at line 177 of file MediumMagboltz.cc.

177 {
178 m_useOpalBeaty = true;
179 m_useGreenSawada = false;
180}

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