Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
MediumMagboltz.hh
Go to the documentation of this file.
1#ifndef G_MEDIUM_MAGBOLTZ_9
2#define G_MEDIUM_MAGBOLTZ_9
3
4#include <array>
5
7#include "MediumGas.hh"
8
9namespace Garfield {
10
11/// Interface to %Magboltz (version 11).
12/// - http://magboltz.web.cern.ch/magboltz/
13
14class MediumMagboltz : public MediumGas {
15 public:
16 /// Constructor
18 /// Destructor
19 virtual ~MediumMagboltz() {}
20
21 /// Set the highest electron energy to be included
22 /// in the table of scattering rates.
23 bool SetMaxElectronEnergy(const double e);
24 /// Get the highest electron energy in the table of scattering rates.
25 double GetMaxElectronEnergy() const { return m_eMax; }
26
27 /// Set the highest photon energy to be included
28 /// in the table of scattering rates.
29 bool SetMaxPhotonEnergy(const double e);
30 /// Get the highest photon energy in the table of scattering rates.
31 double GetMaxPhotonEnergy() const { return m_eFinalGamma; }
32
33 /// Switch on/off the automatic adjustment of the max. energy when an
34 /// energy exceeding the present range is requested
35 void EnableEnergyRangeAdjustment(const bool on) { m_useAutoAdjust = on; }
36
37 /// Switch on/off anisotropic scattering (enabled by default)
38 void EnableAnisotropicScattering(const bool on = true) {
39 m_useAnisotropic = on;
40 m_isChanged = true;
41 }
42
43 /// Sample the secondary electron energy according to the Opal-Beaty model.
45 /// Sample the secondary electron energy according to the Green-Sawada model.
47 /// Sample the secondary electron energy from a flat distribution.
49
50 /// Switch on (microscopic) de-excitation handling.
51 void EnableDeexcitation();
52 /// Switch off (microscopic) de-excitation handling.
53 void DisableDeexcitation() { m_useDeexcitation = false; }
54 /// Switch on discrete photoabsorption levels.
56 /// Switch off discrete photoabsorption levels.
57 void DisableRadiationTrapping() { m_useRadTrap = false; }
58
59 bool EnablePenningTransfer(const double r, const double lambda) override;
60 bool EnablePenningTransfer(const double r, const double lambda,
61 std::string gasname) override;
62 void DisablePenningTransfer() override;
63 bool DisablePenningTransfer(std::string gasname) override;
64
65 /// Write the gas cross-section table to a file during the initialisation.
66 void EnableCrossSectionOutput(const bool on = true) { m_useCsOutput = on; }
67
68 /// Multiply all excitation cross-sections by a uniform scaling factor.
69 void SetExcitationScaling(const double r, std::string gasname);
70
71 /// Initialise the table of scattering rates (called internally when a
72 /// collision rate is requested and the gas mixture or other parameters
73 /// have changed).
74 bool Initialise(const bool verbose = false);
75
76 void PrintGas() override;
77
78 /// Get the overall null-collision rate [ns-1].
79 double GetElectronNullCollisionRate(const int band) override;
80 /// Get the (real) collision rate [ns-1] at a given electron energy e [eV].
81 double GetElectronCollisionRate(const double e, const int band) override;
82 /// Get the collision rate [ns-1] for a specific level.
83 double GetElectronCollisionRate(const double e, const unsigned int level,
84 const int band);
85 /// Sample the collision type.
86 bool GetElectronCollision(const double e, int& type, int& level, double& e1,
87 double& dx, double& dy, double& dz,
88 std::vector<std::pair<int, double> >& secondaries,
89 int& ndxc, int& band) override;
90 void ComputeDeexcitation(int iLevel, int& fLevel);
91 unsigned int GetNumberOfDeexcitationProducts() const override {
92 return m_dxcProducts.size();
93 }
94 bool GetDeexcitationProduct(const unsigned int i, double& t, double& s,
95 int& type, double& energy) const override;
96
97 double GetPhotonCollisionRate(const double e) override;
98 bool GetPhotonCollision(const double e, int& type, int& level, double& e1,
99 double& ctheta, int& nsec, double& esec) override;
100
101 /// Reset the collision counters.
103 /// Get the total number of electron collisions.
104 unsigned int GetNumberOfElectronCollisions() const;
105 /// Get the number of collisions broken down by cross-section type.
106 unsigned int GetNumberOfElectronCollisions(unsigned int& nElastic,
107 unsigned int& nIonising,
108 unsigned int& nAttachment,
109 unsigned int& nInelastic,
110 unsigned int& nExcitation,
111 unsigned int& nSuperelastic) const;
112 /// Get the number of cross-section terms.
113 unsigned int GetNumberOfLevels();
114 /// Get detailed information about a given cross-section term i
115 bool GetLevel(const unsigned int i, int& ngas, int& type, std::string& descr,
116 double& e);
117 /// Get the number of collisions for a specific cross-section term.
118 unsigned int GetNumberOfElectronCollisions(const unsigned int level) const;
119
120 /// Get the number of Penning transfers that occured since the last reset.
121 unsigned int GetNumberOfPenningTransfers() const { return m_nPenning; }
122
123 /// Get the total number of photon collisions.
124 unsigned int GetNumberOfPhotonCollisions() const;
125 /// Get number of photon collisions by collision type.
126 unsigned int GetNumberOfPhotonCollisions(unsigned int& nElastic,
127 unsigned int& nIonising,
128 unsigned int& nInelastic) const;
129
130 /// Take the thermal motion of the gas at the selected temperature
131 /// into account in the calculations done by Magboltz.
132 /// By the default, this feature is off (static gas at 0 K).
133 void EnableThermalMotion(const bool on = true) { m_useGasMotion = on; }
134 /// Let Magboltz determine the upper energy limit (this is the default)
135 /// or use the energy limit specified using SetMaxElectronEnergy).
136 void EnableAutoEnergyLimit(const bool on = true) { m_autoEnergyLimit = on; }
137
138 /** Run Magboltz for a given electric field, magnetic field and angle.
139 * \param[in] e electric field
140 * \param[in] b magnetic field
141 * \param[in] btheta angle between electric and magnetic field
142 * \param[in] ncoll number of collisions (in multiples of 10<sup>7</sup>)
143 to be simulated
144 * \param[in] verbose verbosity flag
145 * \param[out] vx,vy,vz drift velocity vector
146 * \param[out] dl,dt diffusion cofficients
147 * \param[out] alpha Townsend cofficient
148 * \param[out] eta attachment cofficient
149 * \param[out] lor Lorentz angle
150 * \param[out] vxerr,vyerr,vzerr errors on drift velocity
151 * \param[out] dlerr,dterr errors on diffusion coefficients
152 * \param[out] alphaerr,etaerr errors on Townsend/attachment coefficients
153 * \param[out] lorerr error on Lorentz angle
154 * \param[out] alphatof effective Townsend coefficient \f$(\alpha - \eta)\f$
155 * calculated using time-of-flight method
156 * \param[out] difftens components of the diffusion tensor (zz, xx, yy, xz, yz, xy)
157 */
158 void RunMagboltz(const double e, const double b, const double btheta,
159 const int ncoll, bool verbose, double& vx, double& vy,
160 double& vz, double& dl, double& dt, double& alpha,
161 double& eta, double& lor, double& vxerr, double& vyerr,
162 double& vzerr, double& dlerr, double& dterr,
163 double& alphaerr, double& etaerr, double& lorerr,
164 double& alphatof, std::array<double, 6>& difftens);
165
166 /// Generate a new gas table (can later be saved to file) by running
167 /// Magboltz for all electric fields, magnetic fields, and
168 /// angles in the currently set grid.
169 void GenerateGasTable(const int numCollisions = 10,
170 const bool verbose = true);
171
172 private:
173 static constexpr int nEnergyStepsLog = 1000;
174 static constexpr int nEnergyStepsGamma = 5000;
175 static constexpr int nCsTypes = 7;
176 static constexpr int nCsTypesGamma = 4;
177
178 static const int DxcTypeRad;
179 static const int DxcTypeCollIon;
180 static const int DxcTypeCollNonIon;
181
182 /// Simulate thermal motion of the gas or not (when running Magboltz).
183 bool m_useGasMotion = false;
184 /// Automatic calculation of the energy limit by Magboltz or not.
185 bool m_autoEnergyLimit = true;
186
187 /// Max. electron energy in the collision rate tables.
188 double m_eMax;
189 /// Energy spacing in the linear part of the collision rate tables.
190 double m_eStep;
191 double m_eHigh, m_eHighLog;
192 double m_lnStep;
193 bool m_useAutoAdjust = true;
194
195 /// Flag enabling/disabling output of cross-section table to file
196 bool m_useCsOutput = false;
197 /// Number of different cross-section types in the current gas mixture
198 unsigned int m_nTerms = 0;
199 /// Recoil energy parameter
200 std::array<double, m_nMaxGases> m_rgas;
201 /// Opal-Beaty-Peterson splitting parameter [eV]
202 std::array<double, Magboltz::nMaxLevels> m_wOpalBeaty;
203 /// Green-Sawada splitting parameters [eV]
204 /// (&Gamma;s, &Gamma;b, Ts, Ta, Tb).
205 std::array<std::array<double, 5>, m_nMaxGases> m_parGreenSawada;
206 std::array<bool, m_nMaxGases> m_hasGreenSawada;
207 /// Sample secondary electron energies using Opal-Beaty parameterisation
208 bool m_useOpalBeaty = true;
209 /// Sample secondary electron energies using Green-Sawada parameterisation
210 bool m_useGreenSawada = false;
211
212 /// Energy loss
213 std::array<double, Magboltz::nMaxLevels> m_energyLoss;
214 /// Cross-section type
215 std::array<int, Magboltz::nMaxLevels> m_csType;
216
217 /// Fluorescence yield
218 std::array<double, Magboltz::nMaxLevels> m_yFluorescence;
219 /// Number of Auger electrons produced in a collision
220 std::array<unsigned int, Magboltz::nMaxLevels> m_nAuger1;
221 std::array<unsigned int, Magboltz::nMaxLevels> m_nAuger2;
222 /// Energy imparted to Auger electrons
223 std::array<double, Magboltz::nMaxLevels> m_eAuger1;
224 std::array<double, Magboltz::nMaxLevels> m_eAuger2;
225 std::array<unsigned int, Magboltz::nMaxLevels> m_nFluorescence;
226 std::array<double, Magboltz::nMaxLevels> m_eFluorescence;
227
228 // Parameters for calculation of scattering angles
229 bool m_useAnisotropic = true;
230 std::vector<std::vector<double> > m_scatPar;
231 std::vector<std::vector<double> > m_scatCut;
232 std::vector<std::vector<double> > m_scatParLog;
233 std::vector<std::vector<double> > m_scatCutLog;
234 std::array<int, Magboltz::nMaxLevels> m_scatModel;
235
236 /// Level description
237 std::vector<std::string> m_description;
238
239 // Total collision frequency
240 std::vector<double> m_cfTot;
241 std::vector<double> m_cfTotLog;
242 /// Null-collision frequency
243 double m_cfNull = 0.;
244 // Collision frequencies
245 std::vector<std::vector<double> > m_cf;
246 std::vector<std::vector<double> > m_cfLog;
247
248 /// Collision counters
249 /// 0: elastic
250 /// 1: ionisation
251 /// 2: attachment
252 /// 3: inelastic
253 /// 4: excitation
254 /// 5: super-elastic
255 std::array<unsigned int, nCsTypes> m_nCollisions;
256 /// Number of collisions for each cross-section term
257 std::vector<unsigned int> m_nCollisionsDetailed;
258
259 // Penning transfer
260 /// Penning transfer probability (by level)
261 std::array<double, Magboltz::nMaxLevels> m_rPenning;
262 /// Mean distance of Penning ionisation (by level)
263 std::array<double, Magboltz::nMaxLevels> m_lambdaPenning;
264 /// Number of Penning ionisations
265 unsigned int m_nPenning = 0;
266
267 // Deexcitation
268 /// Flag enabling/disabling detailed simulation of de-excitation process
269 bool m_useDeexcitation = false;
270 /// Flag enabling/disable radiation trapping
271 /// (absorption of photons discrete excitation lines)
272 bool m_useRadTrap = true;
273
274 struct Deexcitation {
275 // Gas component
276 int gas;
277 // Associated cross-section term
278 int level;
279 // Level description
280 std::string label;
281 // Energy
282 double energy;
283 // Branching ratios
284 std::vector<double> p;
285 // Final levels
286 std::vector<int> final;
287 // Type of transition
288 std::vector<int> type;
289 // Oscillator strength
290 double osc;
291 // Total decay rate
292 double rate;
293 // Doppler broadening
294 double sDoppler;
295 // Pressure broadening
296 double gPressure;
297 // Effective width
298 double width;
299 // Integrated absorption collision rate
300 double cf;
301 };
302 std::vector<Deexcitation> m_deexcitations;
303 // Mapping between deexcitations and cross-section terms.
304 std::array<int, Magboltz::nMaxLevels> m_iDeexcitation;
305
306 // List of de-excitation products
307 struct dxcProd {
308 // Radial spread
309 double s;
310 // Time delay
311 double t;
312 // Type of deexcitation product
313 int type;
314 // Energy of the electron or photon
315 double energy;
316 };
317 std::vector<dxcProd> m_dxcProducts;
318
319 /// Ionisation potentials of each component
320 std::array<double, m_nMaxGases> m_ionPot;
321 /// Minimum ionisation potential
322 double m_minIonPot = -1.;
323
324 // Scaling factor for excitation cross-sections
325 std::array<double, m_nMaxGases> m_scaleExc;
326
327 // Energy spacing of photon collision rates table
328 double m_eFinalGamma, m_eStepGamma;
329 // Number of photon collision cross-section terms
330 unsigned int m_nPhotonTerms = 0;
331 // Total photon collision frequencies
332 std::vector<double> m_cfTotGamma;
333 // Photon collision frequencies
334 std::vector<std::vector<double> > m_cfGamma;
335 std::vector<int> csTypeGamma;
336 // Photon collision counters
337 // 0: elastic
338 // 1: ionisation
339 // 2: inelastic
340 // 3: excitation
341 std::array<unsigned int, nCsTypesGamma> m_nPhotonCollisions;
342
343 int GetGasNumberMagboltz(const std::string& input) const;
344 bool Mixer(const bool verbose = false);
345 void SetupGreenSawada();
346 void SetScatteringParameters(const int model, const double parIn, double& cut,
347 double& parOut) const;
348 void ComputeAngularCut(const double parIn, double& cut, double& parOut) const;
349 void ComputeDeexcitationTable(const bool verbose);
350 void AddPenningDeexcitation(Deexcitation& dxc, const double rate,
351 const double pPenning) {
352 dxc.p.push_back(rate * pPenning);
353 dxc.p.push_back(rate * (1. - pPenning));
354 dxc.type.push_back(DxcTypeCollIon);
355 dxc.type.push_back(DxcTypeCollNonIon);
356 }
357 double RateConstantWK(const double energy, const double osc,
358 const double pacs, const int igas1,
359 const int igas2) const;
360 double RateConstantHardSphere(const double r1, const double r2,
361 const int igas1, const int igas2) const;
362 void ComputeDeexcitationInternal(int iLevel, int& fLevel);
363 bool ComputePhotonCollisionTable(const bool verbose);
364};
365}
366#endif
Base class for gas media.
Definition: MediumGas.hh:15
static constexpr unsigned int m_nMaxGases
Definition: MediumGas.hh:126
bool EnablePenningTransfer(const double r, const double lambda) override
void SetSplittingFunctionGreenSawada()
Sample the secondary electron energy according to the Green-Sawada model.
void EnableAnisotropicScattering(const bool on=true)
Switch on/off anisotropic scattering (enabled by default)
unsigned int GetNumberOfLevels()
Get the number of cross-section terms.
void SetExcitationScaling(const double r, std::string gasname)
Multiply all excitation cross-sections by a uniform scaling factor.
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.
double GetElectronNullCollisionRate(const int band) override
Get the overall null-collision rate [ns-1].
virtual ~MediumMagboltz()
Destructor.
void ResetCollisionCounters()
Reset the collision counters.
void EnableRadiationTrapping()
Switch on discrete photoabsorption levels.
unsigned int GetNumberOfElectronCollisions() const
Get the total number of electron collisions.
double GetPhotonCollisionRate(const double e) override
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 ComputeDeexcitation(int iLevel, int &fLevel)
void SetSplittingFunctionOpalBeaty()
Sample the secondary electron energy according to the Opal-Beaty model.
void GenerateGasTable(const int numCollisions=10, const bool verbose=true)
double GetMaxElectronEnergy() const
Get the highest electron energy in the table of scattering rates.
bool GetDeexcitationProduct(const unsigned int i, double &t, double &s, int &type, double &energy) const override
void SetSplittingFunctionFlat()
Sample the secondary electron energy from a flat distribution.
MediumMagboltz()
Constructor.
unsigned int GetNumberOfPenningTransfers() const
Get the number of Penning transfers that occured since the last reset.
double GetMaxPhotonEnergy() const
Get the highest photon energy in the table of scattering rates.
unsigned int GetNumberOfPhotonCollisions() const
Get the total number of photon collisions.
bool SetMaxPhotonEnergy(const double e)
void DisablePenningTransfer() override
Switch the simulation of Penning transfers off globally.
bool GetPhotonCollision(const double e, int &type, int &level, double &e1, double &ctheta, int &nsec, double &esec) override
unsigned int GetNumberOfDeexcitationProducts() const override
void DisableDeexcitation()
Switch off (microscopic) de-excitation handling.
void EnableCrossSectionOutput(const bool on=true)
Write the gas cross-section table to a file during the initialisation.
void DisableRadiationTrapping()
Switch off discrete photoabsorption levels.
bool GetElectronCollision(const double e, int &type, int &level, double &e1, double &dx, double &dy, double &dz, std::vector< std::pair< int, double > > &secondaries, int &ndxc, int &band) override
Sample the collision type.
bool Initialise(const bool verbose=false)
double GetElectronCollisionRate(const double e, const int band) override
Get the (real) collision rate [ns-1] at a given electron energy e [eV].
bool SetMaxElectronEnergy(const double e)
void EnableThermalMotion(const bool on=true)
void PrintGas() override
Print information about the present gas mixture and available data.
void EnableDeexcitation()
Switch on (microscopic) de-excitation handling.
void EnableEnergyRangeAdjustment(const bool on)
bool m_isChanged
Definition: Medium.hh:540