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

#include <ViewDrift.hh>

Public Member Functions

 ViewDrift ()
 
 ~ViewDrift ()
 
void SetCanvas (TCanvas *c)
 
void SetArea (const double &xmin, const double &ymin, const double &zmin, const double &xmax, const double &ymax, const double &zmax)
 
void Clear ()
 
void Plot (const bool twod=false, const bool axis=true)
 
void SetClusterMarkerSize (const double &size)
 
void SetCollisionMarkerSize (const double &size)
 
void NewElectronDriftLine (const unsigned int np, int &id, const double x0, const double y0, const double z0)
 
void NewHoleDriftLine (const unsigned int np, int &id, const double x0, const double y0, const double z0)
 
void NewIonDriftLine (const unsigned int np, int &id, const double x0, const double y0, const double z0)
 
void NewPhotonTrack (const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
 
void NewChargedParticleTrack (const unsigned int np, int &id, const double x0, const double y0, const double z0)
 
void SetDriftLinePoint (const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
 
void AddDriftLinePoint (const unsigned int iL, const double x, const double y, const double z)
 
void SetTrackPoint (const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
 
void AddTrackPoint (const unsigned int iL, const double x, const double y, const double z)
 
void AddExcitationMarker (const double x, const double y, const double z)
 
void AddIonisationMarker (const double x, const double y, const double z)
 
void AddAttachmentMarker (const double x, const double y, const double z)
 
void EnableDebugging ()
 
void DisableDebugging ()
 

Friends

class ViewFEMesh
 

Detailed Description

Definition at line 13 of file ViewDrift.hh.

Constructor & Destructor Documentation

◆ ViewDrift()

Garfield::ViewDrift::ViewDrift ( )

Definition at line 11 of file ViewDrift.cc.

12 : m_className("ViewDrift"),
13 m_debug(false),
14 m_label("Drift Lines"),
15 m_canvas(NULL),
16 m_hasExternalCanvas(false),
17 m_xMin(-1.),
18 m_yMin(-1.),
19 m_zMin(-1.),
20 m_xMax(1.),
21 m_yMax(1.),
22 m_zMax(1.),
23 m_view(NULL),
24 m_nDriftLines(0),
25 m_nTracks(0),
26 m_nExcMarkers(0),
27 m_excPlot(NULL),
28 m_nIonMarkers(0),
29 m_ionPlot(NULL),
30 m_nAttMarkers(0),
31 m_attPlot(NULL),
32 m_markerSizeCluster(1.),
33 m_markerSizeCollision(1.) {
34
35 m_driftLines.reserve(1000);
36 m_driftLinePlots.reserve(1000);
37 m_tracks.reserve(100);
38 m_trackPlots.reserve(100);
39 m_excMarkers.reserve(1000);
40 m_ionMarkers.reserve(1000);
41 m_attMarkers.reserve(1000);
42}

◆ ~ViewDrift()

Garfield::ViewDrift::~ViewDrift ( )

Definition at line 44 of file ViewDrift.cc.

44 {
45
46 if (!m_hasExternalCanvas && m_canvas) delete m_canvas;
47 if (m_view) delete m_view;
48 Clear();
49}

Member Function Documentation

◆ AddAttachmentMarker()

void Garfield::ViewDrift::AddAttachmentMarker ( const double  x,
const double  y,
const double  z 
)

Definition at line 334 of file ViewDrift.cc.

335 {
336
337 marker newMarker;
338 newMarker.x = x;
339 newMarker.y = y;
340 newMarker.z = z;
341 m_attMarkers.push_back(newMarker);
342 ++m_nAttMarkers;
343}

◆ AddDriftLinePoint()

void Garfield::ViewDrift::AddDriftLinePoint ( const unsigned int  iL,
const double  x,
const double  y,
const double  z 
)

Definition at line 269 of file ViewDrift.cc.

270 {
271
272 if (iL >= m_nDriftLines) {
273 std::cerr << m_className << "::AddDriftLinePoint:\n";
274 std::cerr << " Drift line index " << iL << " is out of range.\n";
275 return;
276 }
277 marker m;
278 m.x = x;
279 m.y = y;
280 m.z = z;
281 m_driftLines[iL].vect.push_back(m);
282}

◆ AddExcitationMarker()

void Garfield::ViewDrift::AddExcitationMarker ( const double  x,
const double  y,
const double  z 
)

Definition at line 312 of file ViewDrift.cc.

313 {
314
315 marker newMarker;
316 newMarker.x = x;
317 newMarker.y = y;
318 newMarker.z = z;
319 m_excMarkers.push_back(newMarker);
320 ++m_nExcMarkers;
321}

◆ AddIonisationMarker()

void Garfield::ViewDrift::AddIonisationMarker ( const double  x,
const double  y,
const double  z 
)

Definition at line 323 of file ViewDrift.cc.

324 {
325
326 marker newMarker;
327 newMarker.x = x;
328 newMarker.y = y;
329 newMarker.z = z;
330 m_ionMarkers.push_back(newMarker);
331 ++m_nIonMarkers;
332}

◆ AddTrackPoint()

void Garfield::ViewDrift::AddTrackPoint ( const unsigned int  iL,
const double  x,
const double  y,
const double  z 
)

Definition at line 297 of file ViewDrift.cc.

298 {
299
300 if (iL >= m_nTracks) {
301 std::cerr << m_className << "::AddTrackPoint:\n";
302 std::cerr << " Track index " << iL << " is out of range.\n";
303 return;
304 }
305 marker newPoint;
306 newPoint.x = x;
307 newPoint.y = y;
308 newPoint.z = z;
309 m_tracks[iL].vect.push_back(newPoint);
310}

Referenced by Garfield::Track::PlotCluster().

◆ Clear()

void Garfield::ViewDrift::Clear ( )

Definition at line 81 of file ViewDrift.cc.

81 {
82
83 m_driftLines.clear();
84 m_nDriftLines = 0;
85
86 m_tracks.clear();
87 m_nTracks = 0;
88
89 m_nExcMarkers = m_nIonMarkers = m_nAttMarkers = 0;
90 m_excMarkers.clear();
91 m_ionMarkers.clear();
92 m_attMarkers.clear();
93
94 if (m_excPlot) {
95 delete m_excPlot;
96 m_excPlot = NULL;
97 }
98 if (m_ionPlot) {
99 delete m_ionPlot;
100 m_ionPlot = NULL;
101 }
102 if (m_attPlot) {
103 delete m_attPlot;
104 m_attPlot = NULL;
105 }
106 const unsigned int nTrackPlots = m_trackPlots.size();
107 for (unsigned int i = 0; i < nTrackPlots; ++i) {
108 if (m_trackPlots[i]) delete m_trackPlots[i];
109 }
110 const unsigned int nTrackLinePlots = m_trackLinePlots.size();
111 for (unsigned int i = 0; i < nTrackLinePlots; ++i) {
112 if (m_trackLinePlots[i]) delete m_trackLinePlots[i];
113 }
114 const unsigned int nDriftLinePlots = m_driftLinePlots.size();
115 for (unsigned int i = 0; i < nDriftLinePlots; ++i) {
116 if (m_driftLinePlots[i]) delete m_driftLinePlots[i];
117 }
118}

Referenced by ~ViewDrift().

◆ DisableDebugging()

void Garfield::ViewDrift::DisableDebugging ( )
inline

Definition at line 57 of file ViewDrift.hh.

57{ m_debug = false; }

◆ EnableDebugging()

void Garfield::ViewDrift::EnableDebugging ( )
inline

Definition at line 56 of file ViewDrift.hh.

56{ m_debug = true; }

◆ NewChargedParticleTrack()

void Garfield::ViewDrift::NewChargedParticleTrack ( const unsigned int  np,
int &  id,
const double  x0,
const double  y0,
const double  z0 
)

Definition at line 233 of file ViewDrift.cc.

235 {
236
237
238 // Create a new track and add it to the list.
239 track newTrack;
240 if (np <= 0) {
241 // Number of points is not yet known.
242 newTrack.vect.resize(1);
243 } else {
244 newTrack.vect.resize(np);
245 }
246 newTrack.vect[0].x = x0;
247 newTrack.vect[0].y = y0;
248 newTrack.vect[0].z = z0;
249 m_tracks.push_back(newTrack);
250 // Return the index of this drift line.
251 id = m_nTracks;
252 ++m_nTracks;
253}

Referenced by Garfield::Track::PlotNewTrack().

◆ NewElectronDriftLine()

void Garfield::ViewDrift::NewElectronDriftLine ( const unsigned int  np,
int &  id,
const double  x0,
const double  y0,
const double  z0 
)

Definition at line 140 of file ViewDrift.cc.

142 {
143
144 const int col = plottingEngine.GetRootColorElectron();
145 // Create a new electron drift line and add it to the list.
146 driftLine d;
147 if (np <= 0) {
148 // Number of points is not yet known.
149 std::vector<marker> p(1);
150 p[0].x = x0;
151 p[0].y = y0;
152 p[0].z = z0;
153 d.vect = p;
154 } else {
155 std::vector<marker> p(np);
156 for (unsigned int i = 0; i < np; ++i) {
157 p[i].x = x0;
158 p[i].y = y0;
159 p[i].z = z0;
160 }
161 d.vect = p;
162 }
163 d.n = col;
164 m_driftLines.push_back(d);
165 // Return the index of this drift line.
166 id = m_nDriftLines;
167 ++m_nDriftLines;
168}
PlottingEngineRoot plottingEngine

◆ NewHoleDriftLine()

void Garfield::ViewDrift::NewHoleDriftLine ( const unsigned int  np,
int &  id,
const double  x0,
const double  y0,
const double  z0 
)

Definition at line 170 of file ViewDrift.cc.

172 {
173
174 const int col = plottingEngine.GetRootColorHole();
175 driftLine d;
176 marker m0;
177 m0.x = x0;
178 m0.y = y0;
179 m0.z = z0;
180 if (np <= 0) {
181 // Number of points is not yet known.
182 std::vector<marker> p(1, m0);
183 d.vect = p;
184 } else {
185 std::vector<marker> p(np, m0);
186 d.vect = p;
187 }
188 d.n = col;
189 m_driftLines.push_back(d);
190 // Return the index of this drift line.
191 id = m_nDriftLines;
192 ++m_nDriftLines;
193}

◆ NewIonDriftLine()

void Garfield::ViewDrift::NewIonDriftLine ( const unsigned int  np,
int &  id,
const double  x0,
const double  y0,
const double  z0 
)

Definition at line 195 of file ViewDrift.cc.

196 {
197
198 const int col = 47;
199 driftLine d;
200 marker m0;
201 m0.x = x0;
202 m0.y = y0;
203 m0.z = z0;
204 if (np <= 0) {
205 // Number of points is not yet known.
206 std::vector<marker> p(1, m0);
207 d.vect = p;
208 } else {
209 std::vector<marker> p(np, m0);
210 d.vect = p;
211 }
212 d.n = col;
213 m_driftLines.push_back(d);
214 // Return the index of this drift line.
215 id = m_nDriftLines;
216 ++m_nDriftLines;
217}

◆ NewPhotonTrack()

void Garfield::ViewDrift::NewPhotonTrack ( const double  x0,
const double  y0,
const double  z0,
const double  x1,
const double  y1,
const double  z1 
)

Definition at line 219 of file ViewDrift.cc.

221 {
222
223 const int col = plottingEngine.GetRootColorPhoton();
224 // Create a new photon track (line between start and end point).
225 TPolyLine3D p(2);
226 p.SetLineColor(col);
227 p.SetLineStyle(7);
228 p.SetNextPoint(x0, y0, z0);
229 p.SetNextPoint(x1, y1, z1);
230 ++m_nDriftLines;
231}

◆ Plot()

void Garfield::ViewDrift::Plot ( const bool  twod = false,
const bool  axis = true 
)

Definition at line 345 of file ViewDrift.cc.

345 {
346
347 if (twod) {
348 Plot2d(axis);
349 } else {
350 Plot3d(axis);
351 }
352}

◆ SetArea()

void Garfield::ViewDrift::SetArea ( const double &  xmin,
const double &  ymin,
const double &  zmin,
const double &  xmax,
const double &  ymax,
const double &  zmax 
)

Definition at line 62 of file ViewDrift.cc.

65 {
66
67 // Check range, assign if non-null
68 if (xmin == xmax || ymin == ymax || zmin == zmax) {
69 std::cout << m_className << "::SetArea:\n";
70 std::cout << " Null area range not permitted.\n";
71 return;
72 }
73 m_xMin = std::min(xmin, xmax);
74 m_yMin = std::min(ymin, ymax);
75 m_zMin = std::min(zmin, zmax);
76 m_xMax = std::max(xmin, xmax);
77 m_yMax = std::max(ymin, ymax);
78 m_zMax = std::max(zmin, zmax);
79}

◆ SetCanvas()

void Garfield::ViewDrift::SetCanvas ( TCanvas *  c)

Definition at line 51 of file ViewDrift.cc.

51 {
52
53 if (c == NULL) return;
54 if (!m_hasExternalCanvas && m_canvas) {
55 delete m_canvas;
56 m_canvas = NULL;
57 }
58 m_canvas = c;
59 m_hasExternalCanvas = true;
60}

◆ SetClusterMarkerSize()

void Garfield::ViewDrift::SetClusterMarkerSize ( const double &  size)

Definition at line 120 of file ViewDrift.cc.

120 {
121
122 if (size > 0.) {
123 m_markerSizeCluster = size;
124 } else {
125 std::cerr << m_className << "::SetClusterMarkerSize:\n";
126 std::cerr << " Size must be positive.\n";
127 }
128}

◆ SetCollisionMarkerSize()

void Garfield::ViewDrift::SetCollisionMarkerSize ( const double &  size)

Definition at line 130 of file ViewDrift.cc.

130 {
131
132 if (size > 0.) {
133 m_markerSizeCollision = size;
134 } else {
135 std::cerr << m_className << "::SetCollisionMarkerSize:\n";
136 std::cerr << " Size must be positive.\n";
137 }
138}

◆ SetDriftLinePoint()

void Garfield::ViewDrift::SetDriftLinePoint ( const unsigned int  iL,
const unsigned int  iP,
const double  x,
const double  y,
const double  z 
)

Definition at line 255 of file ViewDrift.cc.

257 {
258
259 if (iL >= m_nDriftLines) {
260 std::cerr << m_className << "::SetDriftLinePoint:\n";
261 std::cerr << " Drift line index " << iL << " is out of range.\n";
262 return;
263 }
264 m_driftLines[iL].vect[iP].x = x;
265 m_driftLines[iL].vect[iP].y = y;
266 m_driftLines[iL].vect[iP].z = z;
267}

◆ SetTrackPoint()

void Garfield::ViewDrift::SetTrackPoint ( const unsigned int  iL,
const unsigned int  iP,
const double  x,
const double  y,
const double  z 
)

Definition at line 284 of file ViewDrift.cc.

285 {
286
287 if (iL >= m_nTracks) {
288 std::cerr << m_className << "::SetTrackPoint:\n";
289 std::cerr << " Track index " << iL << " is out of range.\n";
290 return;
291 }
292 m_tracks[iL].vect[iP].x = x;
293 m_tracks[iL].vect[iP].y = y;
294 m_tracks[iL].vect[iP].z = z;
295}

Friends And Related Function Documentation

◆ ViewFEMesh

friend class ViewFEMesh
friend

Definition at line 59 of file ViewDrift.hh.


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