-
Notifications
You must be signed in to change notification settings - Fork 4
/
Copy pathTRestDetectorExperimentalReadoutPixel.h
151 lines (124 loc) · 4.94 KB
/
TRestDetectorExperimentalReadoutPixel.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
/**
* @file TRestDetectorExperimentalReadoutPixel.h
* @brief Defines the TRestDetectorExperimentalReadoutPixel class.
*/
#ifndef REST_TRESTDETECTOREXPERIMENTALREADOUTPIXEL_H
#define REST_TRESTDETECTOREXPERIMENTALREADOUTPIXEL_H
#include <TVector2.h>
#include <TVector3.h>
#include <iostream>
#include <map>
#include <stdexcept>
#include <vector>
/**
* @class TRestDetectorExperimentalReadoutPixel
* @brief Represents a pixel in an experimental readout detector.
* @author Luis Obis
*/
class TRestDetectorExperimentalReadoutPixel {
private:
std::vector<TVector2> fVertices; /**< The vertices of the pixel. */
TVector2 fCenter; /**< The center of the pixel in mm. A pixel is guaranteed to be inside a circle with
center fCenter and radius fRadius. */
double fRadius; /**< The radius of the pixel in mm. A pixel is guaranteed to be inside a circle with
center fCenter and radius fRadius. */
unsigned short fChannel; /**< The channel of the pixel. */
/**
@brief Initializes the vertices of the pixel.
@param vertices The vertices of the pixel.
*/
void InitializeVertices(const std::vector<TVector2>& vertices);
/**
@brief Calculates the vertices of a rectangular pixel.
@param center The center of the pixel.
@param size The size of the pixel (width, height).
@return The vertices of the rectangular pixel.
*/
static std::vector<TVector2> GetRectangularVertices(const TVector2& center,
const std::pair<double, double>& size);
public:
/**
@brief Constructs a TRestDetectorExperimentalReadoutPixel object with given vertices.
@param vertices The vertices of the pixel.
*/
explicit TRestDetectorExperimentalReadoutPixel(const std::vector<TVector2>& vertices,
unsigned short channel = 0);
/**
@brief Constructs a rectangular TRestDetectorExperimentalReadoutPixel object.
@param center The center of the pixel.
@param size The size of the pixel (width, height).
*/
TRestDetectorExperimentalReadoutPixel(const TVector2& center, const std::pair<double, double>& size,
unsigned short channel = 0);
/**
@brief Constructs a square TRestDetectorExperimentalReadoutPixel object.
@param center The center of the pixel.
@param size The size of the pixel (width and height).
*/
TRestDetectorExperimentalReadoutPixel(const TVector2& center, double size, unsigned short channel = 0);
// needed for root
TRestDetectorExperimentalReadoutPixel() = default;
virtual ~TRestDetectorExperimentalReadoutPixel() = default;
/**
@brief Returns the center of the pixel.
@return The center of the pixel.
*/
TVector2 GetCenter() const { return fCenter; }
/**
@brief Returns the radius of the pixel.
@return The radius of the pixel.
*/
double GetRadius() const { return fRadius; }
/**
@brief Returns the vertices of the pixel.
@return The vertices of the pixel.
*/
std::vector<TVector2> GetVertices() const { return fVertices; }
/**
@brief Checks if a given point is inside the pixel.
@param point The point to check.
@return True if the point is inside the pixel, false otherwise.
*/
bool IsInside(const TVector2& point) const;
unsigned short GetChannel() const { return fChannel; }
ClassDef(TRestDetectorExperimentalReadoutPixel, 1);
};
namespace readout {
/**
@brief Calculates the triple product of three vectors.
@param A The first vector.
@param B The second vector.
@param C The third vector.
@return The cross product of the vectors.
*/
double crossProduct(const TVector2& A, const TVector2& B, const TVector2& C);
/**
@brief Calculates the squared distance between two points.
@param A The first point.
@param B The second point.
@return The squared distance between the points.
*/
double squaredDistance(const TVector2& A, const TVector2& B);
/**
@brief Comparator function to compare two points based on their polar angle with respect to the anchor point.
@param A The first point.
@param B The second point.
@param anchor The anchor point.
@return True if point A is smaller than point B, false otherwise.
*/
bool comparePoints(const TVector2& A, const TVector2& B, const TVector2& anchor);
/**
@brief Finds the convex hull of a set of points using the Graham's scan algorithm.
@param _points The input points.
@return The convex hull of the points.
*/
std::vector<TVector2> findConvexHull(const std::vector<TVector2>& _points);
/**
@brief Checks if a given point is inside a convex hull.
@param point The point to check.
@param convexHull The convex hull.
@return True if the point is inside the convex hull, false otherwise.
*/
bool IsPointInsideConvexHull(const TVector2& point, const std::vector<TVector2>& convexHull);
} // namespace readout
#endif // REST_TRESTDETECTOREXPERIMENTALREADOUTPIXEL_H