Garfield++ 3.0
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 <memory>
5#include <string>
6#include <vector>
7
8#include <TGraph.h>
9#include <TPolyLine3D.h>
10#include <TPolyMarker3D.h>
11#include <TView.h>
12
13#include "ViewBase.hh"
14
15namespace Garfield {
16
17/// Visualize drift lines and tracks.
18
19class ViewDrift : public ViewBase {
20 public:
21 /// Constructor.
22 ViewDrift();
23 /// Destructor.
24 ~ViewDrift();
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 friend class ViewFEMesh;
65
66 private:
67
68 struct Marker {
69 double x;
70 double y;
71 double z;
72 };
73
74 // Box dimensions
75 double m_xMin = -1., m_yMin = -1., m_zMin = -1.;
76 double m_xMax = 1., m_yMax = 1., m_zMax = 1.;
77
78 // View
79 std::unique_ptr<TView> m_view;
80
81 struct DriftLine {
82 std::vector<Marker> points;
83 int n; // what kind of particle?
84 };
85 std::vector<DriftLine> m_driftLines;
86 std::vector<TPolyLine3D> m_driftLinePlots;
87
88 std::vector<std::vector<Marker> > m_tracks;
89 std::vector<TPolyMarker3D> m_trackPlots;
90 std::vector<TPolyLine3D> m_trackLinePlots;
91
92 std::vector<Marker> m_excMarkers;
93 std::unique_ptr<TPolyMarker3D> m_excPlot;
94 std::vector<Marker> m_ionMarkers;
95 std::unique_ptr<TPolyMarker3D> m_ionPlot;
96 std::vector<Marker> m_attMarkers;
97 std::unique_ptr<TPolyMarker3D> m_attPlot;
98
99 double m_markerSizeCluster = 1.;
100 double m_markerSizeCollision = 1.;
101
102 void Plot2d(const bool axis);
103 void Plot3d(const bool axis);
104};
105}
106#endif
Base class for visualization classes.
Definition: ViewBase.hh:10
Visualize drift lines and tracks.
Definition: ViewDrift.hh:19
void NewElectronDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:74
void AddExcitationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:211
void AddIonisationMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:220
void Plot(const bool twod=false, const bool axis=true)
Draw the drift lines.
Definition: ViewDrift.cc:238
void SetDriftLinePoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:162
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:26
void AddAttachmentMarker(const double x, const double y, const double z)
Definition: ViewDrift.cc:229
void AddTrackPoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:198
void SetClusterMarkerSize(const double size)
Set the size of the cluster markers (see TAttMarker).
Definition: ViewDrift.cc:58
void NewChargedParticleTrack(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:149
void AddDriftLinePoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:174
ViewDrift()
Constructor.
Definition: ViewDrift.cc:12
void SetTrackPoint(const unsigned int iL, const unsigned int iP, const double x, const double y, const double z)
Definition: ViewDrift.cc:187
void Clear()
Delete existing drift lines, tracks and markers.
Definition: ViewDrift.cc:42
void NewPhotonTrack(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Definition: ViewDrift.cc:138
void NewHoleDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:99
~ViewDrift()
Destructor.
Definition: ViewDrift.cc:22
void NewIonDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:119
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
Definition: ViewDrift.cc:66
Draw the mesh of a field-map component.
Definition: ViewFEMesh.hh:26