Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
HeedDeltaElectronCS.h
Go to the documentation of this file.
1#ifndef HEEDDELTAELECTRONCS_H
2#define HEEDDELTAELECTRONCS_H
3
4#include <vector>
5
10
11#define USE_MEAN_COEF // new variant, means that used mean(1-cos(theta))
12// for low angle scattering. It seems to be more straightforwad and
13// may be very slightly more precise than the old variant
14// that is use of sqrt( mean ( square (1-cos(theta)) ) )
15
16namespace Heed {
17
18/// Cross sections and various parameters for delta-electron transport.
19/// 2003, I. Smirnov
20
22 public:
23 /// Default constructor
25 /// Constructor
27 ElElasticScatLowSigma* feesls, PairProd* fpairprod,
28 int fsruth = 2, double fmlambda = 0.001 * 4.0e-3,
29 double fmthetac = 0.1);
30
31 double get_sigma(double energy, double nscat) const;
32 // copy of similar thing from ElElasticScatLowSigma
33
34 void print(std::ostream& file, int l) const;
35 HeedDeltaElectronCS* copy() const { return new HeedDeltaElectronCS(*this); }
36
37 static constexpr long q_angular_mesh = 50;
38 static constexpr double low_cut_angle_deg = 20.;
39
40 HeedMatterDef* hmd = nullptr;
41 ElElasticScat* ees = nullptr;
43 PairProd* pairprod = nullptr; // in eV
44 /// Table of velocities
45 std::vector<double> beta;
46 /// Table of momenta [MeV/c]
47 std::vector<double> momentum;
48
49 /// Energy losses [MeV/cm] according to very advanced formula with
50 /// Bethe-Bloch and density effect as in GEANT3,
51 /// corresponding to centers of intervals of the common mesh.
52 std::vector<double> eLoss;
53
54 /// Mminimum mean length of range, multiplied by density.
55 /// sm*gr/sm**3 = gr/sm**2
56 double mlambda;
57 /// Path length [cm] at the energy values of the mesh.
58 /// For sruth == 2 mean free path.
59 std::vector<double> lambda;
60 /// Path length for low angle scatterings.
61 /// This is without multiplication used for coef_low_sigma (cm).
62 std::vector<double> low_lambda;
63
64 /// Formula to use.
65 /// 0 - usual multiple scattering formula
66 /// 1 - Rutherford cross-section
67 /// 2 - precise theory?
68 int sruth;
69 /// Minimum threshold turn angle.
70 /// For Rutherford: the interactions with less angle will not be taken
71 /// into account. The actual threshold angle can be larger.
72 /// The second restriction is going from restriction of atomic shell.
73 /// The third one is from mlamBdel.
74 /// For usual multiple scattering:
75 /// Assuming that sigma = mTetacBdel the path lengt is calculated.
76 /// If mlamBdel/density is less then the last is used.
77 double mthetac;
78
79 /// Angular mesh, centers, angles in degrees.
80 std::vector<double> angular_mesh_c;
81
82 // index - energy mesh, angles in degrees
83 std::vector<PointsRan> angular_points_ran;
84 // index - energy mesh
85 std::vector<PointsRan> low_angular_points_ran;
86
87// introduce new name mean_... in order to avoid errors at replacement
88#ifdef USE_MEAN_COEF
89 std::vector<double> mean_coef_low_sigma;
90#else
91 std::vector<double> coef_low_sigma;
92// index - energy, but (attention!) for mesh defined in ElElasticScat
93// this is sigma of cos(theta) - 1.0, supposing that center is 1.0
94#endif
95};
96}
97
98#endif
void print(std::ostream &file, int l) const
HeedDeltaElectronCS * copy() const
double get_sigma(double energy, double nscat) const
std::vector< PointsRan > low_angular_points_ran
ElElasticScatLowSigma * eesls
HeedDeltaElectronCS()
Default constructor.
std::vector< double > beta
Table of velocities.
static constexpr long q_angular_mesh
static constexpr double low_cut_angle_deg
std::vector< double > momentum
Table of momenta [MeV/c].
std::vector< double > eLoss
std::vector< double > angular_mesh_c
Angular mesh, centers, angles in degrees.
std::vector< double > lambda
std::vector< PointsRan > angular_points_ran
std::vector< double > mean_coef_low_sigma
std::vector< double > low_lambda
Definition: BGMesh.cpp:6