Garfield++ 3.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 ()
 Destructor.
 
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)
 
- Public Member Functions inherited from Garfield::ViewBase
 ViewBase ()=delete
 Default constructor.
 
 ViewBase (const std::string &name)
 Constructor.
 
virtual ~ViewBase ()
 Destructor.
 
void SetCanvas (TCanvas *c)
 Set the canvas to be painted on.
 
TCanvas * GetCanvas ()
 Retrieve the canvas.
 
void EnableDebugging (const bool on=true)
 Switch on/off debugging output.
 

Friends

class ViewFEMesh
 

Additional Inherited Members

- Protected Member Functions inherited from Garfield::ViewBase
std::string FindUnusedFunctionName (const std::string &s) const
 
std::string FindUnusedHistogramName (const std::string &s) const
 
- Protected Attributes inherited from Garfield::ViewBase
std::string m_className = "ViewBase"
 
bool m_debug = false
 
TCanvas * m_canvas = nullptr
 
bool m_hasExternalCanvas = false
 
double m_proj [3][3]
 

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

12 : ViewBase("ViewDrift") {
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}
ViewBase()=delete
Default constructor.

◆ ~ViewDrift()

Garfield::ViewDrift::~ViewDrift ( )

Destructor.

Definition at line 22 of file ViewDrift.cc.

22 {
23 Clear();
24}
void Clear()
Delete existing drift lines, tracks and markers.
Definition: ViewDrift.cc:42

Member Function Documentation

◆ AddAttachmentMarker()

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

Definition at line 229 of file ViewDrift.cc.

230 {
231 Marker newMarker;
232 newMarker.x = x;
233 newMarker.y = y;
234 newMarker.z = z;
235 m_attMarkers.push_back(std::move(newMarker));
236}

◆ AddDriftLinePoint()

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

Definition at line 174 of file ViewDrift.cc.

175 {
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}
std::string m_className
Definition: ViewBase.hh:28

◆ AddExcitationMarker()

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

Definition at line 211 of file ViewDrift.cc.

212 {
213 Marker newMarker;
214 newMarker.x = x;
215 newMarker.y = y;
216 newMarker.z = z;
217 m_excMarkers.push_back(std::move(newMarker));
218}

◆ AddIonisationMarker()

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

Definition at line 220 of file ViewDrift.cc.

221 {
222 Marker newMarker;
223 newMarker.x = x;
224 newMarker.y = y;
225 newMarker.z = z;
226 m_ionMarkers.push_back(std::move(newMarker));
227}

◆ AddTrackPoint()

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

Definition at line 198 of file ViewDrift.cc.

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

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

◆ Clear()

void Garfield::ViewDrift::Clear ( )

Delete existing drift lines, tracks and markers.

Definition at line 42 of file ViewDrift.cc.

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

Referenced by ~ViewDrift().

◆ NewChargedParticleTrack()

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

Definition at line 149 of file ViewDrift.cc.

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

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

76 {
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}
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 99 of file ViewDrift.cc.

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

◆ NewIonDriftLine()

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

Definition at line 119 of file ViewDrift.cc.

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

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

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

◆ Plot()

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

Draw the drift lines.

Definition at line 238 of file ViewDrift.cc.

238 {
239 if (twod) {
240 Plot2d(axis);
241 } else {
242 Plot3d(axis);
243 }
244}

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

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

Referenced by main().

◆ SetClusterMarkerSize()

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

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

Definition at line 58 of file ViewDrift.cc.

58 {
59 if (size > 0.) {
60 m_markerSizeCluster = size;
61 } else {
62 std::cerr << m_className << "::SetClusterMarkerSize: Size must be > 0.\n";
63 }
64}

◆ SetCollisionMarkerSize()

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

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

Definition at line 66 of file ViewDrift.cc.

66 {
67 if (size > 0.) {
68 m_markerSizeCollision = size;
69 } else {
70 std::cerr << m_className << "::SetCollisionMarkerSize: Size must be > 0.\n";
71 }
72}

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

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

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

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

Friends And Related Function Documentation

◆ ViewFEMesh

friend class ViewFEMesh
friend

Definition at line 64 of file ViewDrift.hh.


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