REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorReadoutModule.h
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see http://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
23#ifndef RestCore_TRestDetectorReadoutModule
24#define RestCore_TRestDetectorReadoutModule
25
26#include <TMath.h>
27#include <TVector2.h>
28
29#include <iostream>
30
31#include "TRestDetectorReadoutChannel.h"
32#include "TRestDetectorReadoutMapping.h"
33
37 private:
38 Int_t fId = -1;
39
40 TVector2 fOrigin = {0, 0};
41
42 TVector2 fSize = {0, 0};
43
45 Double_t fRotation = 0; //<
46
47 std::pair<Int_t, Int_t> fDaqIdRange = {
48 -1, -1};
49
50 std::vector<TRestDetectorReadoutChannel>
53
55
56 Double_t fTolerance;
58
59 Bool_t showWarnings;
61
62 Int_t fFirstDaqChannel = 0;
63
64 std::string fDecodingFile = "";
65
66 Int_t fMappingNodes = 0;
67
68 Bool_t fDecoding;
70
71 std::string fName; //<
72 std::string fType; //<
73
74 void Initialize();
75
78 inline TVector2 TransformToModuleCoordinates(const TVector2& coords) const {
79 return TVector2(coords - fOrigin).Rotate(-1.0 * fRotation);
80 }
81
84 inline TVector2 TransformToPlaneCoordinates(const TVector2& coords) const {
85 return coords.Rotate(fRotation) + fOrigin;
86 }
87
88 protected:
89 public:
90 // Setters
91
93 inline void SetModuleID(Int_t modID) { fId = modID; }
94
96 inline void SetSize(const TVector2& size) { fSize = size; }
97
99 inline void SetOrigin(const TVector2& origin) { fOrigin = origin; }
100
102 inline void SetRotation(Double_t rotation) { fRotation = rotation; }
103
105 inline void SetName(const std::string& name) { fName = name; }
106
108 inline void SetType(const std::string& type) { fType = type; }
109
111 inline void SetTolerance(Double_t tolerance) { fTolerance = tolerance; }
112
114 inline void SetFirstDaqChannel(Int_t firstDaqChannel) { fFirstDaqChannel = firstDaqChannel; }
115
117 inline void SetMappingNodes(Int_t nodes) { fMappingNodes = nodes; }
118
120 inline Double_t GetTolerance() const { return fTolerance; }
121
123 inline Int_t GetMinDaqID() const { return fDaqIdRange.first; }
124
126 inline Int_t GetMaxDaqID() const { return fDaqIdRange.second; }
127
128 inline Int_t GetMappingNodes() const { return fMappingNodes; }
129
131 inline Int_t DaqToReadoutChannel(Int_t daqChannel) {
132 for (size_t n = 0; n < GetNumberOfChannels(); n++) {
133 if (GetChannel(n)->GetDaqID() == daqChannel) {
134 return n;
135 }
136 }
137 return -1;
138 }
139
141 inline Int_t GetModuleID() const { return fId; }
142
144 inline TVector2 GetOrigin() const { return fOrigin; }
145
147 inline TVector2 GetSize() const { return fSize; }
148
150 inline Double_t GetRotation() const { return fRotation; }
151
154 TVector2 GetModuleCoordinates(const TVector2& p) const { return TransformToModuleCoordinates(p); }
155
158 TVector2 GetPlaneCoordinates(const TVector2& p) { return TransformToPlaneCoordinates(p); }
159
161 inline std::string GetName() const { return fName; }
162 inline std::string GetType() const { return fType; }
163
166
167 inline TRestDetectorReadoutChannel& operator[](int n) { return fReadoutChannel[n]; }
168
171 if (n >= GetNumberOfChannels()) {
172 return nullptr;
173 }
174 return &fReadoutChannel[n];
175 }
176
178 inline size_t GetNumberOfChannels() const { return fReadoutChannel.size(); }
179
181 inline void EnableWarnings() { showWarnings = true; }
182
184 inline void DisableWarnings() { showWarnings = false; }
185
186 void DoReadoutMapping();
187
188 void SetDecodingFile(const std::string& decodingFile);
189
194 Bool_t IsInside(const TVector2& position) const;
195
196 Bool_t IsInsideChannel(Int_t channel, const TVector2& position);
197
198 Bool_t IsInsidePixel(Int_t channel, Int_t pixel, const TVector2& position);
199
200 Bool_t IsDaqIDInside(Int_t daqID);
201 Int_t FindChannel(const TVector2& position);
202 TVector2 GetDistanceToModule(const TVector2& position);
203
204 TVector2 GetPixelOrigin(Int_t channel, Int_t pixel);
205 TVector2 GetPixelVertex(Int_t channel, Int_t pixel, Int_t vertex);
206 TVector2 GetPixelCenter(Int_t channel, Int_t pixel);
207 Bool_t GetPixelTriangle(Int_t channel, Int_t pixel);
208
210 TVector2 GetPixelVertex(TRestDetectorReadoutPixel* pix, Int_t vertex);
213
214 TVector2 GetVertex(int n) const;
215
217
218 void SetMinMaxDaqIDs();
219
220 void Draw();
221
222 void Print(Int_t DetailLevel = 0);
223
224 // Constructor
226 // Destructor
228
229 ClassDef(TRestDetectorReadoutModule, 5);
230};
231#endif
Int_t GetDaqID() const
Returns the corresponding daq channel id.
TRestDetectorReadoutModule()
Default TRestDetectorReadoutModule constructor.
std::vector< TRestDetectorReadoutChannel > fReadoutChannel
void DisableWarnings()
Disables warning output.
void SetSize(const TVector2 &size)
Sets the module size by definition using TVector2 input.
void AddChannel(TRestDetectorReadoutChannel &channel)
Adds a new channel to the module.
void SetType(const std::string &type)
Sets the type of the readout module.
TVector2 TransformToPlaneCoordinates(const TVector2 &coords) const
Int_t GetModuleID() const
Returns the module id.
Bool_t IsInsidePixel(Int_t channel, Int_t pixel, const TVector2 &position)
Determines if the position TVector2 pos is found at a specific pixel id inside the readout channel gi...
void SetMappingNodes(Int_t nodes)
Sets number of nodes.
TVector2 GetPixelCenter(Int_t channel, Int_t pixel)
Returns the center pixel position for a given channel and pixel indexes.
Int_t GetMinDaqID() const
Returns the minimum daq id number.
void Initialize()
TRestDetectorReadoutModule initialization.
*default REST_Warning in TRestDetectorReadout will enable it *Int_t fFirstDaqChannel
First DAQ channel.
Double_t GetTolerance() const
Gets the tolerance for independent pixel overlaps.
TVector2 GetSize() const
Returns the module size (x, y) in mm.
void DoReadoutMapping()
Starts the readout mapping initialization. This process is computationally expensive but it greatly o...
TVector2 fOrigin
The module (x, y) position relative to the readout plane position.
std::string fDecodingFile
Decoding file.
TVector2 GetModuleCoordinates(const TVector2 &p) const
void SetDecodingFile(const std::string &decodingFile)
Set the decoding file in the readout module.
TRestDetectorReadoutChannel * GetChannel(size_t n)
Returns a pointer to a readout channel by index.
TVector2 GetPixelOrigin(Int_t channel, Int_t pixel)
Returns the pixel origin (left-bottom) position for a given channel and pixel indexes.
TVector2 fSize
The module (x, y) size. All pixels should be contained within this size.
void SetModuleID(Int_t modID)
Sets the module by id definition.
TVector2 TransformToModuleCoordinates(const TVector2 &coords) const
Bool_t GetPixelTriangle(Int_t channel, Int_t pixel)
Returns the pixel type for a given channel and pixel indexes.
void SetFirstDaqChannel(Int_t firstDaqChannel)
Sets first DAQ channel.
Bool_t IsInside(const TVector2 &position) const
Determines if the position x,y relative to the readout plane are inside this readout module.
Int_t FindChannel(const TVector2 &position)
Returns the channel index corresponding to the absolute coordinates (absX, absY), but relative to the...
Bool_t IsInsideChannel(Int_t channel, const TVector2 &position)
Determines if the position TVector2 pos is found in any of the pixels of the readout channel index gi...
TVector2 GetOrigin() const
Returns the module origin position.
std::string GetName() const
Returns the module name.
Double_t GetRotation() const
Returns the module rotation in degrees.
void EnableWarnings()
Enables warning output.
TVector2 GetPixelVertex(Int_t channel, Int_t pixel, Int_t vertex)
Returns any of the pixel vertex position for a given channel and pixel indexes.
TVector2 GetPlaneCoordinates(const TVector2 &p)
TRestDetectorReadoutMapping * GetMapping()
Returns a pointer to the readout mapping.
TVector2 GetVertex(int n) const
Returns the coordinates of the specified vertex index n. The physical coordinates relative to the rea...
Bool_t IsDaqIDInside(Int_t daqID)
Determines if a given daqID number is in the range of the module.
void SetTolerance(Double_t tolerance)
Sets the tolerance for independent pixel overlaps.
void SetRotation(Double_t rotation)
Sets the module rotation in degrees.
Double_t fRotation
The rotation of the module around the module origin (fModuleOriginX, fModuleOriginY) in radians.
Int_t GetMaxDaqID() const
Returns the maximum daq id number.
TRestDetectorReadoutMapping fMapping
The readout module uniform grid mapping.
void SetMinMaxDaqIDs()
Initializes the max and min values for the daq channel number.
Int_t fId
The module id given by the readout definition.
void SetName(const std::string &name)
Sets the name of the readout module.
virtual ~TRestDetectorReadoutModule()
Default TRestDetectorReadoutModule destructor.
void Print(Int_t DetailLevel=0)
Prints the module details and channels if fullDetail is enabled.
void SetOrigin(const TVector2 &origin)
Sets the module origin by definition using TVector2 input.
size_t GetNumberOfChannels() const
Returns the total number of channels defined inside the module.
std::pair< Int_t, Int_t > fDaqIdRange
The minimum and maximum daq channel ids associated to the module.
Int_t DaqToReadoutChannel(Int_t daqChannel)
Returns the physical readout channel index for a given daq id channel number.
TVector2 GetDistanceToModule(const TVector2 &position)
Creates a TVector2 with shortest norm going from the given position pos to the module....
A class to store the readout pixel definition used in TRestDetectorReadoutChannel.