Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 "GarfieldConstants.hh"
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() = default;
25
26 /// Delete existing drift lines, tracks and markers.
27 void Clear();
28
29 /// Draw the drift lines.
30 void Plot(const bool twod = false, const bool axis = true,
31 const bool snapshot = false);
32 /// Make a 2D plot of the drift lines in the current viewing plane.
33 void Plot2d(const bool axis = true, const bool snapshot = false);
34 /// Make a 3D plot of the drift lines.
35 void Plot3d(const bool axis = true, const bool ogl = true,
36 const bool snapshot = false);
37
38 /// Draw markers (or not) at every collision along a track.
39 void EnableClusterMarkers(const bool on = true) { m_drawClusters = on; }
40 /// Set the size of the cluster markers (see TAttMarker).
41 void SetClusterMarkerSize(const double size);
42 /// Set the size of the collision markers (see TAttMarker).
43 void SetCollisionMarkerSize(const double size);
44
45 /// Set the colour with which to draw electron drift lines.
46 void SetColourElectrons(const short col) { m_colElectron = col; }
47 /// Set the colour with which to draw hole drift lines.
48 void SetColourHoles(const short col) { m_colHole = col; }
49 /// Set the colour with which to draw ion drift lines.
50 void SetColourIons(const short col) { m_colIon = col; }
51 /// Set the colour with which to draw charged particle tracks.
52 void SetColourTracks(const short col) { m_colTrack = col; }
53 /// Set the colour with which to draw photons.
54 void SetColourPhotons(const short col) { m_colPhoton = col; }
55 /// Set the colour with which to draw excitation markers.
56 void SetColourExcitations(const short col) { m_colExcitation = col; }
57 /// Set the colour with which to draw ionisation markers.
58 void SetColourIonisations(const short col) { m_colIonisation = col; }
59 /// Set the colour with which to draw attachment markers.
60 void SetColourAttachments(const short col) { m_colAttachment = col; }
61
62 /// Get the number of drift lines stored.
63 size_t GetNumberOfDriftLines() const { return m_driftLines.size(); }
64 /// Retrieve the coordinates of a given drift line.
65 void GetDriftLine(const size_t i,
66 std::vector<std::array<float, 3> >& driftLine,
67 bool& electron) const;
68
69 // Functions used by the transport classes.
70 size_t NewDriftLine(const Particle particle, const size_t np,
71 const float x0, const float y0, const float z0);
72 void NewChargedParticleTrack(const size_t np, size_t& id, const float x0,
73 const float y0, const float z0);
74
75 void SetDriftLinePoint(const size_t iL, const size_t iP,
76 const float x, const float y, const float z);
77 void AddDriftLinePoint(const size_t iL, const float x, const float y,
78 const float z);
79 void SetTrackPoint(const size_t iL, const size_t iP,
80 const float x, const float y, const float z);
81 void AddTrackPoint(const size_t iL, const float x, const float y,
82 const float z);
83 void AddExcitation(const float x, const float y, const float z);
84 void AddIonisation(const float x, const float y, const float z);
85 void AddAttachment(const float x, const float y, const float z);
86
87 void AddPhoton(const float x0, const float y0, const float z0,
88 const float x1, const float y1, const float z1);
89
90 friend class ViewFEMesh;
91
92 private:
93 std::mutex m_mutex;
94
95 std::vector<std::pair<std::vector<std::array<float, 3> >,
96 Particle> > m_driftLines;
97
98 std::vector<std::vector<std::array<float, 3> > > m_tracks;
99 std::vector<std::array<std::array<float, 3>, 2> > m_photons;
100
101 std::vector<std::array<float, 3> > m_exc;
102 std::vector<std::array<float, 3> > m_ion;
103 std::vector<std::array<float, 3> > m_att;
104
105 double m_markerSizeCluster = 0.01;
106 double m_markerSizeCollision = 0.5;
107
108 short m_colTrack = kGreen + 3;
109 short m_colPhoton = kBlue + 1;
110 short m_colElectron = kOrange - 3;
111 short m_colHole = kRed + 1;
112 short m_colIon = kRed + 1;
113 short m_colExcitation = kGreen + 3;
114 short m_colIonisation = kOrange - 3;
115 short m_colAttachment = kCyan + 3;
116
117 bool m_drawClusters = false;
118
119 bool SetPlotLimits2d();
120 bool SetPlotLimits3d();
121 void DrawMarkers2d(const std::vector<std::array<float, 3> >& points,
122 const short col, const double size);
123 void DrawMarkers3d(const std::vector<std::array<float, 3> >& points,
124 const short col, const double size);
125};
126}
127#endif
ViewBase()=delete
Default constructor.
void AddIonisation(const float x, const float y, const float z)
Definition ViewDrift.cc:150
void Plot2d(const bool axis=true, const bool snapshot=false)
Make a 2D plot of the drift lines in the current viewing plane.
Definition ViewDrift.cc:170
void EnableClusterMarkers(const bool on=true)
Draw markers (or not) at every collision along a track.
Definition ViewDrift.hh:39
void SetColourPhotons(const short col)
Set the colour with which to draw photons.
Definition ViewDrift.hh:54
void AddDriftLinePoint(const size_t iL, const float x, const float y, const float z)
Definition ViewDrift.cc:112
void AddTrackPoint(const size_t iL, const float x, const float y, const float z)
Definition ViewDrift.cc:133
void SetColourExcitations(const short col)
Set the colour with which to draw excitation markers.
Definition ViewDrift.hh:56
void NewChargedParticleTrack(const size_t np, size_t &id, const float x0, const float y0, const float z0)
Definition ViewDrift.cc:88
void SetColourIonisations(const short col)
Set the colour with which to draw ionisation markers.
Definition ViewDrift.hh:58
void SetColourElectrons(const short col)
Set the colour with which to draw electron drift lines.
Definition ViewDrift.hh:46
void SetTrackPoint(const size_t iL, const size_t iP, const float x, const float y, const float z)
Definition ViewDrift.cc:123
friend class ViewFEMesh
Definition ViewDrift.hh:90
void SetColourHoles(const short col)
Set the colour with which to draw hole drift lines.
Definition ViewDrift.hh:48
void Plot(const bool twod=false, const bool axis=true, const bool snapshot=false)
Draw the drift lines.
Definition ViewDrift.cc:162
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:101
void Plot3d(const bool axis=true, const bool ogl=true, const bool snapshot=false)
Make a 3D plot of the drift lines.
Definition ViewDrift.cc:276
void SetColourTracks(const short col)
Set the colour with which to draw charged particle tracks.
Definition ViewDrift.hh:52
void SetColourIons(const short col)
Set the colour with which to draw ion drift lines.
Definition ViewDrift.hh:50
void AddAttachment(const float x, const float y, const float z)
Definition ViewDrift.cc:156
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:60
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:80
size_t NewDriftLine(const Particle particle, const size_t np, 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:63
void AddExcitation(const float x, const float y, const float z)
Definition ViewDrift.cc:144
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
Definition ViewDrift.cc:46