Garfield++ 4.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 <string>
5#include <vector>
6#include <array>
7#include <utility>
8#include <mutex>
9
10#include <Rtypes.h>
11
12#include "ViewBase.hh"
13
14namespace Garfield {
15
16/// Visualize drift lines and tracks.
17
18class ViewDrift : public ViewBase {
19 public:
20 /// Constructor.
21 ViewDrift();
22 /// Destructor.
23 ~ViewDrift() = default;
24
25 /// Delete existing drift lines, tracks and markers.
26 void Clear();
27
28 /// Draw the drift lines.
29 void Plot(const bool twod = false, const bool axis = true);
30 /// Make a 2D plot of the drift lines in the current viewing plane.
31 void Plot2d(const bool axis);
32 /// Make a 3D plot of the drift lines.
33 void Plot3d(const bool axis, const bool ogl);
34
35 /// Draw markers (or not) at every collision along a track.
36 void EnableClusterMarkers(const bool on = true) { m_drawClusters = on; }
37 /// Set the size of the cluster markers (see TAttMarker).
38 void SetClusterMarkerSize(const double size);
39 /// Set the size of the collision markers (see TAttMarker).
40 void SetCollisionMarkerSize(const double size);
41
42 /// Set the colour with which to draw electron drift lines.
43 void SetColourElectrons(const short col) { m_colElectron = col; }
44 /// Set the colour with which to draw hole drift lines.
45 void SetColourHoles(const short col) { m_colHole = col; }
46 /// Set the colour with which to draw ion drift lines.
47 void SetColourIons(const short col) { m_colIon = col; }
48 /// Set the colour with which to draw charged particle tracks.
49 void SetColourTracks(const short col) { m_colTrack = col; }
50 /// Set the colour with which to draw photons.
51 void SetColourPhotons(const short col) { m_colPhoton = col; }
52 /// Set the colour with which to draw excitation markers.
53 void SetColourExcitations(const short col) { m_colExcitation = col; }
54 /// Set the colour with which to draw ionisation markers.
55 void SetColourIonisations(const short col) { m_colIonisation = col; }
56 /// Set the colour with which to draw attachment markers.
57 void SetColourAttachments(const short col) { m_colAttachment = col; }
58
59 /// Get the number of drift lines stored.
60 size_t GetNumberOfDriftLines() const { return m_driftLines.size(); }
61 /// Retrieve the coordinates of a given drift line.
62 void GetDriftLine(const size_t i,
63 std::vector<std::array<float, 3> >& driftLine,
64 bool& electron) const;
65
66 // Functions used by the transport classes.
67 void NewElectronDriftLine(const size_t np, size_t& id, const float x0,
68 const float y0, const float z0);
69 void NewHoleDriftLine(const size_t np, size_t& id, const float x0,
70 const float y0, const float z0);
71 void NewIonDriftLine(const size_t np, size_t& id, const float x0,
72 const float y0, const float z0);
73 void NewChargedParticleTrack(const size_t np, size_t& id, const float x0,
74 const float y0, const float z0);
75
76 void SetDriftLinePoint(const size_t iL, const size_t iP,
77 const float x, const float y, const float z);
78 void AddDriftLinePoint(const size_t iL, const float x, const float y,
79 const float z);
80 void SetTrackPoint(const size_t iL, const size_t iP,
81 const float x, const float y, const float z);
82 void AddTrackPoint(const size_t iL, const float x, const float y,
83 const float z);
84 void AddExcitation(const float x, const float y, const float z);
85 void AddIonisation(const float x, const float y, const float z);
86 void AddAttachment(const float x, const float y, const float z);
87
88 void AddPhoton(const float x0, const float y0, const float z0,
89 const float x1, const float y1, const float z1);
90
91 friend class ViewFEMesh;
92
93 private:
94 std::mutex m_mutex;
95
96 enum class Particle {
97 Electron,
98 Hole,
99 Ion
100 };
101 std::vector<std::pair<std::vector<std::array<float, 3> >,
102 Particle> > m_driftLines;
103
104 std::vector<std::vector<std::array<float, 3> > > m_tracks;
105 std::vector<std::array<std::array<float, 3>, 2> > m_photons;
106
107 std::vector<std::array<float, 3> > m_exc;
108 std::vector<std::array<float, 3> > m_ion;
109 std::vector<std::array<float, 3> > m_att;
110
111 double m_markerSizeCluster = 0.01;
112 double m_markerSizeCollision = 0.5;
113
114 short m_colTrack = kGreen + 3;
115 short m_colPhoton = kBlue + 1;
116 short m_colElectron = kOrange - 3;
117 short m_colHole = kRed + 1;
118 short m_colIon = kRed + 1;
119 short m_colExcitation = kGreen + 3;
120 short m_colIonisation = kOrange - 3;
121 short m_colAttachment = kCyan + 3;
122
123 bool m_drawClusters = false;
124
125 bool SetPlotLimits2d();
126 bool SetPlotLimits3d();
127};
128}
129#endif
Base class for visualization classes.
Definition: ViewBase.hh:18
Visualize drift lines and tracks.
Definition: ViewDrift.hh:18
void AddIonisation(const float x, const float y, const float z)
Definition: ViewDrift.cc:173
void EnableClusterMarkers(const bool on=true)
Draw markers (or not) at every collision along a track.
Definition: ViewDrift.hh:36
void NewHoleDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
Definition: ViewDrift.cc:81
void SetColourPhotons(const short col)
Set the colour with which to draw photons.
Definition: ViewDrift.hh:51
void AddDriftLinePoint(const size_t iL, const float x, const float y, const float z)
Definition: ViewDrift.cc:135
void Plot(const bool twod=false, const bool axis=true)
Draw the drift lines.
Definition: ViewDrift.cc:185
void AddTrackPoint(const size_t iL, const float x, const float y, const float z)
Definition: ViewDrift.cc:156
void SetColourExcitations(const short col)
Set the colour with which to draw excitation markers.
Definition: ViewDrift.hh:53
void NewChargedParticleTrack(const size_t np, size_t &id, const float x0, const float y0, const float z0)
Definition: ViewDrift.cc:111
void NewIonDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
Definition: ViewDrift.cc:92
void SetColourIonisations(const short col)
Set the colour with which to draw ionisation markers.
Definition: ViewDrift.hh:55
void SetColourElectrons(const short col)
Set the colour with which to draw electron drift lines.
Definition: ViewDrift.hh:43
void SetTrackPoint(const size_t iL, const size_t iP, const float x, const float y, const float z)
Definition: ViewDrift.cc:146
void SetColourHoles(const short col)
Set the colour with which to draw hole drift lines.
Definition: ViewDrift.hh:45
void SetClusterMarkerSize(const double size)
Set the size of the cluster markers (see TAttMarker).
Definition: ViewDrift.cc:38
void SetDriftLinePoint(const size_t iL, const size_t iP, const float x, const float y, const float z)
Definition: ViewDrift.cc:124
void Plot3d(const bool axis, const bool ogl)
Make a 3D plot of the drift lines.
Definition: ViewDrift.cc:309
void Plot2d(const bool axis)
Make a 2D plot of the drift lines in the current viewing plane.
Definition: ViewDrift.cc:193
void SetColourTracks(const short col)
Set the colour with which to draw charged particle tracks.
Definition: ViewDrift.hh:49
void SetColourIons(const short col)
Set the colour with which to draw ion drift lines.
Definition: ViewDrift.hh:47
void AddAttachment(const float x, const float y, const float z)
Definition: ViewDrift.cc:179
ViewDrift()
Constructor.
Definition: ViewDrift.cc:21
~ViewDrift()=default
Destructor.
void GetDriftLine(const size_t i, std::vector< std::array< float, 3 > > &driftLine, bool &electron) const
Retrieve the coordinates of a given drift line.
Definition: ViewDrift.cc:54
void SetColourAttachments(const short col)
Set the colour with which to draw attachment markers.
Definition: ViewDrift.hh:57
void Clear()
Delete existing drift lines, tracks and markers.
Definition: ViewDrift.cc:29
void AddPhoton(const float x0, const float y0, const float z0, const float x1, const float y1, const float z1)
Definition: ViewDrift.cc:103
void NewElectronDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
Definition: ViewDrift.cc:68
size_t GetNumberOfDriftLines() const
Get the number of drift lines stored.
Definition: ViewDrift.hh:60
void AddExcitation(const float x, const float y, const float z)
Definition: ViewDrift.cc:167
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
Definition: ViewDrift.cc:46
Draw the mesh of a field-map component.
Definition: ViewFEMesh.hh:21