Garfield++ 4.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 Component;
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.
23 void SetComponent(Component* c);
24
25 /** Draw equal time contour lines.
26 * \param tstep Time interval.
27 * \param points List of starting points.
28 * \param rev If true, the drift time is measured from the end
29 * points of the drift lines.
30 * \param colour Draw contour lines using colours.
31 * \param markers Draw markers (as opposed to lines).
32 * \param plotDriftLines Draw drift lines together with the isochrons.
33 */
34 void PlotIsochrons(const double tstep,
35 const std::vector<std::array<double, 3> >& points,
36 const bool rev = false,
37 const bool colour = false, const bool markers = false,
38 const bool plotDriftLines = true);
39
40 /// Request electron drift lines with negative (default) or
41 /// positive charge.
42 void DriftElectrons(const bool positive = false) {
43 m_particle = Particle::Electron;
44 m_positive = positive;
45 }
46 /// Request ion drift lines with positive (default) or negative charge.
47 void DriftIons(const bool negative = false) {
48 m_particle = Particle::Ion;
49 m_positive = !negative;
50 }
51 /// Sort (or not) the points on a contour line (default: sorting is done).
52 void EnableSorting(const bool on = true) { m_sortContours = on; }
53 /// Check (or not) that drift-lines do not cross isochrons
54 /// (default: check is done).
55 void CheckCrossings(const bool on = true) { m_checkCrossings = on; }
56 /// Set the aspect ratio above which an isochron is considered
57 /// linear (as opposed to circular).
58 void SetAspectRatioSwitch(const double ar);
59 /// Fractional distance between two points for closing a
60 /// circular isochron (default: 0.2).
61 void SetLoopThreshold(const double thr);
62 /// Fractional distance over which isochron segments are connected
63 /// (default: 0.2).
64 void SetConnectionThreshold(const double thr);
65
66 private:
67 Sensor* m_sensor = nullptr;
68 Component* m_component = nullptr;
69
70 enum class Particle { Electron = 0, Ion };
71 // Type of particle to be used for computing drift lines.
72 Particle m_particle = Particle::Electron;
73 bool m_positive = false;
74
75 short m_markerStyle = 5;
76 short m_lineStyle = 2;
77
78 bool m_sortContours = true;
79 double m_aspectRatio = 3.;
80 double m_loopThreshold = 0.2;
81 double m_connectionThreshold = 0.2;
82 bool m_checkCrossings = true;
83
84 bool SetPlotLimits();
85
86 void ComputeDriftLines(const double tstep,
87 const std::vector<std::array<double, 3> >& points,
88 std::vector<std::vector<std::array<double, 3> > >& driftLines,
89 std::vector<std::array<double, 3> >& startPoints,
90 std::vector<std::array<double, 3> >& endPoints,
91 std::vector<int>& statusCodes, const bool rev = false);
92 void SortContour(
93 std::vector<std::pair<std::array<double, 4>, unsigned int> >& contour,
94 bool& circle);
95
96};
97}
98#endif
Abstract base class for components.
Definition: Component.hh:13
Base class for visualization classes.
Definition: ViewBase.hh:18
Draw equal time contour lines.
void SetComponent(Component *c)
Set the component.
void CheckCrossings(const bool on=true)
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 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 SetLoopThreshold(const double thr)
void SetSensor(Sensor *s)
Set the sensor.