REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorReadoutEventViewer.cxx
1
14
15#include "TRestDetectorReadoutEventViewer.h"
16
17using namespace std;
18
19Int_t planeId = 0;
20
22
23 TRestDetectorReadoutEventViewer::TRestDetectorReadoutEventViewer() {
24 Initialize();
25}
26
27TRestDetectorReadoutEventViewer::~TRestDetectorReadoutEventViewer() {}
28
29void TRestDetectorReadoutEventViewer::Initialize() {
30 TRestEventViewer::Initialize();
31
32 fCanvasXY = new TCanvas("ReadoutMap", "ReadoutMap");
33 fCanvasXZYZ = new TCanvas("XZYZ", "XZYZ");
34 fCanvasXZYZ->Divide(2, 1);
35
36 fHistoXZ = nullptr;
37 fHistoYZ = nullptr;
38
39 fSignalEvent = new TRestDetectorSignalEvent();
40 SetEvent(fSignalEvent);
41}
42
43void TRestDetectorReadoutEventViewer::SetReadout(TRestDetectorReadout* readout) {
44 // Finalize the instantiation based on argument TRestDetectorReadout
45 fReadout = readout;
46 cout << "WARNING : Only plane 0 is drawn. Implementation to draw several "
47 "planes or to choose the plane must be implemented."
48 << endl;
49 fReadout->PrintMetadata();
50 TRestDetectorReadoutPlane* plane = &(*fReadout)[planeId];
51 fHistoXY = plane->GetReadoutHistogram();
52 plane->GetBoundaries(xmin, xmax, ymin, ymax);
53}
54
55void TRestDetectorReadoutEventViewer::AddEvent(TRestEvent* ev) {
56 // Finalize the drawing of current event, adding to the per-channel-signal
57 // vs. time drawn by the generic event viewer, the three 2D histograms of
58 // the XY, XZ and YZ projections.
59 TRestEventViewer::AddEvent(ev);
60
61 if (fPad == nullptr) return;
62
63 fSignalEvent = (TRestDetectorSignalEvent*)ev;
64
65 // XY histo is expected to have always same binning => reset it.
66 // (X|Y)Z may change from event to event => delete (and later on re-create).
67 fHistoXY->Reset(0);
68 if (fHistoXZ != nullptr) {
69 delete fHistoXZ;
70 fHistoXZ = nullptr;
71 }
72 if (fHistoYZ != nullptr) {
73 delete fHistoYZ;
74 fHistoYZ = nullptr;
75 }
76
77 DrawReadoutPulses();
78 fCanvasXY->cd();
79 fHistoXY->Draw("COLZ");
80 fCanvasXY->Update();
81
82 fCanvasXZYZ->cd(1);
83 fHistoXZ->Draw("CONTZ");
84 fCanvasXZYZ->cd(2);
85 fHistoYZ->Draw("CONTZ");
86 fCanvasXZYZ->Update();
87}
88
89void TRestDetectorReadoutEventViewer::DrawReadoutPulses() {
90 int readoutChannel, daqChannel;
91 double charge;
92
93 Int_t modId;
94 TRestDetectorReadoutModule* module = nullptr;
96
97 int maxIndex;
98
99 zmin = std::numeric_limits<Double_t>::max();
100 zmax = 1E-9;
101
102 double z;
103 for (int i = 0; i < fSignalEvent->GetNumberOfSignals(); i++)
104 for (int j = 0; j < fSignalEvent->GetSignal(i)->GetNumberOfPoints(); j++) {
105 z = fSignalEvent->GetSignal(i)->GetTime(j);
106 if (z < zmin) zmin = z;
107 if (z > zmax) zmax = z;
108 }
109
110 zmax++;
111 zmin--;
112
113 fHistoXZ = new TH2D("XZ", "XZ", 100, xmin, xmax, 100, zmin, zmax);
114 fHistoYZ = new TH2D("YZ", "YZ", 100, ymin, ymax, 100, zmin, zmax);
115
116 for (int i = 0; i < fSignalEvent->GetNumberOfSignals(); i++) {
117 daqChannel = fSignalEvent->GetSignal(i)->GetSignalID();
118
119 TRestDetectorReadoutPlane* plane = &(*fReadout)[planeId];
120 for (size_t m = 0; m < plane->GetNumberOfModules(); m++) {
121 module = &(*plane)[m];
122
123 if (module->IsDaqIDInside(daqChannel)) break;
124 }
125 modId = module->GetModuleID();
126
127 readoutChannel = module->DaqToReadoutChannel(daqChannel);
128 cout << "daqChannel " << daqChannel << " readoutChannel " << readoutChannel << endl;
129
130 if ((channel = GetChannel(readoutChannel)) == nullptr) {
131 continue;
132 }
133
134 int nPixels = channel->GetNumberOfPixels();
135
136 Double_t xRead = plane->GetX(modId, readoutChannel);
137 Double_t yRead = plane->GetY(modId, readoutChannel);
138
139 // Pixel readout
140 Int_t xStrip = 0;
141 Int_t yStrip = 0;
142 if (TMath::IsNaN(yRead)) xStrip = 1;
143 if (TMath::IsNaN(xRead)) yStrip = 1;
144 if (xStrip == 0 && yStrip == 0) {
145 for (int j = 0; j < fSignalEvent->GetSignal(i)->GetNumberOfPoints(); j++)
146 fHistoXZ->Fill(xRead, fSignalEvent->GetSignal(i)->GetTime(j),
147 fSignalEvent->GetSignal(i)->GetData(j));
148 for (int j = 0; j < fSignalEvent->GetSignal(i)->GetNumberOfPoints(); j++)
149 fHistoYZ->Fill(yRead, fSignalEvent->GetSignal(i)->GetTime(j),
150 fSignalEvent->GetSignal(i)->GetData(j));
151 }
152 // Strips readout
153 else {
154 if (yStrip == 0) {
155 for (int j = 0; j < fSignalEvent->GetSignal(i)->GetNumberOfPoints(); j++)
156 fHistoXZ->Fill(xRead, fSignalEvent->GetSignal(i)->GetTime(j),
157 fSignalEvent->GetSignal(i)->GetData(j));
158 } else if (xStrip == 0) {
159 for (int j = 0; j < fSignalEvent->GetSignal(i)->GetNumberOfPoints(); j++) {
160 fHistoYZ->Fill(yRead, fSignalEvent->GetSignal(i)->GetTime(j),
161 fSignalEvent->GetSignal(i)->GetData(j));
162 }
163 }
164 }
165
166 // Only the maximum is plotted for the XY readout
167 maxIndex = fSignalEvent->GetSignal(i)->GetMaxIndex();
168 charge = fSignalEvent->GetSignal(i)->GetData(maxIndex);
169
170 for (int px = 0; px < nPixels; px++) {
171 // cout<< "X : " << module->GetPixelCenter( readoutChannel, px ).X()<<" Y
172 // : "<<module->GetPixelCenter( readoutChannel, px ).Y()<<" Ch:
173 // "<<charge<<endl;
174 fHistoXY->Fill(module->GetPixelCenter(readoutChannel, px).X(),
175 module->GetPixelCenter(readoutChannel, px).Y(), charge);
176 }
177 }
178}
179
180TRestDetectorReadoutChannel* TRestDetectorReadoutEventViewer::GetChannel(int readoutChannel) {
181 TRestDetectorReadoutPlane* plane = &(*fReadout)[0];
182 for (size_t n = 0; n < plane->GetNumberOfModules(); n++) {
183 if ((*plane)[n].GetChannel(readoutChannel) == nullptr) continue;
184 return (*plane)[n].GetChannel(readoutChannel);
185 }
186
187 cout << "Readout channel " << readoutChannel << " not found" << endl;
188 return nullptr;
189}
190
191TRestDetectorReadoutModule* TRestDetectorReadoutEventViewer::GetModule(int readoutChannel) {
192 TRestDetectorReadoutPlane* plane = &(*fReadout)[0];
193 for (int n = 0; n < fReadout->GetNumberOfModules(); n++) {
194 if ((*plane)[n].GetChannel(readoutChannel) == nullptr) continue;
195 return &(*plane)[n];
196 }
197
198 cout << "Readout channel " << readoutChannel << " not found" << endl;
199 return nullptr;
200}
Int_t GetNumberOfPixels()
Returns the total number of pixels inside the readout channel.
Int_t GetModuleID() const
Returns the module id.
TVector2 GetPixelCenter(Int_t channel, Int_t pixel)
Returns the center pixel position for a given channel and pixel indexes.
Bool_t IsDaqIDInside(Int_t daqID)
Determines if a given daqID number is in the range of the module.
Int_t DaqToReadoutChannel(Int_t daqChannel)
Returns the physical readout channel index for a given daq id channel number.
TH2Poly * GetReadoutHistogram()
Creates and returns a TH2Poly object with the readout pixel description.
size_t GetNumberOfModules() const
Returns the total number of modules in the readout plane.
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...
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...
A metadata class to generate/store a readout description.
void PrintMetadata() override
Implemented it in the derived metadata class to print out specific metadata information.
Int_t GetNumberOfModules()
Returns the total number of modules implemented in all the readout planes.
A base class for any REST event.
Definition: TRestEvent.h:38