Garfield++ v2r0
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

Visualize drift lines and tracks. More...

#include <ViewDrift.hh>

Public Member Functions

 ViewDrift ()
 Constructor.
 
 ~ViewDrift ()
 Destructor.
 
void SetCanvas (TCanvas *c)
 Set the canvas to be painted on.
 
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.
 
void Clear ()
 Delete existing drift lines, tracks and markers.
 
void Plot (const bool twod=false, const bool axis=true)
 Draw the drift lines.
 
void SetClusterMarkerSize (const double size)
 Set the size of the cluster markers (see TAttMarker).
 
void SetCollisionMarkerSize (const double size)
 Set the size of the collision markers (see TAttMarker).
 
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 (const bool on=true)
 Switch on/off debugging output.
 

Friends

class ViewFEMesh
 

Detailed Description

Visualize drift lines and tracks.

Definition at line 15 of file ViewDrift.hh.

Constructor & Destructor Documentation

◆ ViewDrift()

Garfield::ViewDrift::ViewDrift ( )

Constructor.

Definition at line 11 of file ViewDrift.cc.

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}

◆ ~ViewDrift()

Garfield::ViewDrift::~ViewDrift ( )

Destructor.

Definition at line 38 of file ViewDrift.cc.

38 {
39
40 if (!m_hasExternalCanvas && m_canvas) delete m_canvas;
41 if (m_view) delete m_view;
42 Clear();
43}
void Clear()
Delete existing drift lines, tracks and markers.
Definition: ViewDrift.cc:74

Member Function Documentation

◆ AddAttachmentMarker()

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

Definition at line 303 of file ViewDrift.cc.

304 {
305
306 Marker newMarker;
307 newMarker.x = x;
308 newMarker.y = y;
309 newMarker.z = z;
310 m_attMarkers.push_back(newMarker);
311}

◆ AddDriftLinePoint()

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

Definition at line 243 of file ViewDrift.cc.

244 {
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}

◆ AddExcitationMarker()

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

Definition at line 283 of file ViewDrift.cc.

284 {
285
286 Marker newMarker;
287 newMarker.x = x;
288 newMarker.y = y;
289 newMarker.z = z;
290 m_excMarkers.push_back(newMarker);
291}

◆ AddIonisationMarker()

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

Definition at line 293 of file ViewDrift.cc.

294 {
295
296 Marker newMarker;
297 newMarker.x = x;
298 newMarker.y = y;
299 newMarker.z = z;
300 m_ionMarkers.push_back(newMarker);
301}

◆ AddTrackPoint()

void Garfield::ViewDrift::AddTrackPoint ( 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_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}

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

◆ Clear()

void Garfield::ViewDrift::Clear ( )

Delete existing drift lines, tracks and markers.

Definition at line 74 of file ViewDrift.cc.

74 {
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}

Referenced by ~ViewDrift().

◆ EnableDebugging()

void Garfield::ViewDrift::EnableDebugging ( const bool  on = true)
inline

Switch on/off debugging output.

Definition at line 65 of file ViewDrift.hh.

65{ m_debug = on; }

◆ NewChargedParticleTrack()

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

Definition at line 216 of file ViewDrift.cc.

218 {
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}

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 127 of file ViewDrift.cc.

129 {
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}
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 156 of file ViewDrift.cc.

158 {
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}

◆ NewIonDriftLine()

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

Definition at line 180 of file ViewDrift.cc.

181 {
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}

◆ 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 203 of file ViewDrift.cc.

205 {
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}

◆ Plot()

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

Draw the drift lines.

Definition at line 313 of file ViewDrift.cc.

313 {
314
315 if (twod) {
316 Plot2d(axis);
317 } else {
318 Plot3d(axis);
319 }
320}

◆ SetArea()

void Garfield::ViewDrift::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 at line 56 of file ViewDrift.cc.

59 {
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}

◆ SetCanvas()

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

Set the canvas to be painted on.

Definition at line 45 of file ViewDrift.cc.

45 {
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}

◆ SetClusterMarkerSize()

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

Set the size of the cluster markers (see TAttMarker).

Definition at line 109 of file ViewDrift.cc.

109 {
110
111 if (size > 0.) {
112 m_markerSizeCluster = size;
113 } else {
114 std::cerr << m_className << "::SetClusterMarkerSize: Size must be > 0.\n";
115 }
116}

◆ SetCollisionMarkerSize()

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

Set the size of the collision markers (see TAttMarker).

Definition at line 118 of file ViewDrift.cc.

118 {
119
120 if (size > 0.) {
121 m_markerSizeCollision = size;
122 } else {
123 std::cerr << m_className << "::SetCollisionMarkerSize: Size must be > 0.\n";
124 }
125}

◆ 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 230 of file ViewDrift.cc.

232 {
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}

◆ 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 257 of file ViewDrift.cc.

258 {
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}

Friends And Related Function Documentation

◆ ViewFEMesh

friend class ViewFEMesh
friend

Definition at line 67 of file ViewDrift.hh.


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