Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewDrift.hh
Go to the documentation of this file.
1#ifndef G_VIEW_DRIFT
2#define G_VIEW_DRIFT
3
4#include <string>
5#include <TGraph.h>
6#include <TCanvas.h>
7#include <TPolyLine3D.h>
8#include <TPolyMarker3D.h>
9#include <TView.h>
10
11namespace Garfield {
12
13/// Visualize drift lines and tracks.
14
15class ViewDrift {
16
17 public:
18 /// Constructor
19 ViewDrift();
20 /// Destructor
21 ~ViewDrift();
22
23 /// Set the canvas to be painted on.
24 void SetCanvas(TCanvas* c);
25
26 /// Set the region to be plotted.
27 void SetArea(const double xmin, const double ymin, const double zmin,
28 const double xmax, const double ymax, const double zmax);
29 /// Delete existing drift lines, tracks and markers.
30 void Clear();
31
32 /// Draw the drift lines.
33 void Plot(const bool twod = false, const bool axis = true);
34
35 /// Set the size of the cluster markers (see TAttMarker).
36 void SetClusterMarkerSize(const double size);
37 /// Set the size of the collision markers (see TAttMarker).
38 void SetCollisionMarkerSize(const double size);
39
40 // Functions used by the transport classes.
41 void NewElectronDriftLine(const unsigned int np, int& id, const double x0,
42 const double y0, const double z0);
43 void NewHoleDriftLine(const unsigned int np, int& id, const double x0,
44 const double y0, const double z0);
45 void NewIonDriftLine(const unsigned int np, int& id, const double x0,
46 const double y0, const double z0);
47 void NewPhotonTrack(const double x0, const double y0, const double z0,
48 const double x1, const double y1, const double z1);
49 void NewChargedParticleTrack(const unsigned int np, int& id, const double x0,
50 const double y0, const double z0);
51
52 void SetDriftLinePoint(const unsigned int iL, const unsigned int iP,
53 const double x, const double y, const double z);
54 void AddDriftLinePoint(const unsigned int iL, const double x, const double y,
55 const double z);
56 void SetTrackPoint(const unsigned int iL, const unsigned int iP,
57 const double x, const double y, const double z);
58 void AddTrackPoint(const unsigned int iL, const double x, const double y,
59 const double z);
60 void AddExcitationMarker(const double x, const double y, const double z);
61 void AddIonisationMarker(const double x, const double y, const double z);
62 void AddAttachmentMarker(const double x, const double y, const double z);
63
64 /// Switch on/off debugging output.
65 void EnableDebugging(const bool on = true) { m_debug = on; }
66
67 friend class ViewFEMesh;
68
69 private:
70 std::string m_className;
71
72 // Options
73 bool m_debug;
74
75 struct Marker {
76 double x;
77 double y;
78 double z;
79 };
80 // Canvas
81 TCanvas* m_canvas;
82 bool m_hasExternalCanvas;
83
84 // Box dimensions
85 double m_xMin, m_yMin, m_zMin;
86 double m_xMax, m_yMax, m_zMax;
87 // View
88 TView* m_view;
89
90 struct driftLine {
91 std::vector<Marker> vect;
92 int n; // what kind of particle?
93 };
94 std::vector<driftLine> m_driftLines;
95 std::vector<TPolyLine3D*> m_driftLinePlots;
96
97 std::vector<std::vector<Marker> > m_tracks;
98 std::vector<TPolyMarker3D*> m_trackPlots;
99 std::vector<TPolyLine3D*> m_trackLinePlots;
100
101 std::vector<Marker> m_excMarkers;
102 TPolyMarker3D* m_excPlot;
103 std::vector<Marker> m_ionMarkers;
104 TPolyMarker3D* m_ionPlot;
105 std::vector<Marker> m_attMarkers;
106 TPolyMarker3D* m_attPlot;
107
108 double m_markerSizeCluster;
109 double m_markerSizeCollision;
110
111 void Plot2d(const bool axis);
112 void Plot3d(const bool axis);
113};
114}
115#endif
Visualize drift lines and tracks.
Definition: ViewDrift.hh:15
void NewElectronDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:127
void AddExcitationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:283
void EnableDebugging(const bool on=true)
Switch on/off debugging output.
Definition: ViewDrift.hh:65
void AddIonisationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:293
void Plot(const bool twod=false, const bool axis=true)
Draw the drift lines.
Definition: ViewDrift.cc:313
void SetDriftLinePoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:230
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:56
void AddAttachmentMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:303
void AddTrackPoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:269
void SetClusterMarkerSize(const double size)
Set the size of the cluster markers (see TAttMarker).
Definition: ViewDrift.cc:109
void NewChargedParticleTrack(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:216
void AddDriftLinePoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:243
ViewDrift()
Constructor.
Definition: ViewDrift.cc:11
void SetTrackPoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:257
void Clear()
Delete existing drift lines, tracks and markers.
Definition: ViewDrift.cc:74
void SetCanvas(TCanvas *c)
Set the canvas to be painted on.
Definition: ViewDrift.cc:45
void NewPhotonTrack(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Definition: ViewDrift.cc:203
void NewHoleDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:156
~ViewDrift()
Destructor.
Definition: ViewDrift.cc:38
void NewIonDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:180
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
Definition: ViewDrift.cc:118
Draw the mesh of a field-map component.
Definition: ViewFEMesh.hh:26