REST-for-Physics
v2.3
Rare Event Searches ToolKit for Physics
|
A class to store the readout module definition used in TRestDetectorReadoutPlane. It allows to integrate any number of independent readout channels.
This class stores the readout module geometrical description, module position, orientation, and size. It contains a vector of TRestDetectorReadoutChannel with the definition of the readout channels existing in the readout module.
RESTsoft - Software for Rare Event Searches with TPCs
History of developments:
2015-aug: First concept. Javier Galan
Definition at line 36 of file TRestDetectorReadoutModule.h.
#include <TRestDetectorReadoutModule.h>
Public Member Functions | |
void | AddChannel (TRestDetectorReadoutChannel &channel) |
Adds a new channel to the module. | |
Int_t | DaqToReadoutChannel (Int_t daqChannel) |
Returns the physical readout channel index for a given daq id channel number. | |
void | DisableWarnings () |
Disables warning output. | |
void | DoReadoutMapping () |
Starts the readout mapping initialization. This process is computationally expensive but it greatly optimizes the FindChannel process later on. | |
void | Draw () |
Not implemented. | |
void | EnableWarnings () |
Enables warning output. | |
Int_t | FindChannel (const TVector2 &position) |
Returns the channel index corresponding to the absolute coordinates (absX, absY), but relative to the readout plane coordinate system. More... | |
TRestDetectorReadoutChannel * | GetChannel (size_t n) |
Returns a pointer to a readout channel by index. | |
TVector2 | GetDistanceToModule (const TVector2 &position) |
Creates a TVector2 with shortest norm going from the given position pos to the module. It can be seen as the vector to add to move from the position to the closest border of the module. | |
TRestDetectorReadoutMapping * | GetMapping () |
Returns a pointer to the readout mapping. | |
Int_t | GetMappingNodes () const |
Int_t | GetMaxDaqID () const |
Returns the maximum daq id number. | |
Int_t | GetMinDaqID () const |
Returns the minimum daq id number. | |
TVector2 | GetModuleCoordinates (const TVector2 &p) const |
Int_t | GetModuleID () const |
Returns the module id. | |
std::string | GetName () const |
Returns the module name. | |
size_t | GetNumberOfChannels () const |
Returns the total number of channels defined inside the module. | |
TVector2 | GetOrigin () const |
Returns the module origin position. | |
TVector2 | GetPixelCenter (Int_t channel, Int_t pixel) |
Returns the center pixel position for a given channel and pixel indexes. More... | |
TVector2 | GetPixelCenter (TRestDetectorReadoutPixel *pix) |
TVector2 | GetPixelOrigin (Int_t channel, Int_t pixel) |
Returns the pixel origin (left-bottom) position for a given channel and pixel indexes. | |
TVector2 | GetPixelOrigin (TRestDetectorReadoutPixel *pix) |
Bool_t | GetPixelTriangle (Int_t channel, Int_t pixel) |
Returns the pixel type for a given channel and pixel indexes. More... | |
Bool_t | GetPixelTriangle (TRestDetectorReadoutPixel *pix) |
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. More... | |
TVector2 | GetPixelVertex (TRestDetectorReadoutPixel *pix, Int_t vertex) |
TVector2 | GetPlaneCoordinates (const TVector2 &p) |
Double_t | GetRotation () const |
Returns the module rotation in degrees. | |
TVector2 | GetSize () const |
Returns the module size (x, y) in mm. | |
Double_t | GetTolerance () const |
Gets the tolerance for independent pixel overlaps. | |
std::string | GetType () const |
TVector2 | GetVertex (int n) const |
Returns the coordinates of the specified vertex index n. The physical coordinates relative to the readout plane are returned, including rotation. More... | |
Bool_t | IsDaqIDInside (Int_t daqID) |
Determines if a given daqID number is in the range of the module. | |
Bool_t | IsInside (const TVector2 &position) const |
Determines if the position x,y relative to the readout plane are inside this readout module. More... | |
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 given. | |
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 given. | |
TRestDetectorReadoutChannel & | operator[] (int n) |
void | Print (Int_t DetailLevel=0) |
Prints the module details and channels if fullDetail is enabled. | |
void | SetDecodingFile (const std::string &decodingFile) |
Set the decoding file in the readout module. | |
void | SetFirstDaqChannel (Int_t firstDaqChannel) |
Sets first DAQ channel. | |
void | SetMappingNodes (Int_t nodes) |
Sets number of nodes. | |
void | SetMinMaxDaqIDs () |
Initializes the max and min values for the daq channel number. | |
void | SetModuleID (Int_t modID) |
Sets the module by id definition. | |
void | SetName (const std::string &name) |
Sets the name of the readout module. | |
void | SetOrigin (const TVector2 &origin) |
Sets the module origin by definition using TVector2 input. | |
void | SetRotation (Double_t rotation) |
Sets the module rotation in degrees. | |
void | SetSize (const TVector2 &size) |
Sets the module size by definition using TVector2 input. | |
void | SetTolerance (Double_t tolerance) |
Sets the tolerance for independent pixel overlaps. | |
void | SetType (const std::string &type) |
Sets the type of the readout module. | |
TRestDetectorReadoutModule () | |
Default TRestDetectorReadoutModule constructor. | |
virtual | ~TRestDetectorReadoutModule () |
Default TRestDetectorReadoutModule destructor. | |
Private Member Functions | |
void | Initialize () |
TRestDetectorReadoutModule initialization. | |
TVector2 | TransformToModuleCoordinates (const TVector2 &coords) const |
TVector2 | TransformToPlaneCoordinates (const TVector2 &coords) const |
Private Attributes | |
std::pair< Int_t, Int_t > | fDaqIdRange |
The minimum and maximum daq channel ids associated to the module. More... | |
Bool_t | fDecoding |
std::string | fDecodingFile = "" |
Decoding file. | |
*default REST_Warning in TRestDetectorReadout will enable it *Int_t | fFirstDaqChannel = 0 |
First DAQ channel. | |
Int_t | fId = -1 |
The module id given by the readout definition. | |
TRestDetectorReadoutMapping | fMapping |
The readout module uniform grid mapping. | |
Int_t | fMappingNodes = 0 |
Number of nodes. | |
std::string | fName |
TVector2 | fOrigin = {0, 0} |
The module (x, y) position relative to the readout plane position. | |
std::vector< TRestDetectorReadoutChannel > | fReadoutChannel |
Double_t | fRotation = 0 |
The rotation of the module around the module origin (fModuleOriginX, fModuleOriginY) in radians. | |
TVector2 | fSize = {0, 0} |
The module (x, y) size. All pixels should be contained within this size. | |
Double_t | fTolerance |
std::string | fType |
Bool_t | showWarnings |
Int_t TRestDetectorReadoutModule::FindChannel | ( | const TVector2 & | position | ) |
Returns the channel index corresponding to the absolute coordinates (absX, absY), but relative to the readout plane coordinate system.
The readout mapping (see TRestDetectorReadoutMapping) is used to help finding the pixel where coordinates absX and absY fall in.
Definition at line 322 of file TRestDetectorReadoutModule.cxx.
|
inline |
Converts the coordinates given by TVector2 in the readout plane reference system to the readout module reference system.
Definition at line 154 of file TRestDetectorReadoutModule.h.
TVector2 TRestDetectorReadoutModule::GetPixelCenter | ( | Int_t | channel, |
Int_t | pixel | ||
) |
Returns the center pixel position for a given channel and pixel indexes.
vertex | A value between 0-3 defining the vertex position to be returned |
Definition at line 505 of file TRestDetectorReadoutModule.cxx.
Bool_t TRestDetectorReadoutModule::GetPixelTriangle | ( | Int_t | channel, |
Int_t | pixel | ||
) |
Returns the pixel type for a given channel and pixel indexes.
vertex | A boolean that is true if the pixel is triangular, false otherwise |
Definition at line 520 of file TRestDetectorReadoutModule.cxx.
TVector2 TRestDetectorReadoutModule::GetPixelVertex | ( | Int_t | channel, |
Int_t | pixel, | ||
Int_t | vertex | ||
) |
Returns any of the pixel vertex position for a given channel and pixel indexes.
vertex | A value between 0-3 defining the vertex position to be returned |
Definition at line 490 of file TRestDetectorReadoutModule.cxx.
|
inline |
Converts the coordinates given by TVector2 in the readout module reference system to the readout plane reference system.
Definition at line 158 of file TRestDetectorReadoutModule.h.
TVector2 TRestDetectorReadoutModule::GetVertex | ( | int | n | ) | const |
Returns the coordinates of the specified vertex index n. The physical coordinates relative to the readout plane are returned, including rotation.
n | A value between 0-3 defining the vertex position to be returned |
Definition at line 555 of file TRestDetectorReadoutModule.cxx.
Bool_t TRestDetectorReadoutModule::IsInside | ( | const TVector2 & | position | ) | const |
Determines if the position x,y relative to the readout plane are inside this readout module.
Determines if the position TVector2 pos relative to the readout plane are inside this readout module.
Definition at line 415 of file TRestDetectorReadoutModule.cxx.
|
inlineprivate |
Converts the coordinates (xPhys,yPhys) in the readout plane reference system to the readout module reference system.
Definition at line 78 of file TRestDetectorReadoutModule.h.
|
inlineprivate |
Converts the coordinates (xMod,yMod) in the readout module reference system to the readout plane reference system.
Definition at line 84 of file TRestDetectorReadoutModule.h.
|
private |
The minimum and maximum daq channel ids associated to the module.
Definition at line 47 of file TRestDetectorReadoutModule.h.
|
private |
Defines if a decoding file was used to set the relation between a physical readout channel id and a signal daq id
Definition at line 68 of file TRestDetectorReadoutModule.h.
|
private |
A std::vector of the instances of TRestDetectorReadoutChannel contained in the readout module.
Definition at line 51 of file TRestDetectorReadoutModule.h.
|
private |
Tolerance allowed in overlaps at the pixel boundaries in mm.
Definition at line 56 of file TRestDetectorReadoutModule.h.