Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization 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 <Rtypes.h>
5
6#include <array>
7#include <cstddef>
8#include <mutex>
9#include <vector>
10
12#include "Garfield/ViewBase.hh"
13
14namespace Garfield {
15
17
18class ViewDrift : public ViewBase {
19 public:
23 ~ViewDrift() = default;
24
26 void Clear();
27
29 void Plot(const bool twod = false, const bool axis = true,
30 const bool snapshot = false);
32 void Plot2d(const bool axis = true, const bool snapshot = false);
34 void Plot3d(const bool axis = true, const bool ogl = true,
35 const bool snapshot = false);
36
38 void EnableClusterMarkers(const bool on = true) { m_drawClusters = on; }
40 void SetClusterMarkerSize(const double size);
42 void SetCollisionMarkerSize(const double size);
43
45 void SetColourElectrons(const short col) { m_colElectron = col; }
47 void SetColourHoles(const short col) { m_colHole = col; }
49 void SetColourNegativeIons(const short col) { m_colNegativeIon = col; }
51 void SetColourIons(const short col) { m_colIon = col; }
53 void SetColourTracks(const short col) { m_colTrack = col; }
55 void SetColourPhotons(const short col) { m_colPhoton = col; }
57 void SetColourExcitations(const short col) { m_colExcitation = col; }
59 void SetColourIonisations(const short col) { m_colIonisation = col; }
61 void SetColourAttachments(const short col) { m_colAttachment = col; }
63 void SetElectronsToFront() { m_eTop = true; }
64
66 std::size_t GetNumberOfDriftLines() const { return m_driftLines.size(); }
68 void GetDriftLine(const std::size_t i,
69 std::vector<std::array<float, 3> >& driftLine,
70 bool& electron) const;
71
72 // Functions used by the transport classes.
73 std::size_t NewDriftLine(const Particle particle, const std::size_t np,
74 const float x0, const float y0, const float z0);
75 void NewChargedParticleTrack(const std::size_t np, std::size_t& id,
76 const float x0, const float y0, const float z0);
77
78 void SetDriftLinePoint(const std::size_t iL, const std::size_t iP,
79 const float x, const float y, const float z);
80 void AddDriftLinePoint(const std::size_t iL, const float x, const float y,
81 const float z);
82 void SetTrackPoint(const std::size_t iL, const std::size_t iP, const float x,
83 const float y, const float z);
84 void AddTrackPoint(const std::size_t iL, const float x, const float y,
85 const float z);
86 void AddExcitation(const float x, const float y, const float z);
87 void AddIonisation(const float x, const float y, const float z);
88 void AddAttachment(const float x, const float y, const float z);
89
90 void AddPhoton(const float x0, const float y0, const float z0, const float x1,
91 const float y1, const float z1);
92
93 friend class ViewFEMesh;
94
95 private:
96 std::mutex m_mutex;
97
98 std::vector<std::pair<std::vector<std::array<float, 3> >, Particle> >
100
101 std::vector<std::vector<std::array<float, 3> > > m_tracks;
102 std::vector<std::array<std::array<float, 3>, 2> > m_photons;
103
104 std::vector<std::array<float, 3> > m_exc;
105 std::vector<std::array<float, 3> > m_ion;
106 std::vector<std::array<float, 3> > m_att;
107
108 double m_markerSizeCluster = 0.01;
110
111 short m_colTrack = kGreen + 3;
112 short m_colPhoton = kBlue + 1;
113 short m_colElectron = kOrange - 3;
114 short m_colHole = kRed + 1;
115 short m_colIon = kRed + 1;
116 short m_colNegativeIon = kGray + 2;
117 short m_colExcitation = kGreen + 3;
118 short m_colIonisation = kOrange - 3;
119 short m_colAttachment = kCyan + 3;
120
121 bool m_drawClusters = false;
122
123 bool m_eTop = false;
124
127 void DrawMarkers2d(const std::vector<std::array<float, 3> >& points,
128 const short col, const double size);
129 void DrawMarkers3d(const std::vector<std::array<float, 3> >& points,
130 const short col, const double size);
131};
132} // namespace Garfield
133#endif
ViewBase()=delete
Default constructor.
void AddIonisation(const float x, const float y, const float z)
std::vector< std::array< float, 3 > > m_exc
Definition ViewDrift.hh:104
std::size_t NewDriftLine(const Particle particle, const std::size_t np, const float x0, const float y0, const float z0)
void Plot2d(const bool axis=true, const bool snapshot=false)
Make a 2D plot of the drift lines in the current viewing plane.
void EnableClusterMarkers(const bool on=true)
Draw markers (or not) at every collision along a track.
Definition ViewDrift.hh:38
void SetDriftLinePoint(const std::size_t iL, const std::size_t iP, const float x, const float y, const float z)
void AddTrackPoint(const std::size_t iL, const float x, const float y, const float z)
std::mutex m_mutex
Definition ViewDrift.hh:96
void DrawMarkers3d(const std::vector< std::array< float, 3 > > &points, const short col, const double size)
void SetColourNegativeIons(const short col)
Set the colour with which to draw negative ion drift lines.
Definition ViewDrift.hh:49
void SetColourPhotons(const short col)
Set the colour with which to draw photons.
Definition ViewDrift.hh:55
std::vector< std::array< std::array< float, 3 >, 2 > > m_photons
Definition ViewDrift.hh:102
std::vector< std::array< float, 3 > > m_att
Definition ViewDrift.hh:106
std::vector< std::pair< std::vector< std::array< float, 3 > >, Particle > > m_driftLines
Definition ViewDrift.hh:99
void SetColourExcitations(const short col)
Set the colour with which to draw excitation markers.
Definition ViewDrift.hh:57
std::size_t GetNumberOfDriftLines() const
Get the number of drift lines stored.
Definition ViewDrift.hh:66
void SetColourIonisations(const short col)
Set the colour with which to draw ionisation markers.
Definition ViewDrift.hh:59
void SetColourElectrons(const short col)
Set the colour with which to draw electron drift lines.
Definition ViewDrift.hh:45
void SetTrackPoint(const std::size_t iL, const std::size_t iP, const float x, const float y, const float z)
friend class ViewFEMesh
Definition ViewDrift.hh:93
void SetColourHoles(const short col)
Set the colour with which to draw hole drift lines.
Definition ViewDrift.hh:47
void DrawMarkers2d(const std::vector< std::array< float, 3 > > &points, const short col, const double size)
void Plot(const bool twod=false, const bool axis=true, const bool snapshot=false)
Draw the drift lines.
void SetClusterMarkerSize(const double size)
Set the size of the cluster markers (see TAttMarker).
std::vector< std::vector< std::array< float, 3 > > > m_tracks
Definition ViewDrift.hh:101
double m_markerSizeCollision
Definition ViewDrift.hh:109
void SetElectronsToFront()
Put electron drift lines on top of the ion/hole drift lines.
Definition ViewDrift.hh:63
void Plot3d(const bool axis=true, const bool ogl=true, const bool snapshot=false)
Make a 3D plot of the drift lines.
void SetColourTracks(const short col)
Set the colour with which to draw charged particle tracks.
Definition ViewDrift.hh:53
void SetColourIons(const short col)
Set the colour with which to draw ion drift lines.
Definition ViewDrift.hh:51
void AddAttachment(const float x, const float y, const float z)
std::vector< std::array< float, 3 > > m_ion
Definition ViewDrift.hh:105
ViewDrift()
Constructor.
~ViewDrift()=default
Destructor.
void SetColourAttachments(const short col)
Set the colour with which to draw attachment markers.
Definition ViewDrift.hh:61
void Clear()
Delete existing drift lines, tracks and markers.
void AddPhoton(const float x0, const float y0, const float z0, const float x1, const float y1, const float z1)
void GetDriftLine(const std::size_t i, std::vector< std::array< float, 3 > > &driftLine, bool &electron) const
Retrieve the coordinates of a given drift line.
void AddExcitation(const float x, const float y, const float z)
void AddDriftLinePoint(const std::size_t iL, const float x, const float y, const float z)
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
void NewChargedParticleTrack(const std::size_t np, std::size_t &id, const float x0, const float y0, const float z0)