Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewIsochrons.hh
Go to the documentation of this file.
1#ifndef G_VIEW_ISOCHRONS
2#define G_VIEW_ISOCHRONS
3
4#include "ViewBase.hh"
5
6namespace Garfield {
7
8class Sensor;
9class ComponentBase;
10
11/// Draw equal time contour lines.
12
13class ViewIsochrons : public ViewBase {
14 public:
15 /// Constructor.
17 /// Destructor.
18 ~ViewIsochrons() = default;
19
20 /// Set the sensor.
21 void SetSensor(Sensor* s);
22 /// Set the component.
24
25 /// Set the viewing area (in local coordinates of the current viewing plane).
26 void SetArea(const double xmin, const double ymin, const double xmax,
27 const double ymax);
28 /// Set the viewing area based on the bounding box of the sensor/component.
29 void SetArea() { m_hasUserArea = false; }
30
31 /** Set the projection (viewing plane).
32 * \param fx,fy,fz normal vector
33 * \param x0,y0,z0 in-plane point
34 */
35 void SetPlane(const double fx, const double fy, const double fz,
36 const double x0, const double y0, const double z0);
37 /// Set the default viewing plane (\f$x\f$-\f$y\f$ at \f$z = 0\f$).
39 /// Rotate the viewing plane (angle in radian).
40 void Rotate(const double angle);
41
42 /** Draw equal time contour lines.
43 * \param tstep Time interval.
44 * \param points List of starting points.
45 * \param rev If true, the drift time is measured from the end
46 * points of the drift lines.
47 * \param colour Draw contour lines using colours.
48 * \param markers Draw markers (as opposed to lines).
49 * \param plotDriftLines Draw drift lines together with the isochrons.
50 */
51 void PlotIsochrons(const double tstep,
52 const std::vector<std::array<double, 3> >& points,
53 const bool rev = false,
54 const bool colour = false, const bool markers = false,
55 const bool plotDriftLines = true);
56
57 /// Request electron drift lines with negative (default) or
58 /// positive charge.
59 void DriftElectrons(const bool positive = false) {
60 m_particle = Particle::Electron;
61 m_positive = positive;
62 }
63 /// Request ion drift lines with positive (default) or negative charge.
64 void DriftIons(const bool negative = false) {
65 m_particle = Particle::Ion;
66 m_positive = !negative;
67 }
68 /// Sort (or not) the points on a contour line (default: sorting is done).
69 void EnableSorting(const bool on = true) { m_sortContours = on; }
70 /// Check (or not) that drift-lines do not cross isochrons
71 /// (default: check is done).
72 void CheckCrossings(const bool on = true) { m_checkCrossings = on; }
73 /// Set the aspect ratio above which an isochron is considered
74 /// linear (as opposed to circular).
75 void SetAspectRatioSwitch(const double ar);
76 /// Fractional distance between two points for closing a
77 /// circular isochron (default: 0.2).
78 void SetLoopThreshold(const double thr);
79 /// Fractional distance over which isochron segments are connected
80 /// (default: 0.2).
81 void SetConnectionThreshold(const double thr);
82
83 private:
84 Sensor* m_sensor = nullptr;
85 ComponentBase* m_component = nullptr;
86
87 // Projection for viewing
88 double m_proj[3][3];
89 double m_plane[4];
90 char m_xLabel[50], m_yLabel[50], m_description[50];
91
92 // Plot area
93 bool m_hasUserArea = false;
94 double m_xmin = -1., m_ymin = -1.;
95 double m_xmax = 1., m_ymax = 1.;
96
97 enum class Particle { Electron = 0, Ion };
98 // Type of particle to be used for computing drift lines.
99 Particle m_particle = Particle::Electron;
100 bool m_positive = false;
101
102 short m_markerStyle = 5;
103 short m_lineStyle = 2;
104
105 bool m_sortContours = true;
106 double m_aspectRatio = 3.;
107 double m_loopThreshold = 0.2;
108 double m_connectionThreshold = 0.2;
109 bool m_checkCrossings = true;
110
111 void SetupCanvas();
112 bool Range();
113 void Labels();
114
115 void ComputeDriftLines(const double tstep,
116 const std::vector<std::array<double, 3> >& points,
117 std::vector<std::vector<std::array<double, 3> > >& driftLines,
118 std::vector<std::array<double, 3> >& startPoints,
119 std::vector<std::array<double, 3> >& endPoints,
120 std::vector<int>& statusCodes, const bool rev = false);
121 void SortContour(
122 std::vector<std::pair<std::array<double, 4>, unsigned int> >& contour,
123 bool& circle);
124
125};
126}
127#endif
Abstract base class for components.
Base class for visualization classes.
Definition: ViewBase.hh:10
Draw equal time contour lines.
void CheckCrossings(const bool on=true)
void SetArea()
Set the viewing area based on the bounding box of the sensor/component.
void DriftIons(const bool negative=false)
Request ion drift lines with positive (default) or negative charge.
void PlotIsochrons(const double tstep, const std::vector< std::array< double, 3 > > &points, const bool rev=false, const bool colour=false, const bool markers=false, const bool plotDriftLines=true)
void SetConnectionThreshold(const double thr)
~ViewIsochrons()=default
Destructor.
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
void DriftElectrons(const bool positive=false)
void SetAspectRatioSwitch(const double ar)
void EnableSorting(const bool on=true)
Sort (or not) the points on a contour line (default: sorting is done).
ViewIsochrons()
Constructor.
void SetDefaultProjection()
Set the default viewing plane ( - at ).
void SetLoopThreshold(const double thr)
void SetSensor(Sensor *s)
Set the sensor.
void SetComponent(ComponentBase *c)
Set the component.
void Rotate(const double angle)
Rotate the viewing plane (angle in radian).