Garfield++ v2r0
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 virtual void print(std::ostream& file, int l) const;
35 virtual HeedDeltaElectronCS* copy() const {
36 return new HeedDeltaElectronCS(*this);
37 }
38
39 static const long q_angular_mesh = 50;
40 static const double low_cut_angle_deg;
41
46 /// Table of velocities
47 std::vector<double> beta;
48 /// Table of momenta [MeV/c]
49 std::vector<double> momentum;
50
51 /// Energy losses [MeV/cm] according to very advanced formula with
52 /// Bethe-Bloch and density effect as in GEANT3,
53 /// corresponding to centers of intervals of the common mesh.
54 std::vector<double> eLoss;
55
56 /// Mminimum mean length of range, multiplied by density.
57 /// sm*gr/sm**3 = gr/sm**2
58 double mlambda;
59 /// Path length [cm] at the energy values of the mesh.
60 /// For sruth == 2 mean free path.
61 std::vector<double> lambda;
62 /// Path length for low angle scatterings.
63 /// This is without multiplication used for coef_low_sigma (cm).
64 std::vector<double> low_lambda;
65
66 /// Formula to use.
67 /// 0 - usual multiple scattering formula
68 /// 1 - Rutherford cross-section
69 /// 2 - precise theory?
70 int sruth;
71 /// Minimum threshold turn angle.
72 /// For Rutherford: the interactions with less angle will not be taken
73 /// into account. The actual threshold angle can be larger.
74 /// The second restriction is going from restriction of atomic shell.
75 /// The third one is from mlamBdel.
76 /// For usual multiple scattering:
77 /// Assuming that sigma = mTetacBdel the path lengt is calculated.
78 /// If mlamBdel/density is less then the last is used.
79 double mthetac;
80
81 /// Angular mesh, centers, angles in degrees.
82 std::vector<double> angular_mesh_c;
83
84 // index - energy mesh, angles in degrees
85 std::vector<PointsRan> angular_points_ran;
86 // index - energy mesh
87 std::vector<PointsRan> low_angular_points_ran;
88
89// introduce new name mean_... in order to avoid errors at replacement
90#ifdef USE_MEAN_COEF
91 std::vector<double> mean_coef_low_sigma;
92#else
93 std::vector<double> coef_low_sigma;
94// index - energy, but (attention!) for mesh defined in ElElasticScat
95// this is sigma of cos(theta) - 1.0, supposing that center is 1.0
96#endif
97};
98}
99
100#endif
PassivePtr< ElElasticScat > ees
virtual void print(std::ostream &file, int l) const
virtual HeedDeltaElectronCS * copy() const
double get_sigma(double energy, double nscat) const
PassivePtr< HeedMatterDef > hmd
std::vector< PointsRan > low_angular_points_ran
HeedDeltaElectronCS()
Default constructor.
std::vector< double > beta
Table of velocities.
PassivePtr< PairProd > pairprod
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
PassivePtr< ElElasticScatLowSigma > eesls
static const double low_cut_angle_deg
std::vector< double > mean_coef_low_sigma
std::vector< double > low_lambda
static const long q_angular_mesh
Definition: BGMesh.cpp:5