Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
TrackSrim.hh
Go to the documentation of this file.
1#ifndef G_TRACK_SRIM_H
2#define G_TRACK_SRIM_H
3#include <vector>
4#include "Track.hh"
5
6namespace Garfield {
7
8/// Generate tracks based on SRIM energy loss, range and straggling tables.
9/// - http://www.srim.org
10
11class TrackSrim : public Track {
12
13 public:
14 /// Constructor
15 TrackSrim();
16 /// Destructor
17 virtual ~TrackSrim() {}
18
19 /// Set/get the W value [eV].
20 void SetWorkFunction(const double w) { m_work = w; }
21 double GetWorkFunction() const { return m_work; }
22 /// Set/get the Fano factor.
23 void SetFanoFactor(const double f) { m_fano = f; }
24 double GetFanoFactor() const { return m_fano; }
25 /// Set/get the density [g/cm3] of the target medium.
26 void SetDensity(const double density) { m_density = density; }
27 double GetDensity() const { return m_density; }
28 /// Set/get A and Z of the target medium.
29 void SetAtomicMassNumbers(const double a, const double z) {
30 m_a = a;
31 m_z = z;
32 }
33 void GetAtomicMassMumbers(double& a, double& z) const {
34 a = m_a;
35 z = m_z;
36 }
37
38 /// Set/get the fluctuation model
39 /// (0 = none, 1 = Landau, 2 = Vavilov, 3 = Gaussian, 4 = Combined)
40 void SetModel(const int m) { m_model = m; }
41 int GetModel() const { return m_model; }
42
47
50
51 void SetTargetClusterSize(const int n) { m_nsize = n; }
52 int GetTargetClusterSize() const { return m_nsize; }
53
54 void SetClustersMaximum(const int n) { m_maxclusters = n; }
55 int SetClustersMaximum() const { return m_maxclusters; }
56
57 bool ReadFile(const std::string& file);
58 void Print();
59 void PlotEnergyLoss();
60 void PlotRange();
61 void PlotStraggling();
62
63 virtual bool NewTrack(const double x0, const double y0, const double z0,
64 const double t0, const double dx0, const double dy0,
65 const double dz0);
66 virtual bool GetCluster(double& xcls, double& ycls, double& zcls,
67 double& tcls, int& n, double& e, double& extra);
68
69 protected:
70 /// Use precise Vavilov generator
72 /// Include transverse straggling
74 /// Include longitudinal straggling
76
77 /// Charge gas been defined
79 /// Density [g/cm3]
80 double m_density;
81 /// Work function [eV]
82 double m_work;
83 /// Fano factor [-]
84 double m_fano;
85 /// Charge of ion
86 double m_q;
87 /// Mass of ion [MeV]
88 double m_mass;
89 /// A and Z of the gas
90 double m_a;
91 double m_z;
92
93 /// Maximum number of clusters allowed (infinite if 0)
95 /// Energy in energy loss table [MeV]
96 std::vector<double> m_ekin;
97 /// EM energy loss [MeV cm2/g]
98 std::vector<double> m_emloss;
99 /// Hadronic energy loss [MeV cm2/g]
100 std::vector<double> m_hdloss;
101 /// Projected range [cm]
102 std::vector<double> m_range;
103 /// Transverse straggling [cm]
104 std::vector<double> m_transstraggle;
105 /// Longitudinal straggling [cm]
106 std::vector<double> m_longstraggle;
107
108 /// Index of the next cluster to be returned
109 unsigned int m_currcluster;
110 /// Fluctuation model (0 = none, 1 = Landau, 2 = Vavilov,
111 /// 3 = Gaussian, 4 = Combined)
112 unsigned int m_model;
113 /// Targeted cluster size
115 struct cluster {
116 double x, y, z, t; // Cluster location and time
117 double ec; // Energy spent to make the clusterec
118 double kinetic; // Ion energy when cluster was created
119 int electrons; // Number of electrons in this cluster
120 };
121 std::vector<cluster> m_clusters;
122
123 double DedxEM(const double e) const;
124 double DedxHD(const double e) const;
125 bool PreciseLoss(const double step, const double estart, double& deem,
126 double& dehd) const;
127 bool EstimateRange(const double ekin, const double step, double& stpmax);
128 bool SmallestStep(double ekin, double de, double step, double& stpmin);
129
130 double RndmEnergyLoss(const double ekin, const double de,
131 const double step) const;
132};
133}
134
135#endif
std::vector< double > m_longstraggle
Longitudinal straggling [cm].
Definition: TrackSrim.hh:106
double RndmEnergyLoss(const double ekin, const double de, const double step) const
Definition: TrackSrim.cc:1226
double m_q
Charge of ion.
Definition: TrackSrim.hh:86
void SetClustersMaximum(const int n)
Definition: TrackSrim.hh:54
void GetAtomicMassMumbers(double &a, double &z) const
Definition: TrackSrim.hh:33
virtual bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0)
Definition: TrackSrim.cc:698
unsigned int m_currcluster
Index of the next cluster to be returned.
Definition: TrackSrim.hh:109
std::vector< double > m_hdloss
Hadronic energy loss [MeV cm2/g].
Definition: TrackSrim.hh:100
void SetDensity(const double density)
Set/get the density [g/cm3] of the target medium.
Definition: TrackSrim.hh:26
void DisableTransverseStraggling()
Definition: TrackSrim.hh:44
std::vector< double > m_range
Projected range [cm].
Definition: TrackSrim.hh:102
double DedxHD(const double e) const
Definition: TrackSrim.cc:517
std::vector< cluster > m_clusters
Definition: TrackSrim.hh:121
void EnableLongitudinalStraggling()
Definition: TrackSrim.hh:45
double GetFanoFactor() const
Definition: TrackSrim.hh:24
std::vector< double > m_emloss
EM energy loss [MeV cm2/g].
Definition: TrackSrim.hh:98
void SetTargetClusterSize(const int n)
Definition: TrackSrim.hh:51
virtual bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra)
Definition: TrackSrim.cc:1330
double m_mass
Mass of ion [MeV].
Definition: TrackSrim.hh:88
bool m_chargeset
Charge gas been defined.
Definition: TrackSrim.hh:78
bool SmallestStep(double ekin, double de, double step, double &stpmin)
Definition: TrackSrim.cc:1052
int GetModel() const
Definition: TrackSrim.hh:41
int SetClustersMaximum() const
Definition: TrackSrim.hh:55
double m_fano
Fano factor [-].
Definition: TrackSrim.hh:84
void SetModel(const int m)
Definition: TrackSrim.hh:40
TrackSrim()
Constructor.
Definition: TrackSrim.cc:122
bool PreciseLoss(const double step, const double estart, double &deem, double &dehd) const
Definition: TrackSrim.cc:521
bool m_precisevavilov
Use precise Vavilov generator.
Definition: TrackSrim.hh:71
double m_work
Work function [eV].
Definition: TrackSrim.hh:82
unsigned int m_model
Definition: TrackSrim.hh:112
double m_density
Density [g/cm3].
Definition: TrackSrim.hh:80
void EnablePreciseVavilov()
Definition: TrackSrim.hh:48
bool m_useTransStraggle
Include transverse straggling.
Definition: TrackSrim.hh:73
int m_nsize
Targeted cluster size.
Definition: TrackSrim.hh:114
void SetAtomicMassNumbers(const double a, const double z)
Set/get A and Z of the target medium.
Definition: TrackSrim.hh:29
std::vector< double > m_ekin
Energy in energy loss table [MeV].
Definition: TrackSrim.hh:96
std::vector< double > m_transstraggle
Transverse straggling [cm].
Definition: TrackSrim.hh:104
double GetDensity() const
Definition: TrackSrim.hh:27
virtual ~TrackSrim()
Destructor.
Definition: TrackSrim.hh:17
void SetWorkFunction(const double w)
Set/get the W value [eV].
Definition: TrackSrim.hh:20
void DisableLongitudinalStraggling()
Definition: TrackSrim.hh:46
bool EstimateRange(const double ekin, const double step, double &stpmax)
Definition: TrackSrim.cc:602
bool m_useLongStraggle
Include longitudinal straggling.
Definition: TrackSrim.hh:75
double m_a
A and Z of the gas.
Definition: TrackSrim.hh:90
void EnableTransverseStraggling()
Definition: TrackSrim.hh:43
int m_maxclusters
Maximum number of clusters allowed (infinite if 0)
Definition: TrackSrim.hh:94
double DedxEM(const double e) const
Definition: TrackSrim.cc:513
void SetFanoFactor(const double f)
Set/get the Fano factor.
Definition: TrackSrim.hh:23
double GetWorkFunction() const
Definition: TrackSrim.hh:21
void DisablePreciseVavilov()
Definition: TrackSrim.hh:49
int GetTargetClusterSize() const
Definition: TrackSrim.hh:52
bool ReadFile(const std::string &file)
Definition: TrackSrim.cc:140
Abstract base class for track generation.
Definition: Track.hh:14