Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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, const bool snapshot=false)
 Draw the drift lines.
 
void Plot2d (const bool axis=true, const bool snapshot=false)
 Make a 2D plot of the drift lines in the current viewing plane.
 
void Plot3d (const bool axis=true, const bool ogl=true, const bool snapshot=false)
 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.
 
size_t NewDriftLine (const Particle particle, const size_t np, 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 SetPlaneZX ()
 Set the viewing plane to z-x.
 
void SetPlaneZY ()
 Set the viewing plane to z-y.
 
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
 
- Static Protected Member Functions inherited from Garfield::ViewBase
static bool RangeSet (TVirtualPad *)
 
static void SetRange (TVirtualPad *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 19 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 156 of file ViewDrift.cc.

156 {
157 std::lock_guard<std::mutex> guard(m_mutex);
158 std::array<float, 3> p = {x, y, z};
159 m_att.push_back(std::move(p));
160}

◆ AddDriftLinePoint()

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

Definition at line 112 of file ViewDrift.cc.

113 {
114 std::lock_guard<std::mutex> guard(m_mutex);
115 if (iL >= m_driftLines.size()) {
116 std::cerr << m_className << "::AddDriftLinePoint: Index out of range.\n";
117 return;
118 }
119 std::array<float, 3> p = {x, y, z};
120 m_driftLines[iL].first.push_back(std::move(p));
121}
std::string m_className
Definition ViewBase.hh:82

◆ AddExcitation()

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

Definition at line 144 of file ViewDrift.cc.

144 {
145 std::lock_guard<std::mutex> guard(m_mutex);
146 std::array<float, 3> p = {x, y, z};
147 m_exc.push_back(std::move(p));
148}

◆ AddIonisation()

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

Definition at line 150 of file ViewDrift.cc.

150 {
151 std::lock_guard<std::mutex> guard(m_mutex);
152 std::array<float, 3> p = {x, y, z};
153 m_ion.push_back(std::move(p));
154}

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

81 {
82 std::lock_guard<std::mutex> guard(m_mutex);
83 std::array<float, 3> p0 = {x0, y0, z0};
84 std::array<float, 3> p1 = {x1, y1, z1};
85 m_photons.push_back({p0, p1});
86}

◆ AddTrackPoint()

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

Definition at line 133 of file ViewDrift.cc.

134 {
135 std::lock_guard<std::mutex> guard(m_mutex);
136 if (iL >= m_tracks.size()) {
137 std::cerr << m_className << "::AddTrackPoint: Index out of range.\n";
138 return;
139 }
140 std::array<float, 3> p = {x, y, z};
141 m_tracks[iL].push_back(std::move(p));
142}

◆ 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 39 of file ViewDrift.hh.

39{ 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 63 of file ViewDrift.hh.

63{ 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 88 of file ViewDrift.cc.

90 {
91 std::lock_guard<std::mutex> guard(m_mutex);
92 // Create a new track and add it to the list.
93 constexpr size_t minSize = 1;
94 std::vector<std::array<float, 3> > track(std::max(minSize, np));
95 track[0] = {x0, y0, z0};
96 m_tracks.push_back(std::move(track));
97 // Return the index of this track.
98 id = m_tracks.size() - 1;
99}

◆ NewDriftLine()

size_t Garfield::ViewDrift::NewDriftLine ( const Particle particle,
const size_t np,
const float x0,
const float y0,
const float z0 )

Definition at line 68 of file ViewDrift.cc.

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

◆ Plot()

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

Draw the drift lines.

Definition at line 162 of file ViewDrift.cc.

162 {
163 if (twod) {
164 Plot2d(axis, snapshot);
165 } else {
166 Plot3d(axis, false, snapshot);
167 }
168}
void Plot2d(const bool axis=true, const bool snapshot=false)
Make a 2D plot of the drift lines in the current viewing plane.
Definition ViewDrift.cc:170
void Plot3d(const bool axis=true, const bool ogl=true, const bool snapshot=false)
Make a 3D plot of the drift lines.
Definition ViewDrift.cc:276

Referenced by main().

◆ Plot2d()

void Garfield::ViewDrift::Plot2d ( const bool axis = true,
const bool snapshot = false )

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

Definition at line 170 of file ViewDrift.cc.

170 {
171 auto pad = GetCanvas();
172 pad->cd();
173 pad->SetTitle("Drift lines");
174 // Check if the canvas range has already been set.
175 const bool rangeSet = RangeSet(pad);
176 if (axis || !rangeSet) {
177 // Determine the plot limits.
178 if (!SetPlotLimits2d()) {
179 std::cerr << m_className << "::Plot2d:\n"
180 << " Could not determine the plot limits.\n";
181 return;
182 }
183 }
184 if (axis) {
185 auto frame = pad->DrawFrame(m_xMinPlot, m_yMinPlot,
187 frame->GetXaxis()->SetTitle(LabelX().c_str());
188 frame->GetYaxis()->SetTitle(LabelY().c_str());
189 } else if (!rangeSet) {
191 }
192 if (snapshot) {
193 std::vector<std::array<float, 3> > electrons;
194 std::vector<std::array<float, 3> > holes;
195 std::vector<std::array<float, 3> > ions;
196 for (const auto& driftLine : m_driftLines) {
197 if (driftLine.second == Particle::Electron) {
198 electrons.push_back(driftLine.first.back());
199 } else if (driftLine.second == Particle::Hole) {
200 holes.push_back(driftLine.first.back());
201 } else {
202 ions.push_back(driftLine.first.back());
203 }
204 }
205 DrawMarkers2d(electrons, m_colElectron, m_markerSizeCollision);
206 DrawMarkers2d(holes, m_colHole, m_markerSizeCollision);
207 DrawMarkers2d(ions, m_colIon, m_markerSizeCollision);
208 } else {
209 for (const auto& driftLine : m_driftLines) {
210 const short lw = 1;
211 if (driftLine.second == Particle::Electron) {
212 DrawLine(driftLine.first, m_colElectron, lw);
213 } else if (driftLine.second == Particle::Hole) {
214 DrawLine(driftLine.first, m_colHole, lw);
215 } else {
216 DrawLine(driftLine.first, m_colIon, lw);
217 }
218 }
219 }
220 gPad->Update();
221
222 for (const auto& track : m_tracks) {
223 DrawLine(track, m_colTrack, 2);
224 if (!m_drawClusters) continue;
225 DrawMarkers2d(track, m_colTrack, m_markerSizeCluster);
226 }
227
228 TGraph gr;
229 gr.SetLineColor(m_colPhoton);
230 gr.SetLineStyle(2);
231 for (const auto& photon : m_photons) {
232 float xp0 = 0., yp0 = 0.;
233 float xp1 = 0., yp1 = 0.;
234 ToPlane(photon[0][0], photon[0][1], photon[0][2], xp0, yp0);
235 ToPlane(photon[1][0], photon[1][1], photon[1][2], xp1, yp1);
236 std::vector<float> xgr = {xp0, xp1};
237 std::vector<float> ygr = {yp0, yp1};
238 gr.DrawGraph(2, xgr.data(), ygr.data(), "Lsame");
239 }
240
241 if (!m_exc.empty()) {
242 DrawMarkers2d(m_exc, m_colExcitation, m_markerSizeCollision);
243 }
244 if (!m_ion.empty()) {
245 DrawMarkers2d(m_ion, m_colIonisation, m_markerSizeCollision);
246 }
247 if (!m_att.empty()) {
248 DrawMarkers2d(m_att, m_colAttachment, m_markerSizeCollision);
249 }
250
251 gPad->Update();
252}
std::string LabelY()
Definition ViewBase.cc:510
static void SetRange(TVirtualPad *pad, const double x0, const double y0, const double x1, const double y1)
Definition ViewBase.cc:93
std::string LabelX()
Definition ViewBase.cc:447
void DrawLine(const std::vector< std::array< float, 3 > > &xl, const short col, const short lw)
Definition ViewBase.cc:398
TPad * GetCanvas()
Retrieve the canvas.
Definition ViewBase.cc:74
static bool RangeSet(TVirtualPad *)
Definition ViewBase.cc:83
void ToPlane(const T x, const T y, const T z, T &xp, T &yp) const
Definition ViewBase.hh:113

Referenced by main(), and Plot().

◆ Plot3d()

void Garfield::ViewDrift::Plot3d ( const bool axis = true,
const bool ogl = true,
const bool snapshot = false )

Make a 3D plot of the drift lines.

Definition at line 276 of file ViewDrift.cc.

277 {
278 auto pad = GetCanvas();
279 pad->cd();
280 pad->SetTitle("Drift lines");
281 if (!pad->GetView()) {
282 if (!SetPlotLimits3d()) {
283 std::cerr << m_className << "::Plot3d:\n"
284 << " Could not determine the plot limits.\n";
285 }
286 auto view = TView::CreateView(1, 0, 0);
287 view->SetRange(m_xMinBox, m_yMinBox, m_zMinBox,
289 if (axis) view->ShowAxis();
290 pad->SetView(view);
291 if (ogl) pad->GetViewer3D("ogl");
292 }
293
294 if (snapshot) {
295 std::vector<std::array<float, 3> > electrons;
296 std::vector<std::array<float, 3> > holes;
297 std::vector<std::array<float, 3> > ions;
298 for (const auto& driftLine : m_driftLines) {
299 if (driftLine.second == Particle::Electron) {
300 electrons.push_back(driftLine.first.back());
301 } else if (driftLine.second == Particle::Hole) {
302 holes.push_back(driftLine.first.back());
303 } else {
304 ions.push_back(driftLine.first.back());
305 }
306 }
307 DrawMarkers3d(electrons, m_colElectron, m_markerSizeCollision);
308 DrawMarkers3d(holes, m_colHole, m_markerSizeCollision);
309 DrawMarkers3d(ions, m_colIon, m_markerSizeCollision);
310 } else {
311 for (const auto& driftLine : m_driftLines) {
312 std::vector<float> points;
313 for (const auto& p : driftLine.first) {
314 points.push_back(p[0]);
315 points.push_back(p[1]);
316 points.push_back(p[2]);
317 }
318 const int nP = driftLine.first.size();
319 TPolyLine3D pl(nP, points.data());
320 if (driftLine.second == Particle::Electron) {
321 pl.SetLineColor(m_colElectron);
322 } else if (driftLine.second == Particle::Hole) {
323 pl.SetLineColor(m_colHole);
324 } else {
325 pl.SetLineColor(m_colIon);
326 }
327 pl.SetLineWidth(1);
328 pl.DrawPolyLine(nP, points.data(), "same");
329 }
330 }
331 for (const auto& track : m_tracks) {
332 std::vector<float> points;
333 for (const auto& p : track) {
334 points.push_back(p[0]);
335 points.push_back(p[1]);
336 points.push_back(p[2]);
337 }
338 const int nP = track.size();
339 TPolyLine3D pl(nP, points.data());
340 pl.SetLineColor(m_colTrack);
341 pl.SetLineWidth(1);
342 pl.DrawPolyLine(nP, points.data(), "same");
343 if (m_drawClusters) {
344 DrawMarkers3d(track, m_colTrack, m_markerSizeCluster);
345 }
346 }
347
348 if (!m_exc.empty()) {
349 DrawMarkers3d(m_exc, m_colExcitation, m_markerSizeCollision);
350 }
351 if (!m_ion.empty()) {
352 DrawMarkers3d(m_ion, m_colIonisation, m_markerSizeCollision);
353 }
354 if (!m_att.empty()) {
355 DrawMarkers3d(m_att, m_colAttachment, m_markerSizeCollision);
356 }
357 pad->Modified();
358 pad->Update();
359
360 if (axis) {
361 TAxis3D* ax3d = TAxis3D::GetPadAxis();
362 if (ax3d) {
363 ax3d->SetLabelColor(kGray + 2);
364 ax3d->SetAxisColor(kGray + 2);
365 ax3d->SetXTitle("x");
366 ax3d->SetYTitle("y");
367 ax3d->SetZTitle("z");
368 }
369 pad->Update();
370 }
371}

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 60 of file ViewDrift.hh.

60{ 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 46 of file ViewDrift.hh.

46{ m_colElectron = col; }

◆ SetColourExcitations()

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

Set the colour with which to draw excitation markers.

Definition at line 56 of file ViewDrift.hh.

56{ 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 48 of file ViewDrift.hh.

48{ m_colHole = col; }

◆ SetColourIonisations()

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

Set the colour with which to draw ionisation markers.

Definition at line 58 of file ViewDrift.hh.

58{ 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 50 of file ViewDrift.hh.

50{ m_colIon = col; }

◆ SetColourPhotons()

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

Set the colour with which to draw photons.

Definition at line 54 of file ViewDrift.hh.

54{ 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 52 of file ViewDrift.hh.

52{ 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 101 of file ViewDrift.cc.

103 {
104 std::lock_guard<std::mutex> guard(m_mutex);
105 if (iL >= m_driftLines.size() || iP >= m_driftLines[iL].first.size()) {
106 std::cerr << m_className << "::SetDriftLinePoint: Index out of range.\n";
107 return;
108 }
109 m_driftLines[iL].first[iP] = {x, y, z};
110}

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

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

Friends And Related Symbol Documentation

◆ ViewFEMesh

friend class ViewFEMesh
friend

Definition at line 90 of file ViewDrift.hh.

Referenced by ViewFEMesh.


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