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