Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization 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
4#include <array>
5#include <vector>
6
7#include "Garfield/Track.hh"
8
9namespace Garfield {
10class Medium;
13
14class TrackSrim : public Track {
15 public:
17 TrackSrim() : TrackSrim(nullptr) {}
19 TrackSrim(Sensor* sensor);
21 virtual ~TrackSrim() {}
22
24 bool ReadFile(const std::string& file);
26 void Print();
31 void PlotRange();
35
39 void SetModel(const int m) { m_model = m; }
41 int GetModel() const { return m_model; }
42
44 void SetTargetClusterSize(const int n) { m_nsize = n; }
46 int GetTargetClusterSize() const { return m_nsize; }
47
49 void SetClustersMaximum(const int n) { m_maxclusters = n; }
51 int GetClustersMaximum() const { return m_maxclusters; }
52
54 void SetWorkFunction(const double w) { m_work = w; }
56 double GetWorkFunction() const { return m_work; }
58 void SetFanoFactor(const double f) {
59 m_fano = f;
60 m_fset = true;
61 }
62
63 void UnsetFanoFactor() { m_fset = false; }
65 double GetFanoFactor() const { return m_fano; }
67 void SetDensity(const double density) { m_rho = density; }
69 double GetDensity() const { return m_rho; }
71 void SetAtomicMassNumbers(const double a, const double z) {
72 m_a = a;
73 m_z = z;
74 }
75
76 void GetAtomicMassMumbers(double& a, double& z) const {
77 a = m_a;
78 z = m_z;
79 }
80
82 void EnableTransverseStraggling(const bool on = true) {
84 }
85
86 void EnableLongitudinalStraggling(const bool on = true) {
88 }
89
90 struct Cluster {
91 double x, y, z, t;
92 double energy;
93 double kinetic;
94 int n;
95 };
96
97 bool NewTrack(const double x0, const double y0, const double z0,
98 const double t0, const double dx0, const double dy0,
99 const double dz0) override;
100 const std::vector<Cluster>& GetClusters() const { return m_clusters; }
101 bool GetCluster(double& xc, double& yc, double& zc, double& tc, int& nc,
102 double& ec, double& extra);
103
104 protected:
108 bool m_useLongStraggle = false;
109
111 bool m_chargeset = false;
113 double m_qion = 1.;
115 double m_mion = -1.;
117 double m_rho = -1.;
119 double m_work = -1.;
121 bool m_fset = false;
123 double m_fano = -1.;
125 double m_a = -1.;
127 double m_z = -1.;
128
132 std::vector<double> m_ekin;
134 std::vector<double> m_emloss;
136 std::vector<double> m_hdloss;
138 std::vector<double> m_range;
140 std::vector<double> m_transstraggle;
142 std::vector<double> m_longstraggle;
143
145 size_t m_currcluster = 0;
148 unsigned int m_model = 4;
150 int m_nsize = -1;
151 std::vector<Cluster> m_clusters;
152
153 double Xi(const double x, const double beta2, const double edens) const;
154 double DedxEM(const double e) const;
155 double DedxHD(const double e) const;
156 bool PreciseLoss(const double step, const double estart, double& deem,
157 double& dehd) const;
158 bool EstimateRange(const double ekin, const double step,
159 double& stpmax) const;
160 bool SmallestStep(const double ekin, const double edens, double de,
161 double step, double& stpmin);
162 Medium* GetMedium(const std::array<double, 3>& x) const;
163 double Terminate(const std::array<double, 3>& x0,
164 const std::array<double, 3>& v0, const double step0) const;
165 double TerminateBfield(const std::array<double, 3>& x0,
166 const std::array<double, 3>& v0, const double dt0,
167 const double vmag) const;
168 double RndmEnergyLoss(const double ekin, const double de, const double step,
169 const double edens) const;
170};
171} // namespace Garfield
172
173#endif
std::vector< double > m_longstraggle
Longitudinal straggling [cm].
Definition TrackSrim.hh:142
void SetClustersMaximum(const int n)
Set the max. number of clusters on a track.
Definition TrackSrim.hh:49
void GetAtomicMassMumbers(double &a, double &z) const
Get A and Z of the target medium.
Definition TrackSrim.hh:76
void Print()
Print the energy loss table.
double Xi(const double x, const double beta2, const double edens) const
double m_qion
Charge of the projectile.
Definition TrackSrim.hh:113
double m_z
Effective Z of the target.
Definition TrackSrim.hh:127
std::vector< double > m_hdloss
Hadronic energy loss [MeV cm2/g].
Definition TrackSrim.hh:136
void SetDensity(const double density)
Set the density [g/cm3] of the target medium.
Definition TrackSrim.hh:67
const std::vector< Cluster > & GetClusters() const
Definition TrackSrim.hh:100
std::vector< double > m_range
Projected range [cm].
Definition TrackSrim.hh:138
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
Calculate a new track starting from (x0, y0, z0) at time t0 in direction (dx0, dy0,...
double DedxHD(const double e) const
bool EstimateRange(const double ekin, const double step, double &stpmax) const
std::vector< Cluster > m_clusters
Definition TrackSrim.hh:151
size_t m_currcluster
Index of the next cluster to be returned.
Definition TrackSrim.hh:145
double GetFanoFactor() const
Get the Fano factor.
Definition TrackSrim.hh:65
void EnableTransverseStraggling(const bool on=true)
Simulate transverse straggling (default: on).
Definition TrackSrim.hh:82
std::vector< double > m_emloss
EM energy loss [MeV cm2/g].
Definition TrackSrim.hh:134
void SetTargetClusterSize(const int n)
Specify how many electrons should be grouped to a cluster.
Definition TrackSrim.hh:44
void PlotStraggling()
Plot the transverse and longitudinal straggling as function of the projectile energy.
bool m_chargeset
Has the charge been defined?
Definition TrackSrim.hh:111
int GetModel() const
Get the fluctuation model.
Definition TrackSrim.hh:41
double m_fano
Fano factor [-] of the target.
Definition TrackSrim.hh:123
void SetModel(const int m)
Set the fluctuation model (0 = none, 1 = Landau, 2 = Vavilov, 3 = Gaussian, 4 = Combined).
Definition TrackSrim.hh:39
void EnableLongitudinalStraggling(const bool on=true)
Simulate longitudinal straggling (default: off).
Definition TrackSrim.hh:86
TrackSrim()
Default constructor.
Definition TrackSrim.hh:17
Medium * GetMedium(const std::array< double, 3 > &x) const
double RndmEnergyLoss(const double ekin, const double de, const double step, const double edens) const
double m_mion
Mass [MeV] of the projectile.
Definition TrackSrim.hh:115
bool PreciseLoss(const double step, const double estart, double &deem, double &dehd) const
void UnsetFanoFactor()
Use the default Fano factor.
Definition TrackSrim.hh:63
void PlotEnergyLoss()
Plot the electromagnetic, hadronic, and total energy loss as function of the projectile energy.
double m_work
Work function [eV] of the target.
Definition TrackSrim.hh:119
unsigned int m_model
Fluctuation model (0 = none, 1 = Landau, 2 = Vavilov, 3 = Gaussian, 4 = Combined)
Definition TrackSrim.hh:148
bool GetCluster(double &xc, double &yc, double &zc, double &tc, int &nc, double &ec, double &extra)
bool m_useTransStraggle
Include transverse straggling.
Definition TrackSrim.hh:106
int m_nsize
Targeted cluster size.
Definition TrackSrim.hh:150
void SetAtomicMassNumbers(const double a, const double z)
Set A and Z of the target medium.
Definition TrackSrim.hh:71
std::vector< double > m_ekin
Energy in energy loss table [MeV].
Definition TrackSrim.hh:132
std::vector< double > m_transstraggle
Transverse straggling [cm].
Definition TrackSrim.hh:140
double Terminate(const std::array< double, 3 > &x0, const std::array< double, 3 > &v0, const double step0) const
double GetDensity() const
Get the density [g/cm3] of the target medium.
Definition TrackSrim.hh:69
int GetClustersMaximum() const
Retrieve the max. number of clusters on a track.
Definition TrackSrim.hh:51
double m_rho
Mass density [g/cm3] of the target.
Definition TrackSrim.hh:117
virtual ~TrackSrim()
Destructor.
Definition TrackSrim.hh:21
TrackSrim(Sensor *sensor)
Constructor.
void SetWorkFunction(const double w)
Set the W value [eV].
Definition TrackSrim.hh:54
bool m_useLongStraggle
Include longitudinal straggling.
Definition TrackSrim.hh:108
double m_a
Effective A of the target.
Definition TrackSrim.hh:125
void PlotRange()
Plot the projected range as function of the projectile energy.
bool m_fset
Has the Fano factor been set?
Definition TrackSrim.hh:121
int m_maxclusters
Maximum number of clusters allowed (infinite if 0)
Definition TrackSrim.hh:130
double DedxEM(const double e) const
bool SmallestStep(const double ekin, const double edens, double de, double step, double &stpmin)
void SetFanoFactor(const double f)
Set the Fano factor.
Definition TrackSrim.hh:58
double GetWorkFunction() const
Get the W value [eV].
Definition TrackSrim.hh:56
double TerminateBfield(const std::array< double, 3 > &x0, const std::array< double, 3 > &v0, const double dt0, const double vmag) const
int GetTargetClusterSize() const
Retrieve the target cluster size.
Definition TrackSrim.hh:46
bool ReadFile(const std::string &file)
Load data from a SRIM file.
Track()=delete
Default constructor.
int n
Number of electrons in this cluster.
Definition TrackSrim.hh:94
double kinetic
Ion energy when cluster was created.
Definition TrackSrim.hh:93
double t
Cluster location and time.
Definition TrackSrim.hh:91
double energy
Energy spent to make the cluster.
Definition TrackSrim.hh:92