Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization 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 <array>
5
8
9namespace Garfield {
10
11class Sensor;
12class Component;
13
15
16class ViewIsochrons : public ViewBase {
17 public:
21 ~ViewIsochrons() = default;
22
24 void SetSensor(Sensor* s);
26 void SetComponent(Component* c);
27
37 void PlotIsochrons(const double tstep,
38 const std::vector<std::array<double, 3> >& points,
39 const bool rev = false, const bool colour = false,
40 const bool markers = false,
41 const bool plotDriftLines = true);
42
45 void DriftElectrons(const bool positive = false) {
47 m_positive = positive;
48 }
49
50 void DriftIons(const bool negative = false) {
52 m_positive = !negative;
53 }
54
55 void EnableSorting(const bool on = true) { m_sortContours = on; }
58 void CheckCrossings(const bool on = true) { m_checkCrossings = on; }
61 void SetAspectRatioSwitch(const double ar);
64 void SetLoopThreshold(const double thr);
67 void SetConnectionThreshold(const double thr);
68
69 private:
70 Sensor* m_sensor = nullptr;
71 Component* m_component = nullptr;
72
73 // Type of particle to be used for computing drift lines.
75 bool m_positive = false;
76
77 short m_markerStyle = 5;
78 short m_lineStyle = 2;
79
80 bool m_sortContours = true;
81 double m_aspectRatio = 3.;
82 double m_loopThreshold = 0.2;
84 bool m_checkCrossings = true;
85
87
89 const double tstep, const std::vector<std::array<double, 3> >& points,
90 std::vector<std::vector<std::array<double, 3> > >& driftLines,
91 std::vector<std::array<double, 3> >& startPoints,
92 std::vector<std::array<double, 3> >& endPoints,
93 std::vector<int>& statusCodes, const bool rev = false);
95 std::vector<std::pair<std::array<double, 4>, unsigned int> >& contour,
96 bool& circle);
97};
98} // namespace Garfield
99#endif
ViewBase()=delete
Default constructor.
void SetComponent(Component *c)
Set the component.
void CheckCrossings(const bool on=true)
Check (or not) that drift-lines do not cross isochrons (default: check is done).
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)
Draw equal time contour lines.
void ComputeDriftLines(const double tstep, const std::vector< std::array< double, 3 > > &points, std::vector< std::vector< std::array< double, 3 > > > &driftLines, std::vector< std::array< double, 3 > > &startPoints, std::vector< std::array< double, 3 > > &endPoints, std::vector< int > &statusCodes, const bool rev=false)
void SortContour(std::vector< std::pair< std::array< double, 4 >, unsigned int > > &contour, bool &circle)
void SetConnectionThreshold(const double thr)
Fractional distance over which isochron segments are connected (default: 0.2).
~ViewIsochrons()=default
Destructor.
void DriftElectrons(const bool positive=false)
Request electron drift lines with negative (default) or positive charge.
void SetAspectRatioSwitch(const double ar)
Set the aspect ratio above which an isochron is considered linear (as opposed to circular).
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)
Fractional distance between two points for closing a circular isochron (default: 0....
void SetSensor(Sensor *s)
Set the sensor.