Project
Loading...
Searching...
No Matches
SimpleEventDisplay.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
14
15#include "TH1D.h"
16#include "TH2S.h"
17#include "TROOT.h"
18#include "TString.h"
19#include "TH2Poly.h"
20
21#include "DataFormatsTPC/Defs.h"
22#include "TPCBase/CalArray.h"
23#include "TPCBase/PadPos.h"
24#include "TPCBase/CRU.h"
25#include "TPCBase/Mapper.h"
26
28
29using namespace o2::tpc;
30
32 : CalibRawBase(),
33 mPadMax("qMax", PadSubset::ROC),
34 mPadOccupancy("occupancy", PadSubset::ROC),
35 mHSigIROC(nullptr),
36 mHSigOROC(nullptr),
37 mPedestals(nullptr),
38 mCurrentChannel(-1),
39 mCurrentROC(-1),
40 mLastSector(-1),
41 mSelectedSector(-1),
42 mLastSelSector(-1),
43 mCurrentRow(-1),
44 mCurrentPad(-1),
45 mMaxPadSignal(-1),
46 mMaxTimeBin(-1),
47 mSectorLoop(kFALSE),
48 mFirstTimeBin(0),
49 mLastTimeBin(512),
50 mTPCmapper(Mapper::instance()),
51 mSignalThreshold(0),
52 mShowOccupancy(kFALSE)
53{
54 initHistograms();
55}
56
57void SimpleEventDisplay::initHistograms()
58{
59 delete mHSigOROC;
60 delete mHSigIROC;
61
62 const int nPadsIROC = mTPCmapper.getPadsInIROC();
63 const int nPadsOROC = mTPCmapper.getPadsInOROC();
64 const int numberOfTimeBins = mLastTimeBin - mFirstTimeBin;
65
66 mHSigIROC = new TH2D("PadSigIROC", "Pad Signals IROC", numberOfTimeBins, mFirstTimeBin, mLastTimeBin, nPadsIROC, 0, nPadsIROC);
67 mHSigOROC = new TH2D("PadSigOROC", "Pad Signals OROC", numberOfTimeBins, mFirstTimeBin, mLastTimeBin, nPadsOROC, 0, nPadsOROC);
68}
69
70//_____________________________________________________________________
72 const Int_t row,
73 const Int_t pad,
74 const Int_t timeBin,
75 const Float_t signal)
76{
77 //
78 // Signal filling methode on the fly pedestal and time offset correction if necessary.
79 // no extra analysis necessary. Assumes knowledge of the signal shape!
80 // assumes that it is looped over consecutive time bins of one pad
81 //
82 // printf("update called: %d, %d, %d, %d, %.3f\n", roc, row, pad, timeBin, signal);
83 if (row < 0) {
84 return 0;
85 }
86 if (pad < 0) {
87 return 0;
88 }
89 if (timeBin < 0) {
90 return 0;
91 }
92 if ((timeBin > mLastTimeBin) || (timeBin < mFirstTimeBin)) {
93 return 0;
94 }
95 if (mSectorLoop && roc % 36 != mSelectedSector % 36) {
96 return 0;
97 }
98
99 if (row < 0 || pad < 0) {
100 printf("Wrong Pad or Row number, skipping!");
101 return 0;
102 }
103
104 float corrSignal = signal;
105
106 // ===| get pedestal |========================================================
107 if (mPedestals) {
108 corrSignal -= mPedestals->getValue(ROC(roc), row, pad);
109 }
110
111 const int iChannel = mTPCmapper.getPadNumberInROC(PadROCPos(roc, row, pad));
112
113 // init first pad and roc in this event
114 if (mCurrentChannel == -1) {
115 mCurrentChannel = iChannel;
116 mCurrentROC = roc;
117 mCurrentRow = row;
118 mCurrentPad = pad;
119 }
120
121 // process last pad if we change to a new one
122 if (iChannel != mCurrentChannel) {
123 mLastSector = mCurrentROC;
124 mCurrentChannel = iChannel;
125 mCurrentROC = roc;
126 mCurrentRow = row;
127 mCurrentPad = pad;
128 mMaxPadSignal = 0;
129 }
130
131 // fill signals for current pad
132 if (mCurrentROC % 36 == mSelectedSector % 36) {
133 const Int_t nbins = mLastTimeBin - mFirstTimeBin;
134 const Int_t offset = (nbins + 2) * (iChannel + 1) + (timeBin - mFirstTimeBin) + 1;
135
136 if ((UInt_t)roc < mTPCmapper.getNumberOfIROCs()) {
137 mHSigIROC->GetArray()[offset] = corrSignal >= mSignalThreshold ? corrSignal : 0;
138 } else {
139 mHSigOROC->GetArray()[offset] = corrSignal >= mSignalThreshold ? corrSignal : 0;
140 }
141 }
142
143 CalROC& calROC = mPadMax.getCalArray(mCurrentROC);
144 auto val = calROC.getValue(row, pad);
145
146 if (corrSignal > val && corrSignal >= mSignalThreshold) {
147 calROC.setValue(row, pad, corrSignal);
148 mMaxPadSignal = corrSignal;
149 mMaxTimeBin = timeBin;
150 }
151
152 CalROC& calROCOccupancy = mPadOccupancy.getCalArray(mCurrentROC);
153 const auto occupancy = calROCOccupancy.getValue(row, pad);
154
155 if (corrSignal >= mSignalThreshold) {
156 calROCOccupancy.setValue(row, pad, occupancy + 1.0f);
157 }
158
159 return 0;
160}
161
162//_____________________________________________________________________
164{
165 if (mSelectedSector % 36 != mLastSelSector % 36) {
166 mSectorLoop = kTRUE;
168 mLastSelSector = mSelectedSector;
169 mSectorLoop = kFALSE;
170 }
171}
172
173//_____________________________________________________________________
174TH1D* SimpleEventDisplay::makePadSignals(Int_t roc, Int_t row, Int_t pad)
175{
176 const int padOffset = (roc > 35) * Mapper::getPadsInIROC();
177 const int channel = mTPCmapper.getPadNumberInROC(PadROCPos(roc, row, pad));
178
179 const FECInfo& fecInfo = mTPCmapper.getFECInfo(PadROCPos(roc, row, pad));
180 const int cruNumber = mTPCmapper.getCRU(ROC(roc).getSector(), channel + padOffset);
181 const CRU cru(cruNumber);
182 const PartitionInfo& partInfo = mTPCmapper.getMapPartitionInfo()[cru.partition()];
183 const int nFECs = partInfo.getNumberOfFECs();
184 const int fecOffset = (nFECs + 1) / 2;
185 const int fecInPartition = fecInfo.getIndex() - partInfo.getSectorFECOffset();
186 const int dataWrapperID = fecInPartition >= fecOffset;
187 const int globalLinkID = (fecInPartition % fecOffset) + dataWrapperID * 12;
188
189 mSelectedSector = roc;
190
191 // attention change for if event has changed
193
194 const Int_t nbins = mLastTimeBin - mFirstTimeBin;
195 if (nbins <= 0) {
196 return nullptr;
197 }
198
199 TH2D* hPadSignals = nullptr;
200
201 // ===| ADC vs. Pad vs. Time for row |========================================
202 TH2F* hPadTime = nullptr;
203 if (roc < (Int_t)mTPCmapper.getNumberOfIROCs()) {
204 hPadTime = static_cast<TH2F*>(gROOT->FindObject("hPadTimeValsI"));
205 hPadSignals = mHSigIROC;
206 } else {
207 hPadTime = static_cast<TH2F*>(gROOT->FindObject("hPadTimeValsO"));
208 hPadSignals = mHSigOROC;
209 }
210
211 static Int_t lastRoc = -1;
212 static Int_t lastRow = -1;
213 if (hPadTime && ((lastRoc != roc) || (lastRow != row))) {
214 hPadTime->Reset();
215 const auto nPads = mTPCmapper.getNumberOfPadsInRowROC(roc, row);
216 const auto nBins = hPadTime->GetNbinsY();
217 const auto shift = nBins / 2 - nPads / 2;
218 for (int iPad = 0; iPad < nPads; ++iPad) {
219 const int ichannel = mTPCmapper.getPadNumberInROC(PadROCPos(roc, row, iPad));
220 const Int_t offset = (nbins + 2) * (ichannel + 1);
221 const double* arrSig = hPadSignals->GetArray() + offset;
222 for (int iTime = 0; iTime < nbins; ++iTime) {
223 hPadTime->SetBinContent(iTime + 1, iPad + shift + 1, arrSig[iTime + 1]);
224 }
225 }
226 hPadTime->SetEntries(nPads * nbins);
227 hPadTime->SetTitle(fmt::format("Pad row {}", row).data());
228 hPadTime->SetUniqueID(row);
229 }
230
231 // ===| ADC vs. time for single pad |=========================================
232 TH1D* h = nullptr;
233 const Int_t offset = (nbins + 2) * (channel + 1);
234 Double_t* arrP = nullptr;
235
236 TString title("#splitline{#lower[.1]{#scale[.5]{");
237
238 if (roc < (Int_t)mTPCmapper.getNumberOfIROCs()) {
239 h = (TH1D*)gROOT->FindObject("PadSignals_IROC");
240 if (!h) {
241 h = new TH1D("PadSignals_IROC", "PadSignals IROC;time bin (200ns);amplitude (ADC counts)", nbins, mFirstTimeBin, mLastTimeBin);
242 }
243 h->SetFillColor(kBlue - 10);
244 arrP = mHSigIROC->GetArray() + offset;
245 title += "IROC ";
246 } else {
247 h = (TH1D*)gROOT->FindObject("PadSignals_OROC");
248 if (!h) {
249 h = new TH1D("PadSignals_OROC", "PadSignals OROC;time bin (200ns);amplitude (ADC counts)", nbins, mFirstTimeBin, mLastTimeBin);
250 }
251 h->SetFillColor(kBlue - 10);
252 arrP = mHSigOROC->GetArray() + offset;
253 title += "OROC ";
254 }
255 title += (roc / 18 % 2 == 0) ? "A" : "C";
256 title += Form("%02d (%02d) row: %02d, pad: %03d, globalpad: %05d (in roc)}}}{#scale[.5]{FEC: %02d (%02d), Chip: %02d, Chn: %02d, CRU: %d, Link: %02d (%s%d)}}",
257 roc % 18, roc, row, pad, channel, fecInfo.getIndex(), fecInPartition, fecInfo.getSampaChip(), fecInfo.getSampaChannel(), cruNumber % CRU::CRUperSector, globalLinkID, dataWrapperID ? "B" : "A", globalLinkID % 12);
258
259 h->SetTitle(title.Data());
260 Int_t entries = 0;
261 for (Int_t i = 0; i < nbins; i++) {
262 entries += (Int_t)arrP[i + 1];
263 }
264 h->Set(nbins + 2, arrP);
265 h->SetEntries(entries);
266 return h;
267}
268
269//_____________________________________________________________________
270void SimpleEventDisplay::resetEvent()
271{
272 //
273 //
274 //
275 if (!mSectorLoop) {
276 mPadMax.multiply(0.);
277 mPadOccupancy.multiply(0.);
278 }
279 mHSigIROC->Reset();
280 mHSigOROC->Reset();
281}
282
283//______________________________________________________________________________
285{
286 if (!h) {
287 return;
288 }
289 if (timeBin < mFirstTimeBin || timeBin > mLastTimeBin) {
290 return;
291 }
292
293 int ichannel = 0;
294 const int iTimeBin = timeBin - mFirstTimeBin + 1;
295 // IROC loop
296 for (int ipad = 0; ipad < mHSigIROC->GetNbinsY(); ++ipad, ++ichannel) {
297 h->SetBinContent(ichannel + 1, mHSigIROC->GetBinContent(iTimeBin, ipad + 1));
298 }
299
300 // OROC loop
301 for (int ipad = 0; ipad < mHSigOROC->GetNbinsY(); ++ipad, ++ichannel) {
302 h->SetBinContent(ichannel + 1, mHSigOROC->GetBinContent(iTimeBin, ipad + 1));
303 }
304
305 h->SetEntries(ichannel);
306}
int32_t i
uint32_t roc
Definition RawData.h:3
Class for time synchronization of RawReader instances.
@ CRUperSector
Definition CRU.h:30
unsigned char partition() const
Definition CRU.h:63
void setValue(const size_t channel, const T &value)
Definition CalArray.h:95
const T getValue(const size_t channel) const
Definition CalArray.h:96
const CalType & getCalArray(size_t position) const
Definition CalDet.h:63
const CalDet & multiply(const T &val)
Definition CalDet.h:87
const T getValue(const int sec, const int globalPadInSector) const
Definition CalDet.h:154
Base class for raw data calibrations.
size_t getPresentEventNumber() const
get present event number
ProcessStatus processEvent(int eventNumber=-1)
unsigned char getSampaChannel() const
Definition FECInfo.h:45
unsigned char getIndex() const
Definition FECInfo.h:41
unsigned char getSampaChip() const
Definition FECInfo.h:44
const FECInfo & getFECInfo(const PadROCPos &padROC) const
Definition Mapper.h:179
int getNumberOfPadsInRowROC(int roc, int row) const
Definition Mapper.h:342
static constexpr unsigned short getPadsInIROC()
Definition Mapper.h:409
const std::array< PartitionInfo, 5 > & getMapPartitionInfo() const
Definition Mapper.h:390
static constexpr unsigned short getPadsInOROC()
Definition Mapper.h:413
GlobalPadNumber getPadNumberInROC(const PadROCPos &rocPadPosition) const
Definition Mapper.h:96
static constexpr unsigned short getNumberOfIROCs()
Definition Mapper.h:407
int getCRU(const Sector &sec, GlobalPadNumber globalPad) const
Definition Mapper.h:78
Pad and row inside a ROC.
Definition PadROCPos.h:37
unsigned char getNumberOfFECs() const
unsigned char getSectorFECOffset() const
Sector getSector() const
get sector
Definition ROC.h:105
TH1D * makePadSignals(Int_t roc, Int_t row, Int_t pad)
void fillSectorHistSingleTimeBin(TH2Poly *h, Int_t timeBin)
Int_t updateROC(const Int_t roc, const Int_t row, const Int_t pad, const Int_t timeBin, const Float_t signal) final
GLintptr offset
Definition glcorearb.h:660
GLuint GLfloat * val
Definition glcorearb.h:1582
Global TPC definitions and constants.
Definition SimTraits.h:167
PadSubset
Definition of the different pad subsets.
Definition Defs.h:63
std::vector< int > row