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