Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::TrackTrim Class Reference

#include <TrackTrim.hh>

+ Inheritance diagram for Garfield::TrackTrim:

Classes

struct  Cluster
 

Public Member Functions

 TrackTrim ()
 Constructor.
 
virtual ~TrackTrim ()
 Destructor.
 
bool ReadFile (const std::string &file, const unsigned int nIons=0, const unsigned int nSkip=0)
 Load data from an EXYZ.txt file.
 
void Print ()
 Print a summary of the available TRIM data.
 
void SetWorkFunction (const double w)
 Set the W value [eV].
 
double GetWorkFunction () const
 Get the W value [eV].
 
void SetFanoFactor (const double f)
 Set the Fano factor.
 
double GetFanoFactor () const
 Get the Fano factor.
 
bool NewTrack (const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
 
bool GetCluster (double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra) override
 
- Public Member Functions inherited from Garfield::Track
 Track ()
 Constructor.
 
virtual ~Track ()
 Destructor.
 
virtual void SetParticle (const std::string &part)
 
void SetEnergy (const double e)
 Set the particle energy.
 
void SetBetaGamma (const double bg)
 Set the relative momentum of the particle.
 
void SetBeta (const double beta)
 Set the speed ( $\beta = v/c$) of the particle.
 
void SetGamma (const double gamma)
 Set the Lorentz factor of the particle.
 
void SetMomentum (const double p)
 Set the particle momentum.
 
void SetKineticEnergy (const double ekin)
 Set the kinetic energy of the particle.
 
double GetEnergy () const
 Return the particle energy.
 
double GetBetaGamma () const
 Return the $\beta\gamma$ of the projectile.
 
double GetBeta () const
 Return the speed ( $\beta = v/c$) of the projectile.
 
double GetGamma () const
 Return the Lorentz factor of the projectile.
 
double GetMomentum () const
 Return the particle momentum.
 
double GetKineticEnergy () const
 Return the kinetic energy of the projectile.
 
double GetCharge () const
 Get the charge of the projectile.
 
double GetMass () const
 Get the mass [eV / c2] of the projectile.
 
void SetSensor (Sensor *s)
 Set the sensor through which to transport the particle.
 
virtual bool NewTrack (const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0)=0
 
virtual bool GetCluster (double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra)=0
 
virtual double GetClusterDensity ()
 
virtual double GetStoppingPower ()
 Get the stopping power (mean energy loss [eV] per cm).
 
void EnablePlotting (ViewDrift *viewer)
 Switch on plotting.
 
void DisablePlotting ()
 Switch off plotting.
 
void EnableDebugging ()
 Switch on debugging messages.
 
void DisableDebugging ()
 Switch off debugging messages.
 

Protected Member Functions

void AddIon (const std::vector< float > &x, const std::vector< float > &y, const std::vector< float > &z, const std::vector< float > &dedx, const std::vector< float > &ekin)
 
- Protected Member Functions inherited from Garfield::Track
void PlotNewTrack (const double x0, const double y0, const double z0)
 
void PlotCluster (const double x0, const double y0, const double z0)
 

Protected Attributes

double m_work = -1.
 Work function [eV].
 
double m_fano = -1.
 Fano factor [-].
 
double m_ekin = 0.
 Projectile energy [eV].
 
std::vector< std::vector< std::array< float, 5 > > > m_ions
 List of tracks.
 
size_t m_ion = 0
 Index of the current track.
 
std::vector< Clusterm_clusters
 Clusters on the current track.
 
size_t m_cluster = 0
 Index of the next cluster to be returned.
 
- Protected Attributes inherited from Garfield::Track
std::string m_className = "Track"
 
double m_q = -1.
 
int m_spin = 1
 
double m_mass
 
double m_energy = 0.
 
double m_beta2
 
bool m_isElectron = false
 
std::string m_particleName = "mu-"
 
Sensorm_sensor = nullptr
 
bool m_isChanged = true
 
ViewDriftm_viewer = nullptr
 
bool m_debug = false
 
size_t m_plotId = 0
 

Detailed Description

Generate tracks based on TRIM output files.

Definition at line 14 of file TrackTrim.hh.

Constructor & Destructor Documentation

◆ TrackTrim()

Garfield::TrackTrim::TrackTrim ( )

Constructor.

Definition at line 28 of file TrackTrim.cc.

28: Track() { m_className = "TrackTrim"; }
std::string m_className
Definition: Track.hh:102
Track()
Constructor.
Definition: Track.cc:13

◆ ~TrackTrim()

virtual Garfield::TrackTrim::~TrackTrim ( )
inlinevirtual

Destructor.

Definition at line 19 of file TrackTrim.hh.

19{}

Member Function Documentation

◆ AddIon()

void Garfield::TrackTrim::AddIon ( const std::vector< float > &  x,
const std::vector< float > &  y,
const std::vector< float > &  z,
const std::vector< float > &  dedx,
const std::vector< float > &  ekin 
)
protected

Definition at line 117 of file TrackTrim.cc.

121 {
122
123 const size_t nPoints = x.size();
124 if (nPoints < 2) return;
125 std::vector<std::array<float, 5> > path;
126 for (size_t i = 0; i < nPoints; ++i) {
127 float eloss = 0.;
128 if (i < nPoints - 1) {
129 const float dx = x[i + 1] - x[i];
130 const float dy = y[i + 1] - y[i];
131 const float dz = z[i + 1] - z[i];
132 const float step = sqrt(dx * dx + dy * dy + dz * dz);
133 if (i == 0 && dedx[i] > 10. * dedx[i + 1]) {
134 eloss = step * dedx[i + 1];
135 } else {
136 eloss = step * dedx[i];
137 }
138 const float dekin = ekin[i] - ekin[i + 1];
139 if (dekin > 0.) eloss = std::min(eloss, dekin);
140 }
141 path.push_back({x[i], y[i], z[i], eloss, ekin[i]});
142 }
143 m_ions.push_back(path);
144}
std::vector< std::vector< std::array< float, 5 > > > m_ions
List of tracks.
Definition: TrackTrim.hh:51
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

Referenced by ReadFile().

◆ GetCluster()

bool Garfield::TrackTrim::GetCluster ( double &  xcls,
double &  ycls,
double &  zcls,
double &  tcls,
int &  n,
double &  e,
double &  extra 
)
overridevirtual

Get the next "cluster" (ionising collision of the charged particle).

Parameters
xcls,ycls,zclscoordinates of the collision
tclstime of the collision
nnumber of electrons produced
edeposited energy
extraadditional information (not always implemented)

Implements Garfield::Track.

Definition at line 272 of file TrackTrim.cc.

273 {
274 if (m_debug) {
275 std::cout << m_className << "::GetCluster: Cluster " << m_cluster
276 << " of " << m_clusters.size() << "\n";
277 }
278 // Stop if we have exhausted the list of clusters.
279 if (m_cluster >= m_clusters.size()) return false;
280
281 const auto& cluster = m_clusters[m_cluster];
282 xcls = cluster.x;
283 ycls = cluster.y;
284 zcls = cluster.z;
285 tcls = cluster.t;
286
287 n = cluster.electrons;
288 e = cluster.ec;
289 extra = cluster.ekin;
290 // Move to the next cluster.
291 ++m_cluster;
292 return true;
293}
size_t m_cluster
Index of the next cluster to be returned.
Definition: TrackTrim.hh:64
std::vector< Cluster > m_clusters
Clusters on the current track.
Definition: TrackTrim.hh:62
bool m_debug
Definition: Track.hh:118

◆ GetFanoFactor()

double Garfield::TrackTrim::GetFanoFactor ( ) const
inline

Get the Fano factor.

Definition at line 34 of file TrackTrim.hh.

34{ return m_fano; }
double m_fano
Fano factor [-].
Definition: TrackTrim.hh:46

◆ GetWorkFunction()

double Garfield::TrackTrim::GetWorkFunction ( ) const
inline

Get the W value [eV].

Definition at line 30 of file TrackTrim.hh.

30{ return m_work; }
double m_work
Work function [eV].
Definition: TrackTrim.hh:44

◆ NewTrack()

bool Garfield::TrackTrim::NewTrack ( const double  x0,
const double  y0,
const double  z0,
const double  t0,
const double  dx0,
const double  dy0,
const double  dz0 
)
overridevirtual

Calculate a new track starting from (x0, y0, z0) at time t0 in direction (dx0, dy0, dz0).

Implements Garfield::Track.

Definition at line 159 of file TrackTrim.cc.

160 {
161
162 // TRMGEN - Generates TRIM clusters
163
164 if (m_ions.empty()) {
165 std::cerr << m_className << "::NewTrack: No TRIM data present.\n";
166 return false;
167 }
168
169 if (m_ion >= m_ions.size()) {
170 // Rewind.
171 std::cout << m_className << "::NewTrack: Rewinding.\n";
172 m_ion = 0;
173 }
174
175 // Verify that a sensor has been set.
176 if (!m_sensor) {
177 std::cerr << m_className << "::NewTrack: Sensor is not defined.\n";
178 return false;
179 }
180
181 // Normalise and store the direction.
182 const double d0 = sqrt(dx0 * dx0 + dy0 * dy0 + dz0 * dz0);
183 double dx = dx0;
184 double dy = dy0;
185 double dz = dz0;
186 if (d0 < Small) {
187 if (m_debug) {
188 std::cout << m_className << "::NewTrack: Randomizing initial direction.\n";
189 }
190 // Null vector. Sample the direction isotropically.
191 RndmDirection(dx, dy, dz);
192 } else {
193 // Normalise the direction vector.
194 dx /= d0;
195 dy /= d0;
196 dz /= d0;
197 }
198 const double dt = sqrt(dx * dx + dy * dy);
199 double phi = 0.;
200 double theta = 0.;
201 if (dt > 0.) {
202 phi = atan2(dy, dx);
203 theta = atan2(dz, dt);
204 } else {
205 theta = dz < 0. ? -HalfPi : HalfPi;
206 }
207 const double ctheta = cos(theta);
208 const double stheta = sin(theta);
209 const double cphi = cos(phi);
210 const double sphi = sin(phi);
211
212 // Make sure all necessary parameters have been set.
213 if (m_work < Small) {
214 std::cerr << m_className << "::NewTrack: Work function not set.\n";
215 return false;
216 }
217
218 // Plot.
219 if (m_viewer) PlotNewTrack(x0, y0, z0);
220
221 // Reset the cluster count.
222 m_cluster = 0;
223 m_clusters.clear();
224
225 // Pool of unused energy
226 double epool = 0.0;
227
228 const auto& path = m_ions[m_ion];
229 const size_t nPoints = path.size();
230 for (size_t i = 1; i < nPoints; ++i) {
231 const double x = path[i][0];
232 const double y = path[i][1];
233 const double z = path[i][2];
234 Cluster cluster;
235 cluster.x = x0 + cphi * ctheta * x - sphi * y - cphi * stheta * z;
236 cluster.y = y0 + sphi * ctheta * x + cphi * y - sphi * stheta * z;
237 cluster.z = z0 + stheta * x + ctheta * z;
238 // Is this point inside an ionisable medium?
239 Medium* medium = nullptr;
240 if (!m_sensor->GetMedium(cluster.x, cluster.y, cluster.z, medium)) {
241 continue;
242 }
243 if (!medium || !medium->IsIonisable()) continue;
244 cluster.t = t0;
245 double eloss = path[i - 1][3];
246 if (m_fano < Small) {
247 // No fluctuations.
248 cluster.electrons = int((eloss + epool) / m_work);
249 cluster.ec = m_work * cluster.electrons;
250 } else {
251 double ec = eloss + epool;
252 cluster.electrons = 0;
253 cluster.ec = 0.0;
254 while (true) {
255 const double er = RndmHeedWF(m_work, m_fano);
256 if (er > ec) break;
257 cluster.electrons++;
258 cluster.ec += er;
259 ec -= er;
260 }
261 }
262 cluster.ekin = path[i][4];
263 epool += eloss - cluster.ec;
264 m_clusters.push_back(std::move(cluster));
265 if (m_viewer) PlotCluster(cluster.x, cluster.y, cluster.z);
266 }
267 // Move to the next ion in the list.
268 ++m_ion;
269 return true;
270}
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:159
size_t m_ion
Index of the current track.
Definition: TrackTrim.hh:53
Sensor * m_sensor
Definition: Track.hh:112
ViewDrift * m_viewer
Definition: Track.hh:116
void PlotCluster(const double x0, const double y0, const double z0)
Definition: Track.cc:192
void PlotNewTrack(const double x0, const double y0, const double z0)
Definition: Track.cc:187
double RndmHeedWF(const double w, const double f)
Definition: Random.cc:699
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:107
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384

◆ Print()

void Garfield::TrackTrim::Print ( )

Print a summary of the available TRIM data.

Definition at line 146 of file TrackTrim.cc.

146 {
147 std::cout << m_className << "::Print:\n";
148 if (m_ions.empty()) {
149 std::cerr << " No TRIM data present.\n";
150 return;
151 }
152 std::cout << " Projectile: " << m_particleName << ", "
153 << m_ekin * 1.e-3 << " keV\n"
154 << " Number of tracks: " << m_ions.size() << "\n"
155 << " Work function: " << m_work << " eV\n"
156 << " Fano factor: " << m_fano << "\n";
157}
double m_ekin
Projectile energy [eV].
Definition: TrackTrim.hh:49
std::string m_particleName
Definition: Track.hh:110

◆ ReadFile()

bool Garfield::TrackTrim::ReadFile ( const std::string &  file,
const unsigned int  nIons = 0,
const unsigned int  nSkip = 0 
)

Load data from an EXYZ.txt file.

Definition at line 30 of file TrackTrim.cc.

31 {
32
33 // TRMREE - Reads the TRIM EXYZ file.
34
35 // Reset.
36 m_ions.clear();
37 m_ion = 0;
38 m_clusters.clear();
39 m_cluster = 0;
40
41 std::ifstream infile;
42 infile.open(filename.c_str(), std::ios::in);
43 if (infile.fail()) {
44 std::cerr << m_className << "::ReadFile:\n"
45 << " Unable to open the EXYZ file (" << filename << ").\n";
46 return false;
47 }
48
49 constexpr double Angstrom = 1.e-8;
50 unsigned int nRead = 0;
51 unsigned int ionNumber = 0;
52 bool header = true;
53 std::vector<float> x;
54 std::vector<float> y;
55 std::vector<float> z;
56 std::vector<float> dedx;
57 std::vector<float> ekin;
58 for (std::string line; std::getline(infile, line);) {
59 if (line.find("------- ") != std::string::npos) {
60 // Reached the end of the header.
61 header = false;
62 continue;
63 } else if (header) {
64 if (line.find("Ion Data: ") != std::string::npos) {
65 // Read the next line.
66 std::getline(infile, line);
67 auto words = tokenize(line);
68 if (words.size() >= 3) {
69 m_particleName = words[0];
70 auto pos = words[2].find("keV");
71 if (pos != std::string::npos) {
72 m_ekin = 1.e3 * std::stod(words[2].substr(0, pos));
73 }
74 }
75 }
76 // Otherwise, skip the header.
77 continue;
78 }
79 auto words = tokenize(line);
80 if (words.size() < 6) {
81 std::cerr << m_className << "::ReadFile: Unexpected line:\n"
82 << line << "\n";
83 continue;
84 }
85 if (ionNumber != std::stoul(words[0])) {
86 // New ion.
87 if (ionNumber > 0) {
88 if (nRead >= nSkip) AddIon(x, y, z, dedx, ekin);
89 x.clear();
90 y.clear();
91 z.clear();
92 dedx.clear();
93 ekin.clear();
94 ++nRead;
95 // Stop if we are done reading the requested number of ions.
96 if (nIons > 0 && m_ions.size() >= nIons) break;
97 }
98 ionNumber = std::stoi(words[0]);
99 }
100 if (nRead < nSkip) continue;
101 // Convert coordinates to cm.
102 x.push_back(std::stof(words[2]) * Angstrom);
103 y.push_back(std::stof(words[3]) * Angstrom);
104 z.push_back(std::stof(words[4]) * Angstrom);
105 // Convert stopping power from eV/A to eV/cm.
106 dedx.push_back(std::stof(words[5]) / Angstrom);
107 // Convert ion energy from keV to eV.
108 ekin.push_back(std::stof(words[1]) * 1.e3);
109 }
110 infile.close();
111 AddIon(x, y, z, dedx, ekin);
112 std::cout << m_className << "::ReadFile: Read energy vs position for "
113 << m_ions.size() << " ions.\n";
114 return true;
115}
void AddIon(const std::vector< float > &x, const std::vector< float > &y, const std::vector< float > &z, const std::vector< float > &dedx, const std::vector< float > &ekin)
Definition: TrackTrim.cc:117

◆ SetFanoFactor()

void Garfield::TrackTrim::SetFanoFactor ( const double  f)
inline

Set the Fano factor.

Definition at line 32 of file TrackTrim.hh.

32{ m_fano = f; }

◆ SetWorkFunction()

void Garfield::TrackTrim::SetWorkFunction ( const double  w)
inline

Set the W value [eV].

Definition at line 28 of file TrackTrim.hh.

28{ m_work = w; }

Member Data Documentation

◆ m_cluster

size_t Garfield::TrackTrim::m_cluster = 0
protected

Index of the next cluster to be returned.

Definition at line 64 of file TrackTrim.hh.

Referenced by GetCluster(), NewTrack(), and ReadFile().

◆ m_clusters

std::vector<Cluster> Garfield::TrackTrim::m_clusters
protected

Clusters on the current track.

Definition at line 62 of file TrackTrim.hh.

Referenced by GetCluster(), NewTrack(), and ReadFile().

◆ m_ekin

double Garfield::TrackTrim::m_ekin = 0.
protected

Projectile energy [eV].

Definition at line 49 of file TrackTrim.hh.

Referenced by Print(), and ReadFile().

◆ m_fano

double Garfield::TrackTrim::m_fano = -1.
protected

Fano factor [-].

Definition at line 46 of file TrackTrim.hh.

Referenced by GetFanoFactor(), NewTrack(), Print(), and SetFanoFactor().

◆ m_ion

size_t Garfield::TrackTrim::m_ion = 0
protected

Index of the current track.

Definition at line 53 of file TrackTrim.hh.

Referenced by NewTrack(), and ReadFile().

◆ m_ions

std::vector<std::vector<std::array<float, 5> > > Garfield::TrackTrim::m_ions
protected

List of tracks.

Definition at line 51 of file TrackTrim.hh.

Referenced by AddIon(), NewTrack(), Print(), and ReadFile().

◆ m_work

double Garfield::TrackTrim::m_work = -1.
protected

Work function [eV].

Definition at line 44 of file TrackTrim.hh.

Referenced by GetWorkFunction(), NewTrack(), Print(), and SetWorkFunction().


The documentation for this class was generated from the following files: