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

+ Inheritance diagram for Garfield::ViewDrift:

Public Member Functions

 ViewDrift ()
 Constructor.
 
 ~ViewDrift ()=default
 Destructor.
 
void Clear ()
 Delete existing drift lines, tracks and markers.
 
void Plot (const bool twod=false, const bool axis=true)
 Draw the drift lines.
 
void Plot2d (const bool axis)
 Make a 2D plot of the drift lines in the current viewing plane.
 
void Plot3d (const bool axis, const bool ogl)
 Make a 3D plot of the drift lines.
 
void EnableClusterMarkers (const bool on=true)
 Draw markers (or not) at every collision along a track.
 
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 SetColourElectrons (const short col)
 Set the colour with which to draw electron drift lines.
 
void SetColourHoles (const short col)
 Set the colour with which to draw hole drift lines.
 
void SetColourIons (const short col)
 Set the colour with which to draw ion drift lines.
 
void SetColourTracks (const short col)
 Set the colour with which to draw charged particle tracks.
 
void SetColourPhotons (const short col)
 Set the colour with which to draw photons.
 
void SetColourExcitations (const short col)
 Set the colour with which to draw excitation markers.
 
void SetColourIonisations (const short col)
 Set the colour with which to draw ionisation markers.
 
void SetColourAttachments (const short col)
 Set the colour with which to draw attachment markers.
 
size_t GetNumberOfDriftLines () const
 Get the number of drift lines stored.
 
void GetDriftLine (const size_t i, std::vector< std::array< float, 3 > > &driftLine, bool &electron) const
 Retrieve the coordinates of a given drift line.
 
void NewElectronDriftLine (const size_t np, size_t &id, const float x0, const float y0, const float z0)
 
void NewHoleDriftLine (const size_t np, size_t &id, const float x0, const float y0, const float z0)
 
void NewIonDriftLine (const size_t np, size_t &id, const float x0, const float y0, const float z0)
 
void NewChargedParticleTrack (const size_t np, size_t &id, const float x0, const float y0, const float z0)
 
void SetDriftLinePoint (const size_t iL, const size_t iP, const float x, const float y, const float z)
 
void AddDriftLinePoint (const size_t iL, const float x, const float y, const float z)
 
void SetTrackPoint (const size_t iL, const size_t iP, const float x, const float y, const float z)
 
void AddTrackPoint (const size_t iL, const float x, const float y, const float z)
 
void AddExcitation (const float x, const float y, const float z)
 
void AddIonisation (const float x, const float y, const float z)
 
void AddAttachment (const float x, const float y, const float z)
 
void AddPhoton (const float x0, const float y0, const float z0, const float x1, const float y1, const float z1)
 
- Public Member Functions inherited from Garfield::ViewBase
 ViewBase ()=delete
 Default constructor.
 
 ViewBase (const std::string &name)
 Constructor.
 
virtual ~ViewBase ()=default
 Destructor.
 
void SetCanvas (TPad *pad)
 Set the canvas to be painted on.
 
void SetCanvas ()
 Unset an external canvas.
 
TPad * GetCanvas ()
 Retrieve the canvas.
 
void SetArea (const double xmin, const double ymin, const double xmax, const double ymax)
 
virtual void SetArea (const double xmin, const double ymin, const double zmin, const double xmax, const double ymax, const double zmax)
 Set a bounding box (if applicable).
 
void SetArea ()
 
virtual void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
 
virtual void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0, const double hx, const double hy, const double hz)
 Set the projection plane specifying a hint for the in-plane x axis.
 
void Rotate (const double angle)
 Rotate the viewing plane (angle in radian).
 
void SetPlaneXY ()
 Set the viewing plane to x-y.
 
void SetPlaneXZ ()
 Set the viewing plane to x-z.
 
void SetPlaneYZ ()
 Set the viewing plane to y-z.
 
void EnableDebugging (const bool on=true)
 Switch on/off debugging output.
 

Friends

class ViewFEMesh
 

Additional Inherited Members

- Static Public Member Functions inherited from Garfield::ViewBase
static std::string FindUnusedFunctionName (const std::string &s)
 Find an unused function name.
 
static std::string FindUnusedHistogramName (const std::string &s)
 Find an unused histogram name.
 
static std::string FindUnusedCanvasName (const std::string &s)
 Find an unused canvas name.
 
- Protected Member Functions inherited from Garfield::ViewBase
void UpdateProjectionMatrix ()
 
template<typename T >
void ToPlane (const T x, const T y, const T z, T &xp, T &yp) const
 
template<typename T >
bool InBox (const std::array< T, 3 > &x) const
 
void Clip (const std::array< float, 3 > &x0, const std::array< float, 3 > &x1, std::array< float, 3 > &xc) const
 
void DrawLine (const std::vector< std::array< float, 3 > > &xl, const short col, const short lw)
 
std::string LabelX ()
 
std::string LabelY ()
 
std::string PlaneDescription ()
 
bool PlotLimits (Sensor *sensor, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimits (Component *cmp, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimitsFromUserBox (double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimits (std::array< double, 3 > &bbmin, std::array< double, 3 > &bbmax, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool RangeSet (TPad *)
 
void SetRange (TPad *pad, const double x0, const double y0, const double x1, const double y1)
 
- Protected Attributes inherited from Garfield::ViewBase
std::string m_className = "ViewBase"
 
bool m_debug = false
 
bool m_userPlotLimits = false
 
double m_xMinPlot = -1.
 
double m_xMaxPlot = 1.
 
double m_yMinPlot = -1.
 
double m_yMaxPlot = 1.
 
bool m_userBox = false
 
double m_xMinBox = -1.
 
double m_xMaxBox = 1.
 
double m_yMinBox = -1.
 
double m_yMaxBox = 1.
 
double m_zMinBox = -1.
 
double m_zMaxBox = 1.
 
std::array< std::array< double, 3 >, 3 > m_proj
 
std::array< double, 4 > m_plane {{0, 0, 1, 0}}
 
std::array< std::array< double, 3 >, 3 > m_prmat
 

Detailed Description

Visualize drift lines and tracks.

Definition at line 18 of file ViewDrift.hh.

Constructor & Destructor Documentation

◆ ViewDrift()

Garfield::ViewDrift::ViewDrift ( )

Constructor.

Definition at line 21 of file ViewDrift.cc.

21 : ViewBase("ViewDrift") {
22 m_driftLines.reserve(1000);
23 m_tracks.reserve(100);
24 m_exc.reserve(1000);
25 m_ion.reserve(1000);
26 m_att.reserve(1000);
27}
ViewBase()=delete
Default constructor.

◆ ~ViewDrift()

Garfield::ViewDrift::~ViewDrift ( )
default

Destructor.

Member Function Documentation

◆ AddAttachment()

void Garfield::ViewDrift::AddAttachment ( const float  x,
const float  y,
const float  z 
)

Definition at line 179 of file ViewDrift.cc.

179 {
180 std::lock_guard<std::mutex> guard(m_mutex);
181 std::array<float, 3> p = {x, y, z};
182 m_att.push_back(std::move(p));
183}

◆ AddDriftLinePoint()

void Garfield::ViewDrift::AddDriftLinePoint ( const size_t  iL,
const float  x,
const float  y,
const float  z 
)

Definition at line 135 of file ViewDrift.cc.

136 {
137 std::lock_guard<std::mutex> guard(m_mutex);
138 if (iL >= m_driftLines.size()) {
139 std::cerr << m_className << "::AddDriftLinePoint: Index out of range.\n";
140 return;
141 }
142 std::array<float, 3> p = {x, y, z};
143 m_driftLines[iL].first.push_back(std::move(p));
144}
std::string m_className
Definition: ViewBase.hh:78

◆ AddExcitation()

void Garfield::ViewDrift::AddExcitation ( const float  x,
const float  y,
const float  z 
)

Definition at line 167 of file ViewDrift.cc.

167 {
168 std::lock_guard<std::mutex> guard(m_mutex);
169 std::array<float, 3> p = {x, y, z};
170 m_exc.push_back(std::move(p));
171}

◆ AddIonisation()

void Garfield::ViewDrift::AddIonisation ( const float  x,
const float  y,
const float  z 
)

Definition at line 173 of file ViewDrift.cc.

173 {
174 std::lock_guard<std::mutex> guard(m_mutex);
175 std::array<float, 3> p = {x, y, z};
176 m_ion.push_back(std::move(p));
177}

◆ AddPhoton()

void Garfield::ViewDrift::AddPhoton ( const float  x0,
const float  y0,
const float  z0,
const float  x1,
const float  y1,
const float  z1 
)

Definition at line 103 of file ViewDrift.cc.

104 {
105 std::lock_guard<std::mutex> guard(m_mutex);
106 std::array<float, 3> p0 = {x0, y0, z0};
107 std::array<float, 3> p1 = {x1, y1, z1};
108 m_photons.push_back({p0, p1});
109}

◆ AddTrackPoint()

void Garfield::ViewDrift::AddTrackPoint ( const size_t  iL,
const float  x,
const float  y,
const float  z 
)

Definition at line 156 of file ViewDrift.cc.

157 {
158 std::lock_guard<std::mutex> guard(m_mutex);
159 if (iL >= m_tracks.size()) {
160 std::cerr << m_className << "::AddTrackPoint: Index out of range.\n";
161 return;
162 }
163 std::array<float, 3> p = {x, y, z};
164 m_tracks[iL].push_back(std::move(p));
165}

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

◆ Clear()

void Garfield::ViewDrift::Clear ( )

Delete existing drift lines, tracks and markers.

Definition at line 29 of file ViewDrift.cc.

29 {
30 m_driftLines.clear();
31 m_tracks.clear();
32
33 m_exc.clear();
34 m_ion.clear();
35 m_att.clear();
36}

◆ EnableClusterMarkers()

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

Draw markers (or not) at every collision along a track.

Definition at line 36 of file ViewDrift.hh.

36{ m_drawClusters = on; }

◆ GetDriftLine()

void Garfield::ViewDrift::GetDriftLine ( const size_t  i,
std::vector< std::array< float, 3 > > &  driftLine,
bool &  electron 
) const

Retrieve the coordinates of a given drift line.

Definition at line 54 of file ViewDrift.cc.

56 {
57 driftLine.clear();
58 if (i >= m_driftLines.size()) return;
59 std::copy(m_driftLines[i].first.begin(), m_driftLines[i].first.end(),
60 std::back_inserter(driftLine));
61 if (m_driftLines[i].second == Particle::Electron) {
62 electron = true;
63 } else {
64 electron = false;
65 }
66}

◆ GetNumberOfDriftLines()

size_t Garfield::ViewDrift::GetNumberOfDriftLines ( ) const
inline

Get the number of drift lines stored.

Definition at line 60 of file ViewDrift.hh.

60{ return m_driftLines.size(); }

◆ NewChargedParticleTrack()

void Garfield::ViewDrift::NewChargedParticleTrack ( const size_t  np,
size_t &  id,
const float  x0,
const float  y0,
const float  z0 
)

Definition at line 111 of file ViewDrift.cc.

113 {
114 std::lock_guard<std::mutex> guard(m_mutex);
115 // Create a new track and add it to the list.
116 constexpr size_t minSize = 1;
117 std::vector<std::array<float, 3> > track(std::max(minSize, np));
118 track[0] = {x0, y0, z0};
119 m_tracks.push_back(std::move(track));
120 // Return the index of this track.
121 id = m_tracks.size() - 1;
122}

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

◆ NewElectronDriftLine()

void Garfield::ViewDrift::NewElectronDriftLine ( const size_t  np,
size_t &  id,
const float  x0,
const float  y0,
const float  z0 
)

Definition at line 68 of file ViewDrift.cc.

70 {
71 std::lock_guard<std::mutex> guard(m_mutex);
72 // Create a new electron drift line and add it to the list.
73 std::array<float, 3> p = {x0, y0, z0};
74 constexpr size_t minSize = 1;
75 std::vector<std::array<float, 3> > dl(std::max(minSize, np), p);
76 m_driftLines.push_back(std::make_pair(std::move(dl), Particle::Electron));
77 // Return the index of this drift line.
78 id = m_driftLines.size() - 1;
79}

◆ NewHoleDriftLine()

void Garfield::ViewDrift::NewHoleDriftLine ( const size_t  np,
size_t &  id,
const float  x0,
const float  y0,
const float  z0 
)

Definition at line 81 of file ViewDrift.cc.

82 {
83 std::lock_guard<std::mutex> guard(m_mutex);
84 std::array<float, 3> p = {x0, y0, z0};
85 constexpr size_t minSize = 1;
86 std::vector<std::array<float, 3> > dl(std::max(minSize, np), p);
87 m_driftLines.push_back(std::make_pair(std::move(dl), Particle::Hole));
88 // Return the index of this drift line.
89 id = m_driftLines.size() - 1;
90}

◆ NewIonDriftLine()

void Garfield::ViewDrift::NewIonDriftLine ( const size_t  np,
size_t &  id,
const float  x0,
const float  y0,
const float  z0 
)

Definition at line 92 of file ViewDrift.cc.

93 {
94 std::lock_guard<std::mutex> guard(m_mutex);
95 std::array<float, 3> p = {x0, y0, z0};
96 constexpr size_t minSize = 1;
97 std::vector<std::array<float, 3> > dl(std::max(minSize, np), p);
98 m_driftLines.push_back(std::make_pair(std::move(dl), Particle::Ion));
99 // Return the index of this drift line.
100 id = m_driftLines.size() - 1;
101}

◆ Plot()

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

Draw the drift lines.

Definition at line 185 of file ViewDrift.cc.

185 {
186 if (twod) {
187 Plot2d(axis);
188 } else {
189 Plot3d(axis, false);
190 }
191}
void Plot3d(const bool axis, const bool ogl)
Make a 3D plot of the drift lines.
Definition: ViewDrift.cc:309
void Plot2d(const bool axis)
Make a 2D plot of the drift lines in the current viewing plane.
Definition: ViewDrift.cc:193

Referenced by main().

◆ Plot2d()

void Garfield::ViewDrift::Plot2d ( const bool  axis)

Make a 2D plot of the drift lines in the current viewing plane.

Definition at line 193 of file ViewDrift.cc.

193 {
194 auto pad = GetCanvas();
195 pad->cd();
196 pad->SetTitle("Drift lines");
197 // Check if the canvas range has already been set.
198 const bool rangeSet = RangeSet(pad);
199 if (axis || !rangeSet) {
200 // Determine the plot limits.
201 if (!SetPlotLimits2d()) {
202 std::cerr << m_className << "::Plot2d:\n"
203 << " Could not determine the plot limits.\n";
204 return;
205 }
206 }
207 if (axis) {
208 auto frame = pad->DrawFrame(m_xMinPlot, m_yMinPlot,
210 frame->GetXaxis()->SetTitle(LabelX().c_str());
211 frame->GetYaxis()->SetTitle(LabelY().c_str());
212 } else if (!rangeSet) {
214 }
215
216 for (const auto& driftLine : m_driftLines) {
217 const short lw = 1;
218 if (driftLine.second == Particle::Electron) {
219 DrawLine(driftLine.first, m_colElectron, lw);
220 } else if (driftLine.second == Particle::Hole) {
221 DrawLine(driftLine.first, m_colHole, lw);
222 } else {
223 DrawLine(driftLine.first, m_colIon, lw);
224 }
225 }
226 gPad->Update();
227
228 TGraph gr;
229 gr.SetMarkerColor(m_colTrack);
230 gr.SetMarkerSize(m_markerSizeCluster);
231 for (const auto& track : m_tracks) {
232 DrawLine(track, m_colTrack, 2);
233 if (!m_drawClusters) continue;
234 std::vector<float> xgr;
235 std::vector<float> ygr;
236 for (const auto& p : track) {
237 if (!InBox(p)) continue;
238 float xp = 0., yp = 0.;
239 ToPlane(p[0], p[1], p[2], xp, yp);
240 xgr.push_back(xp);
241 ygr.push_back(yp);
242 }
243 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "Psame");
244 }
245
246 gr.SetLineColor(m_colPhoton);
247 gr.SetLineStyle(2);
248 for (const auto& photon : m_photons) {
249 float xp0 = 0., yp0 = 0.;
250 float xp1 = 0., yp1 = 0.;
251 ToPlane(photon[0][0], photon[0][1], photon[0][2], xp0, yp0);
252 ToPlane(photon[1][0], photon[1][1], photon[1][2], xp1, yp1);
253 std::vector<float> xgr = {xp0, xp1};
254 std::vector<float> ygr = {yp0, yp1};
255 gr.DrawGraph(2, xgr.data(), ygr.data(), "Lsame");
256 }
257
258 gr.SetMarkerSize(m_markerSizeCollision);
259 gr.SetMarkerStyle(20);
260 if (!m_exc.empty()) {
261 gr.SetMarkerColor(m_colExcitation);
262 std::vector<float> xgr;
263 std::vector<float> ygr;
264 for (const auto& p : m_exc) {
265 if (!InBox(p)) continue;
266 float xp = 0., yp = 0.;
267 ToPlane(p[0], p[1], p[2], xp, yp);
268 xgr.push_back(xp);
269 ygr.push_back(yp);
270 }
271 if (!xgr.empty()) {
272 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "Psame");
273 }
274 }
275 if (!m_ion.empty()) {
276 gr.SetMarkerColor(m_colIonisation);
277 std::vector<float> xgr;
278 std::vector<float> ygr;
279 for (const auto& p : m_ion) {
280 if (!InBox(p)) continue;
281 float xp = 0., yp = 0.;
282 ToPlane(p[0], p[1], p[2], xp, yp);
283 xgr.push_back(xp);
284 ygr.push_back(yp);
285 }
286 if (!xgr.empty()) {
287 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "Psame");
288 }
289 }
290 if (!m_att.empty()) {
291 gr.SetMarkerColor(m_colAttachment);
292 std::vector<float> xgr;
293 std::vector<float> ygr;
294 for (const auto& p : m_att) {
295 if (!InBox(p)) continue;
296 float xp = 0., yp = 0.;
297 ToPlane(p[0], p[1], p[2], xp, yp);
298 xgr.push_back(xp);
299 ygr.push_back(yp);
300 }
301 if (!xgr.empty()) {
302 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "Psame");
303 }
304 }
305
306 gPad->Update();
307}
std::string LabelY()
Definition: ViewBase.cc:483
std::string LabelX()
Definition: ViewBase.cc:420
bool InBox(const std::array< T, 3 > &x) const
Definition: ViewBase.hh:115
void DrawLine(const std::vector< std::array< float, 3 > > &xl, const short col, const short lw)
Definition: ViewBase.cc:371
bool RangeSet(TPad *)
Definition: ViewBase.cc:83
TPad * GetCanvas()
Retrieve the canvas.
Definition: ViewBase.cc:74
void SetRange(TPad *pad, const double x0, const double y0, const double x1, const double y1)
Definition: ViewBase.cc:93
void ToPlane(const T x, const T y, const T z, T &xp, T &yp) const
Definition: ViewBase.hh:109

Referenced by main(), and Plot().

◆ Plot3d()

void Garfield::ViewDrift::Plot3d ( const bool  axis,
const bool  ogl 
)

Make a 3D plot of the drift lines.

Definition at line 309 of file ViewDrift.cc.

309 {
310 auto pad = GetCanvas();
311 pad->cd();
312 pad->SetTitle("Drift lines");
313 if (!pad->GetView()) {
314 if (!SetPlotLimits3d()) {
315 std::cerr << m_className << "::Plot3d:\n"
316 << " Could not determine the plot limits.\n";
317 }
318 auto view = TView::CreateView(1, 0, 0);
319 view->SetRange(m_xMinBox, m_yMinBox, m_zMinBox,
321 if (axis) view->ShowAxis();
322 pad->SetView(view);
323 if (ogl) pad->GetViewer3D("ogl");
324 }
325 for (const auto& driftLine : m_driftLines) {
326 std::vector<float> points;
327 for (const auto& p : driftLine.first) {
328 points.push_back(p[0]);
329 points.push_back(p[1]);
330 points.push_back(p[2]);
331 }
332 const int nP = driftLine.first.size();
333 TPolyLine3D pl(nP, points.data());
334 if (driftLine.second == Particle::Electron) {
335 pl.SetLineColor(m_colElectron);
336 } else if (driftLine.second == Particle::Hole) {
337 pl.SetLineColor(m_colHole);
338 } else {
339 pl.SetLineColor(m_colIon);
340 }
341 pl.SetLineWidth(1);
342 pl.DrawPolyLine(nP, points.data(), "same");
343 }
344
345 for (const auto& track : m_tracks) {
346 std::vector<float> points;
347 for (const auto& p : track) {
348 points.push_back(p[0]);
349 points.push_back(p[1]);
350 points.push_back(p[2]);
351 }
352 const int nP = track.size();
353 TPolyLine3D pl(nP, points.data());
354 pl.SetLineColor(m_colTrack);
355 pl.SetLineWidth(1);
356 pl.DrawPolyLine(nP, points.data(), "same");
357 if (!m_drawClusters) continue;
358 TPolyMarker3D pm(nP, points.data());
359 pm.SetMarkerColor(m_colTrack);
360 pm.SetMarkerSize(m_markerSizeCluster);
361 pm.DrawPolyMarker(nP, points.data(), 20, "same");
362 }
363
364 if (!m_exc.empty()) {
365 const size_t nP = m_exc.size();
366 std::vector<float> points;
367 for (size_t i = 0; i < nP; ++i) {
368 points.push_back(m_exc[i][0]);
369 points.push_back(m_exc[i][1]);
370 points.push_back(m_exc[i][2]);
371 }
372 TPolyMarker3D pm(nP, points.data());
373 pm.SetMarkerColor(m_colExcitation);
374 pm.SetMarkerSize(m_markerSizeCollision);
375 pm.DrawPolyMarker(nP, points.data(), 20, "same");
376 }
377
378 if (!m_ion.empty()) {
379 const size_t nP = m_ion.size();
380 std::vector<float> points;
381 for (size_t i = 0; i < nP; ++i) {
382 points.push_back(m_ion[i][0]);
383 points.push_back(m_ion[i][1]);
384 points.push_back(m_ion[i][2]);
385 }
386 TPolyMarker3D pm(nP, points.data());
387 pm.SetMarkerColor(m_colIonisation);
388 pm.SetMarkerSize(m_markerSizeCollision);
389 pm.DrawPolyMarker(nP, points.data(), 20, "same");
390 }
391
392 if (!m_att.empty()) {
393 const size_t nP = m_att.size();
394 std::vector<float> points;
395 for (size_t i = 0; i < nP; ++i) {
396 points.push_back(m_att[i][0]);
397 points.push_back(m_att[i][1]);
398 points.push_back(m_att[i][2]);
399 }
400 TPolyMarker3D pm(nP, points.data());
401 pm.SetMarkerColor(m_colAttachment);
402 pm.SetMarkerSize(m_markerSizeCollision);
403 pm.DrawPolyMarker(nP, points.data(), 20, "same");
404 }
405 pad->Modified();
406 pad->Update();
407
408 if (axis) {
409 TAxis3D* ax3d = TAxis3D::GetPadAxis();
410 if (ax3d) {
411 ax3d->SetLabelColor(kGray + 2);
412 ax3d->SetAxisColor(kGray + 2);
413 ax3d->SetXTitle("x");
414 ax3d->SetYTitle("y");
415 ax3d->SetZTitle("z");
416 }
417 pad->Update();
418 }
419}

Referenced by Plot().

◆ SetClusterMarkerSize()

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

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

Definition at line 38 of file ViewDrift.cc.

38 {
39 if (size > 0.) {
40 m_markerSizeCluster = size;
41 } else {
42 std::cerr << m_className << "::SetClusterMarkerSize: Size must be > 0.\n";
43 }
44}

Referenced by main().

◆ SetCollisionMarkerSize()

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

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

Definition at line 46 of file ViewDrift.cc.

46 {
47 if (size > 0.) {
48 m_markerSizeCollision = size;
49 } else {
50 std::cerr << m_className << "::SetCollisionMarkerSize: Size must be > 0.\n";
51 }
52}

Referenced by main().

◆ SetColourAttachments()

void Garfield::ViewDrift::SetColourAttachments ( const short  col)
inline

Set the colour with which to draw attachment markers.

Definition at line 57 of file ViewDrift.hh.

57{ m_colAttachment = col; }

◆ SetColourElectrons()

void Garfield::ViewDrift::SetColourElectrons ( const short  col)
inline

Set the colour with which to draw electron drift lines.

Definition at line 43 of file ViewDrift.hh.

43{ m_colElectron = col; }

◆ SetColourExcitations()

void Garfield::ViewDrift::SetColourExcitations ( const short  col)
inline

Set the colour with which to draw excitation markers.

Definition at line 53 of file ViewDrift.hh.

53{ m_colExcitation = col; }

◆ SetColourHoles()

void Garfield::ViewDrift::SetColourHoles ( const short  col)
inline

Set the colour with which to draw hole drift lines.

Definition at line 45 of file ViewDrift.hh.

45{ m_colHole = col; }

◆ SetColourIonisations()

void Garfield::ViewDrift::SetColourIonisations ( const short  col)
inline

Set the colour with which to draw ionisation markers.

Definition at line 55 of file ViewDrift.hh.

55{ m_colIonisation = col; }

◆ SetColourIons()

void Garfield::ViewDrift::SetColourIons ( const short  col)
inline

Set the colour with which to draw ion drift lines.

Definition at line 47 of file ViewDrift.hh.

47{ m_colIon = col; }

◆ SetColourPhotons()

void Garfield::ViewDrift::SetColourPhotons ( const short  col)
inline

Set the colour with which to draw photons.

Definition at line 51 of file ViewDrift.hh.

51{ m_colPhoton = col; }

◆ SetColourTracks()

void Garfield::ViewDrift::SetColourTracks ( const short  col)
inline

Set the colour with which to draw charged particle tracks.

Definition at line 49 of file ViewDrift.hh.

49{ m_colTrack = col; }

◆ SetDriftLinePoint()

void Garfield::ViewDrift::SetDriftLinePoint ( const size_t  iL,
const size_t  iP,
const float  x,
const float  y,
const float  z 
)

Definition at line 124 of file ViewDrift.cc.

126 {
127 std::lock_guard<std::mutex> guard(m_mutex);
128 if (iL >= m_driftLines.size() || iP >= m_driftLines[iL].first.size()) {
129 std::cerr << m_className << "::SetDriftLinePoint: Index out of range.\n";
130 return;
131 }
132 m_driftLines[iL].first[iP] = {x, y, z};
133}

◆ SetTrackPoint()

void Garfield::ViewDrift::SetTrackPoint ( const size_t  iL,
const size_t  iP,
const float  x,
const float  y,
const float  z 
)

Definition at line 146 of file ViewDrift.cc.

147 {
148 std::lock_guard<std::mutex> guard(m_mutex);
149 if (iL >= m_tracks.size() || iP >= m_tracks[iL].size()) {
150 std::cerr << m_className << "::SetTrackPoint: Index out of range.\n";
151 return;
152 }
153 m_tracks[iL][iP] = {x, y, z};
154}

Friends And Related Function Documentation

◆ ViewFEMesh

friend class ViewFEMesh
friend

Definition at line 91 of file ViewDrift.hh.


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