Garfield++ v1r0
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
13class ViewDrift {
14
15 public:
16 // Constructor
17 ViewDrift();
18 // Destructor
19 ~ViewDrift();
20
21 void SetCanvas(TCanvas* c);
22
23 // Set area to be plotted.
24 void SetArea(const double& xmin, const double& ymin, const double& zmin,
25 const double& xmax, const double& ymax, const double& zmax);
26 void Clear();
27 void Plot(const bool twod = false, const bool axis = true);
28
29 void SetClusterMarkerSize(const double& size);
30 void SetCollisionMarkerSize(const double& size);
31
32 // Functions to be used by transport classes.
33 void NewElectronDriftLine(const unsigned int np, int& id, const double x0,
34 const double y0, const double z0);
35 void NewHoleDriftLine(const unsigned int np, int& id, const double x0,
36 const double y0, const double z0);
37 void NewIonDriftLine(const unsigned int np, int& id, const double x0,
38 const double y0, const double z0);
39 void NewPhotonTrack(const double x0, const double y0, const double z0,
40 const double x1, const double y1, const double z1);
41 void NewChargedParticleTrack(const unsigned int np, int& id, const double x0,
42 const double y0, const double z0);
43
44 void SetDriftLinePoint(const unsigned int iL, const unsigned int iP,
45 const double x, const double y, const double z);
46 void AddDriftLinePoint(const unsigned int iL, const double x, const double y,
47 const double z);
48 void SetTrackPoint(const unsigned int iL, const unsigned int iP,
49 const double x, const double y, const double z);
50 void AddTrackPoint(const unsigned int iL, const double x, const double y,
51 const double z);
52 void AddExcitationMarker(const double x, const double y, const double z);
53 void AddIonisationMarker(const double x, const double y, const double z);
54 void AddAttachmentMarker(const double x, const double y, const double z);
55
56 void EnableDebugging() { m_debug = true; }
57 void DisableDebugging() { m_debug = false; }
58
59 friend class ViewFEMesh;
60
61 private:
62 std::string m_className;
63
64 // Options
65 bool m_debug;
66
67 std::string m_label;
68
69 struct marker {
70 double x;
71 double y;
72 double z;
73 };
74 // Canvas
75 TCanvas* m_canvas;
76 bool m_hasExternalCanvas;
77
78 // Box dimensions
79 double m_xMin, m_yMin, m_zMin;
80 double m_xMax, m_yMax, m_zMax;
81 // View
82 TView* m_view;
83
84 unsigned int m_nDriftLines;
85 struct driftLine {
86 std::vector<marker> vect;
87 int n; // what kind of particle?
88 };
89 std::vector<driftLine> m_driftLines;
90 std::vector<TPolyLine3D*> m_driftLinePlots;
91
92 unsigned int m_nTracks;
93 struct track {
94 std::vector<marker> vect;
95 };
96 std::vector<track> m_tracks;
97 std::vector<TPolyMarker3D*> m_trackPlots;
98 std::vector<TPolyLine3D*> m_trackLinePlots;
99
100 unsigned int m_nExcMarkers;
101 std::vector<marker> m_excMarkers;
102 TPolyMarker3D* m_excPlot;
103 unsigned int m_nIonMarkers;
104 std::vector<marker> m_ionMarkers;
105 TPolyMarker3D* m_ionPlot;
106 unsigned int m_nAttMarkers;
107 std::vector<marker> m_attMarkers;
108 TPolyMarker3D* m_attPlot;
109
110 double m_markerSizeCluster;
111 double m_markerSizeCollision;
112
113 void Plot2d(const bool axis);
114 void Plot3d(const bool axis);
115};
116}
117#endif
void NewElectronDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:140
void AddExcitationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:312
void EnableDebugging()
Definition: ViewDrift.hh:56
void AddIonisationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:323
void SetArea(const double &xmin, const double &ymin, const double &zmin, const double &xmax, const double &ymax, const double &zmax)
Definition: ViewDrift.cc:62
void Plot(const bool twod=false, const bool axis=true)
Definition: ViewDrift.cc:345
void SetDriftLinePoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:255
void AddAttachmentMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:334
void AddTrackPoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:297
void SetCollisionMarkerSize(const double &size)
Definition: ViewDrift.cc:130
void SetClusterMarkerSize(const double &size)
Definition: ViewDrift.cc:120
void NewChargedParticleTrack(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:233
void AddDriftLinePoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:269
void SetTrackPoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:284
void DisableDebugging()
Definition: ViewDrift.hh:57
void SetCanvas(TCanvas *c)
Definition: ViewDrift.cc:51
void NewPhotonTrack(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Definition: ViewDrift.cc:219
void NewHoleDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:170
void NewIonDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:195