REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorReadoutPlane.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_TRestDetectorReadoutPlane
24 #define RestCore_TRestDetectorReadoutPlane
25 
26 #include <TGraph.h>
27 #include <TH2Poly.h>
28 
29 #include <iostream>
30 
31 #include "TRestDetectorReadoutChannel.h"
32 #include "TRestDetectorReadoutModule.h"
33 
37  private:
39  Int_t fId = -1; //<
40 
42  TVector3 fPosition = {0, 0, 0}; //<
43 
45  TVector3 fNormal = {0, 0, 1}; //<
46 
48  TVector3 fAxisX; //<
49 
51  TVector3 fAxisY; //<
52 
54  Double_t fChargeCollection = 1; //<
55 
57  Double_t fHeight = 0; //<
58 
60  Double_t fRotation = 0; //<
61 
62  std::string fName; //<
63  std::string fType; //<
64 
66  std::vector<TRestDetectorReadoutModule> fReadoutModules; //<
67 
68  void UpdateAxes();
69 
70  public:
71  // Setters
73  void SetID(int id) { fId = id; }
74 
76  void SetPosition(const TVector3& position) { fPosition = position; }
77 
78  void SetNormal(const TVector3& normal);
79 
81  void SetChargeCollection(Double_t charge) { fChargeCollection = charge; }
82 
83  void SetHeight(Double_t d);
84 
85  void SetRotation(Double_t radians);
86 
87  void SetAxisX(const TVector3& axis);
88 
89  // Getters
91  inline Int_t GetID() const { return fId; }
92 
94  inline TVector3 GetPosition() const { return fPosition; }
95 
97  Double_t GetRotation() const { return fRotation; }
98 
100  inline TVector3 GetCathodePosition() const { return fPosition + fHeight * fNormal; }
101 
103  inline TVector3 GetNormal() const { return fNormal; }
104 
106  inline TVector3 GetAxisX() const { return fAxisX; }
107 
109  inline TVector3 GetAxisY() const { return fAxisY; }
110 
112  inline Double_t GetChargeCollection() const { return fChargeCollection; }
113 
115  inline Double_t GetDriftDistance() const { return fHeight; }
116 
118  inline Double_t GetHeight() const { return fHeight; }
119 
121  Double_t GetDistanceTo(const TVector3& pos) const;
122 
124  bool IsInside(const TVector3& point) const;
125 
128  TVector2 GetDistanceToModule(Int_t mod, const TVector2& position) {
129  return GetModuleByID(mod)->GetDistanceToModule(position);
130  }
131 
132  TRestDetectorReadoutModule& operator[](int mod) { return fReadoutModules[mod]; }
133 
136  if (mod >= GetNumberOfModules()) {
137  return nullptr;
138  }
139  return &fReadoutModules[mod];
140  }
141 
143  size_t GetNumberOfModules() const { return fReadoutModules.size(); }
144 
146  void AddModule(const TRestDetectorReadoutModule& module);
147 
149  void PrintMetadata() { Print(); }
150 
151  Int_t GetNumberOfChannels();
152 
154 
155  std::string GetType() const { return fType; }
156 
157  std::string GetName() const { return fName; }
158 
159  void SetType(const std::string& type) { fType = type; }
160 
161  void SetName(const std::string& name) { fName = name; }
162 
163  Int_t isZInsideDriftVolume(Double_t z);
164 
165  Int_t isZInsideDriftVolume(const TVector3& position);
166 
167  Bool_t isDaqIDInside(Int_t daqId);
168 
169  Int_t GetModuleIDFromPosition(const TVector3& position) const;
170 
171  TVector2 GetPositionInPlane(const TVector3& point) const;
172  Double_t GetDistanceToPlane(const TVector3& point) const;
173 
174  TVector3 GetPositionInWorld(const TVector2& point, Double_t height = 0) const;
175 
176  void Draw();
177 
178  void Print(Int_t DetailLevel = 0);
179 
180  Int_t FindChannel(Int_t module, const TVector2& position);
181 
182  Double_t GetX(Int_t modID, Int_t chID);
183  Double_t GetY(Int_t modID, Int_t chID);
184 
185  TH2Poly* GetReadoutHistogram();
186  void GetBoundaries(double& xmin, double& xmax, double& ymin, double& ymax);
187 
188  // Constructor
190  // Destructor
192 
193  ClassDef(TRestDetectorReadoutPlane, 6);
194 };
195 #endif
TVector2 GetDistanceToModule(const TVector2 &position)
Creates a TVector2 with shortest norm going from the given position pos to the module....
Bool_t isDaqIDInside(Int_t daqId)
This method determines if the daqId given is associated to any of the readout readout channels in any...
Int_t GetNumberOfChannels()
Returns the total number of channels in the readout plane.
Int_t isZInsideDriftVolume(Double_t z)
This method determines if a given position in z is inside the drift volume drifting distance for this...
TRestDetectorReadoutModule * GetModuleByID(Int_t modID)
Returns a pointer to a module using its internal module id.
void SetPosition(const TVector3 &position)
Sets the readout plane position.
void SetNormal(const TVector3 &normal)
It updates the value of the normal vector and recalculates the corresponding X and Y axis.
TH2Poly * GetReadoutHistogram()
Creates and returns a TH2Poly object with the readout pixel description.
TVector3 GetNormal() const
Returns a TVector3 with a vector normal to the readout plane.
Double_t fRotation
Rotation (in radians) of the readout plane around the normal vector.
Double_t GetHeight() const
Returns the height of the readout volume.
Double_t GetChargeCollection() const
Returns the charge collection ratio at this readout plane.
TVector3 GetAxisY() const
Returns a TVector3 with a vector that defines the Y-axis plane coordinate system.
TVector3 GetPosition() const
Returns a TVector3 with the readout plane position.
Int_t GetModuleIDFromPosition(const TVector3 &position) const
This method returns the module id where pos is found. The z-coordinate must be found in between the c...
Int_t fId
The plane id number imposed by the order of creation. Being the first id=0.
Double_t fChargeCollection
The fraction of charge/energy this readout plane collects from a hit position.
void Print(Int_t DetailLevel=0)
Prints information with details of the readout plane and modules defined inside the readout plane.
void SetChargeCollection(Double_t charge)
Sets the value for the charge collection.
Double_t GetRotation() const
Returns the rotation angle in radians of the reference frame with respect to the normal vector.
TRestDetectorReadoutPlane()
Default TRestDetectorReadoutPlane constructor.
bool IsInside(const TVector3 &point) const
Check if the point is inside any module of the readout plane.
void PrintMetadata()
Prints the readout plane description.
size_t GetNumberOfModules() const
Returns the total number of modules in the readout plane.
void SetID(int id)
Sets the planeId. This is done by TRestDetectorReadout during initialization.
TVector3 GetAxisX() const
Returns a TVector3 with a vector that defines the X-axis plane coordinate system.
void Draw()
Draws the readout plane using GetReadoutHistogram.
Double_t GetX(Int_t modID, Int_t chID)
Returns the X coordinate of a given channel in a given module using their internal module and channel...
Double_t GetDistanceTo(const TVector3 &pos) const
Returns the perpendicular distance to the readout plane from a given position pos.
Int_t FindChannel(Int_t module, const TVector2 &position)
Finds the readout channel index for a given module stored in a given module index stored in the reado...
TRestDetectorReadoutModule * GetModule(size_t mod)
Returns a pointer to a readout module using its std::vector index.
Double_t fHeight
A length in mm that confers a 3rd dimension to the readout plane and defines a volume.
void GetBoundaries(double &xmin, double &xmax, double &ymin, double &ymax)
Finds the xy boundaries of the readout plane delimited by the modules.
Double_t GetY(Int_t modID, Int_t chID)
Returns the Y coordinate of a given channel in a given module using their internal module and channel...
virtual ~TRestDetectorReadoutPlane()
Default TRestDetectorReadoutPlane destructor.
TVector3 fPosition
The position of the readout plane that defines the origin (0,0) in plane coordinates.
TVector3 GetCathodePosition() const
Returns a TVector3 with the cathode position.
void SetHeight(Double_t d)
Used to define the height of the readout volume with sign crosscheck.
TVector3 fAxisX
A vector contained in the plane that defines the plane X-axis.
TVector3 fAxisY
A vector contained in the plane that defines the plane Y-axis.
TVector2 GetDistanceToModule(Int_t mod, const TVector2 &position)
TVector3 fNormal
A vector that defines the plane orientation and the side of the active volume.
Int_t GetID() const
Returns an integer with the plane id number.
Double_t GetDriftDistance() const
Returns the total drift distance. Equivalent to the height of the readout volume.
std::vector< TRestDetectorReadoutModule > fReadoutModules
< A list of TRestDetectorReadoutModule components contained in the readout plane.
void AddModule(const TRestDetectorReadoutModule &module)
Adds a new module to the readout plane.