Garfield++ v2r0
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_canvas(NULL),
15 m_hasExternalCanvas(false),
16 m_xMin(-1.),
17 m_yMin(-1.),
18 m_zMin(-1.),
19 m_xMax(1.),
20 m_yMax(1.),
21 m_zMax(1.),
22 m_view(NULL),
23 m_excPlot(NULL),
24 m_ionPlot(NULL),
25 m_attPlot(NULL),
26 m_markerSizeCluster(1.),
27 m_markerSizeCollision(1.) {
28
29 m_driftLines.reserve(1000);
30 m_driftLinePlots.reserve(1000);
31 m_tracks.reserve(100);
32 m_trackPlots.reserve(100);
33 m_excMarkers.reserve(1000);
34 m_ionMarkers.reserve(1000);
35 m_attMarkers.reserve(1000);
36}
37
39
40 if (!m_hasExternalCanvas && m_canvas) delete m_canvas;
41 if (m_view) delete m_view;
42 Clear();
43}
44
45void ViewDrift::SetCanvas(TCanvas* c) {
46
47 if (!c) return;
48 if (!m_hasExternalCanvas && m_canvas) {
49 delete m_canvas;
50 m_canvas = NULL;
51 }
52 m_canvas = c;
53 m_hasExternalCanvas = true;
54}
55
56void ViewDrift::SetArea(const double xmin, const double ymin,
57 const double zmin,
58 const double xmax, const double ymax,
59 const double zmax) {
60
61 // Check range, assign if non-null
62 if (xmin == xmax || ymin == ymax || zmin == zmax) {
63 std::cerr << m_className << "::SetArea: Null area range not permitted.\n";
64 return;
65 }
66 m_xMin = std::min(xmin, xmax);
67 m_yMin = std::min(ymin, ymax);
68 m_zMin = std::min(zmin, zmax);
69 m_xMax = std::max(xmin, xmax);
70 m_yMax = std::max(ymin, ymax);
71 m_zMax = std::max(zmin, zmax);
72}
73
75
76 m_driftLines.clear();
77 m_tracks.clear();
78
79 m_excMarkers.clear();
80 m_ionMarkers.clear();
81 m_attMarkers.clear();
82
83 if (m_excPlot) {
84 delete m_excPlot;
85 m_excPlot = NULL;
86 }
87 if (m_ionPlot) {
88 delete m_ionPlot;
89 m_ionPlot = NULL;
90 }
91 if (m_attPlot) {
92 delete m_attPlot;
93 m_attPlot = NULL;
94 }
95 const unsigned int nTrackPlots = m_trackPlots.size();
96 for (unsigned int i = 0; i < nTrackPlots; ++i) {
97 if (m_trackPlots[i]) delete m_trackPlots[i];
98 }
99 const unsigned int nTrackLinePlots = m_trackLinePlots.size();
100 for (unsigned int i = 0; i < nTrackLinePlots; ++i) {
101 if (m_trackLinePlots[i]) delete m_trackLinePlots[i];
102 }
103 const unsigned int nDriftLinePlots = m_driftLinePlots.size();
104 for (unsigned int i = 0; i < nDriftLinePlots; ++i) {
105 if (m_driftLinePlots[i]) delete m_driftLinePlots[i];
106 }
107}
108
109void ViewDrift::SetClusterMarkerSize(const double size) {
110
111 if (size > 0.) {
112 m_markerSizeCluster = size;
113 } else {
114 std::cerr << m_className << "::SetClusterMarkerSize: Size must be > 0.\n";
115 }
116}
117
118void ViewDrift::SetCollisionMarkerSize(const double size) {
119
120 if (size > 0.) {
121 m_markerSizeCollision = size;
122 } else {
123 std::cerr << m_className << "::SetCollisionMarkerSize: Size must be > 0.\n";
124 }
125}
126
127void ViewDrift::NewElectronDriftLine(const unsigned int np, int& id,
128 const double x0, const double y0,
129 const double z0) {
130
131 const int col = plottingEngine.GetRootColorElectron();
132 // Create a new electron drift line and add it to the list.
133 driftLine d;
134 if (np <= 0) {
135 // Number of points is not yet known.
136 std::vector<Marker> p(1);
137 p[0].x = x0;
138 p[0].y = y0;
139 p[0].z = z0;
140 d.vect = p;
141 } else {
142 std::vector<Marker> p(np);
143 for (unsigned int i = 0; i < np; ++i) {
144 p[i].x = x0;
145 p[i].y = y0;
146 p[i].z = z0;
147 }
148 d.vect = p;
149 }
150 d.n = col;
151 m_driftLines.push_back(d);
152 // Return the index of this drift line.
153 id = m_driftLines.size() - 1;
154}
155
156void ViewDrift::NewHoleDriftLine(const unsigned int np, int& id,
157 const double x0, const double y0,
158 const double z0) {
159
160 const int col = plottingEngine.GetRootColorHole();
161 driftLine d;
162 Marker m0;
163 m0.x = x0;
164 m0.y = y0;
165 m0.z = z0;
166 if (np <= 0) {
167 // Number of points is not yet known.
168 std::vector<Marker> p(1, m0);
169 d.vect = p;
170 } else {
171 std::vector<Marker> p(np, m0);
172 d.vect = p;
173 }
174 d.n = col;
175 m_driftLines.push_back(d);
176 // Return the index of this drift line.
177 id = m_driftLines.size() - 1;
178}
179
180void ViewDrift::NewIonDriftLine(const unsigned int np, int& id, const double x0,
181 const double y0, const double z0) {
182
183 const int col = 47;
184 driftLine d;
185 Marker m0;
186 m0.x = x0;
187 m0.y = y0;
188 m0.z = z0;
189 if (np <= 0) {
190 // Number of points is not yet known.
191 std::vector<Marker> p(1, m0);
192 d.vect = p;
193 } else {
194 std::vector<Marker> p(np, m0);
195 d.vect = p;
196 }
197 d.n = col;
198 m_driftLines.push_back(d);
199 // Return the index of this drift line.
200 id = m_driftLines.size() - 1;
201}
202
203void ViewDrift::NewPhotonTrack(const double x0, const double y0,
204 const double z0, const double x1,
205 const double y1, const double z1) {
206
207 const int col = plottingEngine.GetRootColorPhoton();
208 // Create a new photon track (line between start and end point).
209 TPolyLine3D p(2);
210 p.SetLineColor(col);
211 p.SetLineStyle(7);
212 p.SetNextPoint(x0, y0, z0);
213 p.SetNextPoint(x1, y1, z1);
214}
215
216void ViewDrift::NewChargedParticleTrack(const unsigned int np, int& id,
217 const double x0, const double y0,
218 const double z0) {
219
220 // Create a new track and add it to the list.
221 std::vector<Marker> track(std::max(1U, np));
222 track[0].x = x0;
223 track[0].y = y0;
224 track[0].z = z0;
225 m_tracks.push_back(track);
226 // Return the index of this track.
227 id = m_tracks.size() - 1;
228}
229
230void ViewDrift::SetDriftLinePoint(const unsigned int iL, const unsigned int iP,
231 const double x, const double y,
232 const double z) {
233
234 if (iL >= m_driftLines.size()) {
235 std::cerr << m_className << "::SetDriftLinePoint: Index out of range.\n";
236 return;
237 }
238 m_driftLines[iL].vect[iP].x = x;
239 m_driftLines[iL].vect[iP].y = y;
240 m_driftLines[iL].vect[iP].z = z;
241}
242
243void ViewDrift::AddDriftLinePoint(const unsigned int iL, const double x,
244 const double y, const double z) {
245
246 if (iL >= m_driftLines.size()) {
247 std::cerr << m_className << "::AddDriftLinePoint: Index out of range.\n";
248 return;
249 }
250 Marker m;
251 m.x = x;
252 m.y = y;
253 m.z = z;
254 m_driftLines[iL].vect.push_back(m);
255}
256
257void ViewDrift::SetTrackPoint(const unsigned int iL, const unsigned int iP,
258 const double x, const double y, const double z) {
259
260 if (iL >= m_tracks.size()) {
261 std::cerr << m_className << "::SetTrackPoint: Index out of range.\n";
262 return;
263 }
264 m_tracks[iL][iP].x = x;
265 m_tracks[iL][iP].y = y;
266 m_tracks[iL][iP].z = z;
267}
268
269void ViewDrift::AddTrackPoint(const unsigned int iL, const double x,
270 const double y, const double z) {
271
272 if (iL >= m_tracks.size()) {
273 std::cerr << m_className << "::AddTrackPoint: Index out of range.\n";
274 return;
275 }
276 Marker newPoint;
277 newPoint.x = x;
278 newPoint.y = y;
279 newPoint.z = z;
280 m_tracks[iL].push_back(newPoint);
281}
282
283void ViewDrift::AddExcitationMarker(const double x, const double y,
284 const double z) {
285
286 Marker newMarker;
287 newMarker.x = x;
288 newMarker.y = y;
289 newMarker.z = z;
290 m_excMarkers.push_back(newMarker);
291}
292
293void ViewDrift::AddIonisationMarker(const double x, const double y,
294 const double z) {
295
296 Marker newMarker;
297 newMarker.x = x;
298 newMarker.y = y;
299 newMarker.z = z;
300 m_ionMarkers.push_back(newMarker);
301}
302
303void ViewDrift::AddAttachmentMarker(const double x, const double y,
304 const double z) {
305
306 Marker newMarker;
307 newMarker.x = x;
308 newMarker.y = y;
309 newMarker.z = z;
310 m_attMarkers.push_back(newMarker);
311}
312
313void ViewDrift::Plot(const bool twod, const bool axis) {
314
315 if (twod) {
316 Plot2d(axis);
317 } else {
318 Plot3d(axis);
319 }
320}
321
322void ViewDrift::Plot2d(const bool axis) {
323
324 std::cout << m_className << "::Plot: Plotting in 2D.\n";
325 if (!m_canvas) {
326 m_canvas = new TCanvas();
327 if (m_hasExternalCanvas) m_hasExternalCanvas = false;
328 }
329 m_canvas->cd();
330 m_canvas->Update();
331
332 const unsigned int nDriftLines = m_driftLines.size();
333 for (unsigned int i = 0; i < nDriftLines; ++i) {
334 const unsigned int nPoints = m_driftLines[i].vect.size();
335 TGraph t(nPoints);
336 for (unsigned int j = 0; j < nPoints; ++j) {
337 t.SetPoint(j, m_driftLines[i].vect[j].x, m_driftLines[i].vect[j].y);
338 }
339 t.SetLineColor(m_driftLines[i].n);
340 t.SetLineWidth(1);
341 if (i == 0) {
342 t.GetXaxis()->SetLimits(m_xMin, m_xMax);
343 t.GetHistogram()->SetMaximum(m_yMax);
344 t.GetHistogram()->SetMinimum(m_yMin);
345 if (axis) {
346 t.DrawClone("ALsame");
347 } else {
348 t.DrawClone("Lsame");
349 }
350 } else {
351 t.DrawClone("Lsame");
352 }
353 m_canvas->Update();
354 }
355
356 const int trackCol = plottingEngine.GetRootColorChargedParticle();
357 const unsigned int nTracks = m_tracks.size();
358 for (unsigned int i = 0; i < nTracks; ++i) {
359 const std::vector<Marker>& track = m_tracks[i];
360 const unsigned int nPoints = track.size();
361 TGraph t(nPoints);
362 for (unsigned int j = 0; j < nPoints; ++j) {
363 t.SetPoint(j, track[j].x, track[j].y);
364 }
365 t.SetLineColor(trackCol);
366 t.SetLineWidth(2);
367 if (i == 0) {
368 t.GetXaxis()->SetLimits(m_xMin, m_xMax);
369 t.GetHistogram()->SetMaximum(m_yMax);
370 t.GetHistogram()->SetMinimum(m_yMin);
371 if (axis && m_driftLines.empty()) {
372 t.DrawClone("ALsame");
373 } else {
374 t.DrawClone("Lsame");
375 }
376 } else {
377 t.DrawClone("Lsame");
378 }
379 m_canvas->Update();
380 }
381
382}
383
384void ViewDrift::Plot3d(const bool axis) {
385
386 std::cout << m_className << "::Plot: Plotting in 3D.\n";
387 if (!m_canvas) {
388 m_canvas = new TCanvas();
389 if (m_hasExternalCanvas) m_hasExternalCanvas = false;
390 }
391 m_canvas->cd();
392 if (axis) {
393 if (!m_canvas->GetView()) {
394 if (!m_view) m_view = TView::CreateView(1, 0, 0);
395 m_view->SetRange(m_xMin, m_yMin, m_zMin, m_xMax, m_yMax, m_zMax);
396 m_view->ShowAxis();
397 m_view->Top();
398 m_canvas->SetView(m_view);
399 }
400 }
401
402 const unsigned int nDriftLines = m_driftLines.size();
403 for (unsigned int i = 0; i < nDriftLines; ++i) {
404 const unsigned int nPoints = m_driftLines[i].vect.size();
405 TPolyLine3D* t = new TPolyLine3D(nPoints);
406 for (unsigned int j = 0; j < nPoints; ++j) {
407 t->SetNextPoint(m_driftLines[i].vect[j].x, m_driftLines[i].vect[j].y,
408 m_driftLines[i].vect[j].z);
409 }
410 t->SetLineColor(m_driftLines[i].n);
411 m_driftLinePlots.push_back(t);
412 t->Draw("same");
413 }
414 const int trackCol = plottingEngine.GetRootColorChargedParticle();
415 const unsigned int nTracks = m_tracks.size();
416 for (unsigned int i = 0; i < nTracks; ++i) {
417 const std::vector<Marker>& track = m_tracks[i];
418 const unsigned int nPoints = track.size();
419 TPolyMarker3D* t = new TPolyMarker3D(nPoints);
420 TPolyLine3D* l = new TPolyLine3D(nPoints);
421 for (unsigned int j = 0; j < nPoints; ++j) {
422 t->SetNextPoint(track[j].x, track[j].y, track[j].z);
423 l->SetNextPoint(track[j].x, track[j].y, track[j].z);
424 }
425 t->SetMarkerStyle(20);
426 t->SetMarkerColor(trackCol);
427 t->SetMarkerSize(m_markerSizeCluster);
428 t->Draw("same");
429 m_trackPlots.push_back(t);
430 l->SetLineColor(trackCol);
431 l->SetLineWidth(1);
432 l->Draw("same");
433 m_trackLinePlots.push_back(l);
434 }
435 if (m_excPlot) {
436 delete m_excPlot;
437 m_excPlot = NULL;
438 }
439 if (!m_excMarkers.empty()) {
440 const unsigned int nExcMarkers = m_excMarkers.size();
441 m_excPlot = new TPolyMarker3D(nExcMarkers);
442 m_excPlot->SetMarkerColor(plottingEngine.GetRootColorLine2());
443 m_excPlot->SetMarkerStyle(20);
444 m_excPlot->SetMarkerSize(m_markerSizeCollision);
445 for (unsigned int i = 0; i < nExcMarkers; ++i) {
446 m_excPlot->SetNextPoint(m_excMarkers[i].x, m_excMarkers[i].y,
447 m_excMarkers[i].z);
448 }
449 m_excPlot->Draw("same");
450 }
451 if (m_ionPlot) {
452 delete m_ionPlot;
453 m_ionPlot = NULL;
454 }
455 if (!m_ionMarkers.empty()) {
456 const unsigned int nIonMarkers = m_ionMarkers.size();
457 m_ionPlot = new TPolyMarker3D(nIonMarkers);
458 m_ionPlot->SetMarkerColor(plottingEngine.GetRootColorIon());
459 m_ionPlot->SetMarkerStyle(20);
460 m_ionPlot->SetMarkerSize(m_markerSizeCollision);
461 for (unsigned int i = 0; i < nIonMarkers; ++i) {
462 m_ionPlot->SetNextPoint(m_ionMarkers[i].x, m_ionMarkers[i].y,
463 m_ionMarkers[i].z);
464 }
465 m_ionPlot->Draw("same");
466 }
467 if (m_attPlot) {
468 delete m_attPlot;
469 m_attPlot = NULL;
470 }
471 if (!m_attMarkers.empty()) {
472 const unsigned int nAttMarkers = m_attMarkers.size();
473 m_attPlot = new TPolyMarker3D(nAttMarkers);
474 m_attPlot->SetMarkerColor(plottingEngine.GetRootColorLine1());
475 m_attPlot->SetMarkerStyle(20);
476 m_attPlot->SetMarkerSize(m_markerSizeCollision);
477 for (unsigned int i = 0; i < nAttMarkers; ++i) {
478 m_attPlot->SetNextPoint(m_attMarkers[i].x, m_attMarkers[i].y,
479 m_attMarkers[i].z);
480 }
481 m_attPlot->Draw("same");
482 }
483 m_canvas->Update();
484}
485}
void NewElectronDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:127
void AddExcitationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:283
void AddIonisationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:293
void Plot(const bool twod=false, const bool axis=true)
Draw the drift lines.
Definition: ViewDrift.cc:313
void SetDriftLinePoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:230
void SetArea(const double xmin, const double ymin, const double zmin, const double xmax, const double ymax, const double zmax)
Set the region to be plotted.
Definition: ViewDrift.cc:56
void AddAttachmentMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:303
void AddTrackPoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:269
void SetClusterMarkerSize(const double size)
Set the size of the cluster markers (see TAttMarker).
Definition: ViewDrift.cc:109
void NewChargedParticleTrack(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:216
void AddDriftLinePoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:243
ViewDrift()
Constructor.
Definition: ViewDrift.cc:11
void SetTrackPoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:257
void Clear()
Delete existing drift lines, tracks and markers.
Definition: ViewDrift.cc:74
void SetCanvas(TCanvas *c)
Set the canvas to be painted on.
Definition: ViewDrift.cc:45
void NewPhotonTrack(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Definition: ViewDrift.cc:203
void NewHoleDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:156
~ViewDrift()
Destructor.
Definition: ViewDrift.cc:38
void NewIonDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:180
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
Definition: ViewDrift.cc:118
PlottingEngineRoot plottingEngine