Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewDrift.cc
Go to the documentation of this file.
1#include <iostream>
2
3#include <TAxis.h>
4#include <TH1F.h>
5
6#include "Plotting.hh"
7#include "ViewDrift.hh"
8
9namespace Garfield {
10
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}
43
45
46 if (!m_hasExternalCanvas && m_canvas) delete m_canvas;
47 if (m_view) delete m_view;
48 Clear();
49}
50
51void ViewDrift::SetCanvas(TCanvas* c) {
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}
61
62void ViewDrift::SetArea(const double& xmin, const double& ymin,
63 const double& zmin,
64 const double& xmax, const double& ymax,
65 const double& zmax) {
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}
80
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}
119
120void ViewDrift::SetClusterMarkerSize(const double& size) {
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}
129
130void ViewDrift::SetCollisionMarkerSize(const double& size) {
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}
139
140void ViewDrift::NewElectronDriftLine(const unsigned int np, int& id,
141 const double x0, const double y0,
142 const double z0) {
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}
169
170void ViewDrift::NewHoleDriftLine(const unsigned int np, int& id,
171 const double x0, const double y0,
172 const double z0) {
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}
194
195void ViewDrift::NewIonDriftLine(const unsigned int np, int& id, const double x0,
196 const double y0, const double z0) {
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}
218
219void ViewDrift::NewPhotonTrack(const double x0, const double y0,
220 const double z0, const double x1,
221 const double y1, const double z1) {
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}
232
233void ViewDrift::NewChargedParticleTrack(const unsigned int np, int& id,
234 const double x0, const double y0,
235 const double z0) {
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}
254
255void ViewDrift::SetDriftLinePoint(const unsigned int iL, const unsigned int iP,
256 const double x, const double y,
257 const double z) {
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}
268
269void ViewDrift::AddDriftLinePoint(const unsigned int iL, const double x,
270 const double y, const double z) {
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}
283
284void ViewDrift::SetTrackPoint(const unsigned int iL, const unsigned int iP,
285 const double x, const double y, const double z) {
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}
296
297void ViewDrift::AddTrackPoint(const unsigned int iL, const double x,
298 const double y, const double z) {
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}
311
312void ViewDrift::AddExcitationMarker(const double x, const double y,
313 const double z) {
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}
322
323void ViewDrift::AddIonisationMarker(const double x, const double y,
324 const double z) {
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}
333
334void ViewDrift::AddAttachmentMarker(const double x, const double y,
335 const double z) {
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}
344
345void ViewDrift::Plot(const bool twod, const bool axis) {
346
347 if (twod) {
348 Plot2d(axis);
349 } else {
350 Plot3d(axis);
351 }
352}
353
354void ViewDrift::Plot2d(const bool axis) {
355
356 std::cout << m_className << "::Plot:\n";
357 std::cout << " Plotting in 2D.\n";
358 if (m_canvas == NULL) {
359 m_canvas = new TCanvas();
360 m_canvas->SetTitle(m_label.c_str());
361 if (m_hasExternalCanvas) m_hasExternalCanvas = false;
362 }
363 m_canvas->cd();
364 m_canvas->Update();
365
366 for (unsigned int i = 0; i < m_nDriftLines; ++i) {
367 const unsigned int nPoints = m_driftLines[i].vect.size();
368 TGraph t(nPoints);
369 for (unsigned int j = 0; j < nPoints; ++j) {
370 t.SetPoint(j, m_driftLines[i].vect[j].x, m_driftLines[i].vect[j].y);
371 }
372 t.SetLineColor(m_driftLines[i].n);
373 t.SetLineWidth(1);
374 if (i == 0) {
375 t.GetXaxis()->SetLimits(m_xMin, m_xMax);
376 t.GetHistogram()->SetMaximum(m_yMax);
377 t.GetHistogram()->SetMinimum(m_yMin);
378 if (axis) {
379 t.DrawClone("ALsame");
380 } else {
381 t.DrawClone("Lsame");
382 }
383 } else {
384 t.DrawClone("Lsame");
385 }
386 m_canvas->Update();
387 }
388
389 const int trackCol = plottingEngine.GetRootColorChargedParticle();
390 for (unsigned int i = 0; i < m_nTracks; ++i) {
391 const unsigned int nPoints = m_tracks[i].vect.size();
392 TGraph t(nPoints);
393 for (unsigned int j = 0; j < nPoints; ++j) {
394 t.SetPoint(j, m_tracks[i].vect[j].x, m_tracks[i].vect[j].y);
395 }
396 t.SetLineColor(trackCol);
397 t.SetLineWidth(2);
398 if (i == 0) {
399 t.GetXaxis()->SetLimits(m_xMin, m_xMax);
400 t.GetHistogram()->SetMaximum(m_yMax);
401 t.GetHistogram()->SetMinimum(m_yMin);
402 if (axis && m_nDriftLines == 0) {
403 t.DrawClone("ALsame");
404 } else {
405 t.DrawClone("Lsame");
406 }
407 } else {
408 t.DrawClone("Lsame");
409 }
410 m_canvas->Update();
411 }
412
413}
414
415void ViewDrift::Plot3d(const bool axis) {
416
417 std::cout << m_className << "::Plot:\n";
418 std::cout << " Plotting in 3D.\n";
419 if (m_canvas == NULL) {
420 m_canvas = new TCanvas();
421 m_canvas->SetTitle(m_label.c_str());
422 if (m_hasExternalCanvas) m_hasExternalCanvas = false;
423 }
424 m_canvas->cd();
425 if (axis) {
426 if (m_canvas->GetView() == NULL) {
427 if (m_view == NULL) m_view = TView::CreateView(1, 0, 0);
428 m_view->SetRange(m_xMin, m_yMin, m_zMin, m_xMax, m_yMax, m_zMax);
429 m_view->ShowAxis();
430 m_view->Top();
431 m_canvas->SetView(m_view);
432 }
433 }
434
435 for (unsigned int i = 0; i < m_nDriftLines; ++i) {
436 const unsigned int nPoints = m_driftLines[i].vect.size();
437 TPolyLine3D* t = new TPolyLine3D(nPoints);
438 for (unsigned int j = 0; j < nPoints; ++j) {
439 t->SetNextPoint(m_driftLines[i].vect[j].x, m_driftLines[i].vect[j].y,
440 m_driftLines[i].vect[j].z);
441 }
442 t->SetLineColor(m_driftLines[i].n);
443 m_driftLinePlots.push_back(t);
444 t->Draw("same");
445 }
446 const int trackCol = plottingEngine.GetRootColorChargedParticle();
447 for (unsigned int i = 0; i < m_nTracks; ++i) {
448 const unsigned int nPoints = m_tracks[i].vect.size();
449 TPolyMarker3D* t = new TPolyMarker3D(nPoints);
450 TPolyLine3D* l = new TPolyLine3D(nPoints);
451 for (unsigned int j = 0; j < nPoints; ++j) {
452 t->SetNextPoint(m_tracks[i].vect[j].x, m_tracks[i].vect[j].y,
453 m_tracks[i].vect[j].z);
454 l->SetNextPoint(m_tracks[i].vect[j].x, m_tracks[i].vect[j].y,
455 m_tracks[i].vect[j].z);
456 }
457 t->SetMarkerStyle(20);
458 t->SetMarkerColor(trackCol);
459 t->SetMarkerSize(m_markerSizeCluster);
460 t->Draw("same");
461 m_trackPlots.push_back(t);
462 l->SetLineColor(trackCol);
463 l->SetLineWidth(1);
464 l->Draw("same");
465 m_trackLinePlots.push_back(l);
466 }
467 if (m_excPlot) {
468 delete m_excPlot;
469 m_excPlot = NULL;
470 }
471 if (m_nExcMarkers > 0) {
472 m_excPlot = new TPolyMarker3D(m_nExcMarkers);
473 m_excPlot->SetMarkerColor(plottingEngine.GetRootColorLine2());
474 m_excPlot->SetMarkerStyle(20);
475 m_excPlot->SetMarkerSize(m_markerSizeCollision);
476 for (unsigned int i = 0; i < m_nExcMarkers; ++i) {
477 m_excPlot->SetNextPoint(m_excMarkers[i].x, m_excMarkers[i].y,
478 m_excMarkers[i].z);
479 }
480 m_excPlot->Draw("same");
481 }
482 if (m_ionPlot) {
483 delete m_ionPlot;
484 m_ionPlot = NULL;
485 }
486 if (m_nIonMarkers > 0) {
487 m_ionPlot = new TPolyMarker3D(m_nIonMarkers);
488 m_ionPlot->SetMarkerColor(plottingEngine.GetRootColorIon());
489 m_ionPlot->SetMarkerStyle(20);
490 m_ionPlot->SetMarkerSize(m_markerSizeCollision);
491 for (unsigned int i = 0; i < m_nIonMarkers; ++i) {
492 m_ionPlot->SetNextPoint(m_ionMarkers[i].x, m_ionMarkers[i].y,
493 m_ionMarkers[i].z);
494 }
495 m_ionPlot->Draw("same");
496 }
497 if (m_attPlot) {
498 delete m_attPlot;
499 m_attPlot = NULL;
500 }
501 if (m_nAttMarkers > 0) {
502 m_attPlot = new TPolyMarker3D(m_nAttMarkers);
503 m_attPlot->SetMarkerColor(plottingEngine.GetRootColorLine1());
504 m_attPlot->SetMarkerStyle(20);
505 m_attPlot->SetMarkerSize(m_markerSizeCollision);
506 for (unsigned int i = 0; i < m_nAttMarkers; ++i) {
507 m_attPlot->SetNextPoint(m_attMarkers[i].x, m_attMarkers[i].y,
508 m_attMarkers[i].z);
509 }
510 m_attPlot->Draw("same");
511 }
512 m_canvas->Update();
513}
514}
void NewElectronDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:140
void AddExcitationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:312
void AddIonisationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:323
void SetArea(const double &xmin, const double &ymin, const double &zmin, const double &xmax, const double &ymax, const double &zmax)
Definition: ViewDrift.cc:62
void Plot(const bool twod=false, const bool axis=true)
Definition: ViewDrift.cc:345
void SetDriftLinePoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:255
void AddAttachmentMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:334
void AddTrackPoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:297
void SetCollisionMarkerSize(const double &size)
Definition: ViewDrift.cc:130
void SetClusterMarkerSize(const double &size)
Definition: ViewDrift.cc:120
void NewChargedParticleTrack(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:233
void AddDriftLinePoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:269
void SetTrackPoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:284
void SetCanvas(TCanvas *c)
Definition: ViewDrift.cc:51
void NewPhotonTrack(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Definition: ViewDrift.cc:219
void NewHoleDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:170
void NewIonDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:195
PlottingEngineRoot plottingEngine