Project
Loading...
Searching...
No Matches
Clusterer.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
16
17#include <TTree.h>
18#include "Framework/Logger.h"
22
23#include <algorithm>
24
25#ifdef WITH_OPENMP
26#include <omp.h>
27#endif
28using namespace o2::itsmft;
29
30namespace o2::its3
31{
32void Clusterer::process(int nThreads, PixelReader& reader, CompClusCont* compClus,
33 PatternCont* patterns, ROFRecCont* vecROFRec, MCTruth* labelsCl)
34{
35#ifdef _PERFORM_TIMING_
36 mTimer.Start(false);
37#endif
38 if (nThreads < 1) {
39 nThreads = 1;
40 }
41 auto autoDecode = reader.getDecodeNextAuto();
42 do {
43 if (autoDecode) {
44 reader.setDecodeNextAuto(false); // internally do not autodecode
45 if (!reader.decodeNextTrigger()) {
46 break; // on the fly decoding was requested, but there were no data left
47 }
48 }
49 if (reader.getInteractionRecord().isDummy()) {
50 continue; // No IR info was found
51 }
52 // pre-fetch all non-empty chips of current ROF
53 ChipPixelData* curChipData = nullptr;
54 mFiredChipsPtr.clear();
55 size_t nPix = 0;
56 while ((curChipData = reader.getNextChipData(mChips))) {
57 mFiredChipsPtr.push_back(curChipData);
58 nPix += curChipData->getData().size();
59 }
60
61 auto& rof = vecROFRec->emplace_back(reader.getInteractionRecord(), vecROFRec->size(), compClus->size(), 0); // create new ROF
62
63 uint16_t nFired = mFiredChipsPtr.size();
64 if (!nFired) {
65 if (autoDecode) {
66 continue;
67 }
68 break; // just 1 ROF was asked to be processed
69 }
70 if (nFired < nThreads) {
71 nThreads = nFired;
72 }
73#ifndef WITH_OPENMP
74 nThreads = 1;
75#endif
76 uint16_t chipStep = nThreads > 1 ? (nThreads == 2 ? 20 : 10) : nFired;
77 int dynGrp = std::min(4, std::max(1, nThreads / 2));
78 if (nThreads > mThreads.size()) {
79 int oldSz = mThreads.size();
80 mThreads.resize(nThreads);
81 for (size_t i = oldSz; i < nThreads; i++) {
82 mThreads[i] = std::make_unique<ClustererThread>(this, i);
83 }
84 }
85#ifdef WITH_OPENMP
86#pragma omp parallel for schedule(dynamic, dynGrp) num_threads(nThreads)
87 //>> start of MT region
88 for (uint16_t ic = 0; ic < nFired; ic += chipStep) {
89 auto ith = omp_get_thread_num();
90 if (nThreads > 1) {
91 mThreads[ith]->process(ic, std::min(chipStep, uint16_t(nFired - ic)),
92 &mThreads[ith]->compClusters,
93 patterns ? &mThreads[ith]->patterns : nullptr,
94 labelsCl ? reader.getDigitsMCTruth() : nullptr,
95 labelsCl ? &mThreads[ith]->labels : nullptr, rof);
96 } else { // put directly to the destination
97 mThreads[0]->process(0, nFired, compClus, patterns, labelsCl ? reader.getDigitsMCTruth() : nullptr, labelsCl, rof);
98 }
99 }
100 //<< end of MT region
101#else
102 mThreads[0]->process(0, nFired, compClus, patterns, labelsCl ? reader.getDigitsMCTruth() : nullptr, labelsCl, rof);
103#endif
104 // copy data of all threads but the 1st one to final destination
105 if (nThreads > 1) {
106#ifdef _PERFORM_TIMING_
107 mTimerMerge.Start(false);
108#endif
109 size_t nClTot = 0, nPattTot = 0;
110 int chid = 0;
111 std::vector<size_t> thrStatIdx(nThreads, 0);
112 for (int ith = 0; ith < nThreads; ith++) {
113 std::sort(mThreads[ith]->stats.begin(), mThreads[ith]->stats.end(), [](const ThreadStat& a, const ThreadStat& b) { return a.firstChip < b.firstChip; });
114 nClTot += mThreads[ith]->compClusters.size();
115 nPattTot += mThreads[ith]->patterns.size();
116 }
117 compClus->reserve(nClTot);
118 if (patterns) {
119 patterns->reserve(nPattTot);
120 }
121 while (chid < nFired) {
122 for (int ith = 0; ith < nThreads; ith++) {
123 if (thrStatIdx[ith] >= mThreads[ith]->stats.size()) {
124 continue;
125 }
126 const auto& stat = mThreads[ith]->stats[thrStatIdx[ith]];
127 if (stat.firstChip == chid) {
128 thrStatIdx[ith]++;
129 chid += stat.nChips; // next chip to look
130 const auto clbeg = mThreads[ith]->compClusters.begin() + stat.firstClus;
131 compClus->insert(compClus->end(), clbeg, clbeg + stat.nClus);
132 if (patterns) {
133 const auto ptbeg = mThreads[ith]->patterns.begin() + stat.firstPatt;
134 patterns->insert(patterns->end(), ptbeg, ptbeg + stat.nPatt);
135 }
136 if (labelsCl) {
137 labelsCl->mergeAtBack(mThreads[ith]->labels, stat.firstClus, stat.nClus);
138 }
139 }
140 }
141 }
142 for (int ith = 0; ith < nThreads; ith++) {
143 mThreads[ith]->patterns.clear();
144 mThreads[ith]->compClusters.clear();
145 mThreads[ith]->labels.clear();
146 mThreads[ith]->stats.clear();
147 }
148#ifdef _PERFORM_TIMING_
149 mTimerMerge.Stop();
150#endif
151 } else {
152 mThreads[0]->stats.clear();
153 }
154 rof.setNEntries(compClus->size() - rof.getFirstEntry()); // update
155 } while (autoDecode);
156 reader.setDecodeNextAuto(autoDecode); // restore setting
157#ifdef _PERFORM_TIMING_
158 mTimer.Stop();
159#endif
160}
161
162//__________________________________________________
163void Clusterer::ClustererThread::process(uint16_t chip, uint16_t nChips, CompClusCont* compClusPtr, PatternCont* patternsPtr,
164 const ConstMCTruth* labelsDigPtr, MCTruth* labelsClPtr, const ROFRecord& rofPtr)
165{
166 if (stats.empty() || stats.back().firstChip + stats.back().nChips != chip) { // there is a jump, register new block
167 stats.emplace_back(ThreadStat{chip, 0, uint32_t(compClusPtr->size()), patternsPtr ? uint32_t(patternsPtr->size()) : 0, 0, 0});
168 }
169 for (int ic = 0; ic < nChips; ic++) {
170 auto* curChipData = parent->mFiredChipsPtr[chip + ic];
171 auto chipID = curChipData->getChipID();
172 if (parent->mMaxBCSeparationToMask > 0) { // mask pixels fired from the previous ROF
173 const auto& chipInPrevROF = parent->mChipsOld[chipID];
174 if (std::abs(rofPtr.getBCData().differenceInBC(chipInPrevROF.getInteractionRecord())) < parent->mMaxBCSeparationToMask) {
175 parent->mMaxRowColDiffToMask != 0 ? curChipData->maskFiredInSample(parent->mChipsOld[chipID], parent->mMaxRowColDiffToMask) : curChipData->maskFiredInSample(parent->mChipsOld[chipID]);
176 }
177 }
178 auto validPixID = curChipData->getFirstUnmasked();
179 auto npix = curChipData->getData().size();
180 if (validPixID < npix) { // chip data may have all of its pixels masked!
181 auto valp = validPixID++;
182 if (validPixID == npix) { // special case of a single pixel fired on the chip
183 finishChipSingleHitFast(valp, curChipData, compClusPtr, patternsPtr, labelsDigPtr, labelsClPtr);
184 } else {
185 initChip(curChipData, valp);
186 for (; validPixID < npix; validPixID++) {
187 if (!curChipData->getData()[validPixID].isMasked()) {
188 updateChip(curChipData, validPixID);
189 }
190 }
191 finishChip(curChipData, compClusPtr, patternsPtr, labelsDigPtr, labelsClPtr);
192 }
193 }
194 if (parent->mMaxBCSeparationToMask > 0) { // current chip data will be used in the next ROF to mask overflow pixels
195 parent->mChipsOld[chipID].swap(*curChipData);
196 }
197 }
198 auto& currStat = stats.back();
199 currStat.nChips += nChips;
200 currStat.nClus = compClusPtr->size() - currStat.firstClus;
201 currStat.nPatt = patternsPtr ? (patternsPtr->size() - currStat.firstPatt) : 0;
202}
203
204//__________________________________________________
206 PatternCont* patternsPtr, const ConstMCTruth* labelsDigPtr, MCTruth* labelsClusPtr)
207{
208 const auto& pixData = curChipData->getData();
209 for (size_t i1 = 0; i1 < preClusterHeads.size(); ++i1) {
210 auto ci = preClusterIndices[i1];
211 if (ci < 0) {
212 continue;
213 }
214 BBox bbox(curChipData->getChipID());
215 int nlab = 0;
216 int next = preClusterHeads[i1];
217 pixArrBuff.clear();
218 while (next >= 0) {
219 const auto& pixEntry = pixels[next];
220 const auto pix = pixData[pixEntry.second];
221 pixArrBuff.push_back(pix); // needed for cluster topology
222 bbox.adjust(pix.getRowDirect(), pix.getCol());
223 if (labelsClusPtr) {
224 if (parent->mSquashingDepth) { // the MCtruth for this pixel is stored in chip data: due to squashing we lose contiguity
225 fetchMCLabels(curChipData->getOrderedPixId(pixEntry.second), labelsDigPtr, nlab);
226 } else { // the MCtruth for this pixel is at curChipData->startID+pixEntry.second
227 fetchMCLabels(pixEntry.second + curChipData->getStartID(), labelsDigPtr, nlab);
228 }
229 }
230 next = pixEntry.first;
231 }
232 preClusterIndices[i1] = -1;
233 for (size_t i2 = i1 + 1; i2 < preClusterHeads.size(); ++i2) {
234 if (preClusterIndices[i2] != ci) {
235 continue;
236 }
237 next = preClusterHeads[i2];
238 while (next >= 0) {
239 const auto& pixEntry = pixels[next];
240 const auto pix = pixData[pixEntry.second]; // PixelData
241 pixArrBuff.push_back(pix); // needed for cluster topology
242 bbox.adjust(pix.getRowDirect(), pix.getCol());
243 if (labelsClusPtr != nullptr) {
244 if (parent->mSquashingDepth) { // the MCtruth for this pixel is stored in chip data: due to squashing we lose contiguity
245 fetchMCLabels(curChipData->getOrderedPixId(pixEntry.second), labelsDigPtr, nlab);
246 } else { // the MCtruth for this pixel is at curChipData->startID+pixEntry.second
247 fetchMCLabels(pixEntry.second + curChipData->getStartID(), labelsDigPtr, nlab);
248 }
249 }
250 next = pixEntry.first;
251 }
252 preClusterIndices[i2] = -1;
253 }
254 if (bbox.isAcceptableSize()) {
255 parent->streamCluster(pixArrBuff, &labelsBuff, bbox, parent->mPattIdConverter, compClusPtr, patternsPtr, labelsClusPtr, nlab);
256 } else {
257 auto warnLeft = MaxHugeClusWarn - parent->mNHugeClus;
258 if (warnLeft > 0) {
259 LOGP(warn, "Splitting a huge cluster: chipID {}, rows {}:{} cols {}:{}{}", bbox.chipID, bbox.rowMin, bbox.rowMax, bbox.colMin, bbox.colMax,
260 warnLeft == 1 ? " (Further warnings will be muted)" : "");
261#ifdef WITH_OPENMP
262#pragma omp critical
263#endif
264 {
265 parent->mNHugeClus++;
266 }
267 }
268 BBox bboxT(bbox); // truncated box
269 std::vector<PixelData> pixbuf;
270 do {
271 bboxT.rowMin = bbox.rowMin;
272 bboxT.colMax = std::min(bbox.colMax, uint16_t(bboxT.colMin + o2::itsmft::ClusterPattern::MaxColSpan - 1));
273 do { // Select a subset of pixels fitting the reduced bounding box
274 bboxT.rowMax = std::min(bbox.rowMax, uint16_t(bboxT.rowMin + o2::itsmft::ClusterPattern::MaxRowSpan - 1));
275 for (const auto& pix : pixArrBuff) {
276 if (bboxT.isInside(pix.getRowDirect(), pix.getCol())) {
277 pixbuf.push_back(pix);
278 }
279 }
280 if (!pixbuf.empty()) { // Stream a piece of cluster only if the reduced bounding box is not empty
281 parent->streamCluster(pixbuf, &labelsBuff, bboxT, parent->mPattIdConverter, compClusPtr, patternsPtr, labelsClusPtr, nlab, true);
282 pixbuf.clear();
283 }
284 bboxT.rowMin = bboxT.rowMax + 1;
285 } while (bboxT.rowMin < bbox.rowMax);
286 bboxT.colMin = bboxT.colMax + 1;
287 } while (bboxT.colMin < bbox.colMax);
288 }
289 }
290}
291
292//__________________________________________________
294 PatternCont* patternsPtr, const ConstMCTruth* labelsDigPtr, MCTruth* labelsClusPtr)
295{
296 auto pix = curChipData->getData()[hit];
297 uint16_t row = pix.getRowDirect(), col = pix.getCol();
298
299 if (labelsClusPtr) { // MC labels were requested
300 int nlab = 0;
301 fetchMCLabels(curChipData->getStartID() + hit, labelsDigPtr, nlab);
302 auto cnt = compClusPtr->size();
303 for (int i = nlab; i--;) {
304 labelsClusPtr->addElement(cnt, labelsBuff[i]);
305 }
306 }
307
308 // add to compact clusters, which must be always filled
309 unsigned char patt[ClusterPattern::MaxPatternBytes]{0x1 << (7 - (0 % 8))}; // unrolled 1 hit version of full loop in finishChip
310 uint16_t pattID = (parent->mPattIdConverter.size() == 0) ? CompCluster::InvalidPatternID : parent->mPattIdConverter.findGroupID(1, 1, patt);
311 if ((pattID == CompCluster::InvalidPatternID || parent->mPattIdConverter.isGroup(pattID)) && patternsPtr) {
312 patternsPtr->emplace_back(1); // rowspan
313 patternsPtr->emplace_back(1); // colspan
314 patternsPtr->insert(patternsPtr->end(), std::begin(patt), std::begin(patt) + 1);
315 }
316 compClusPtr->emplace_back(row, col, pattID, curChipData->getChipID());
317}
318
319//__________________________________________________
321{
322#ifdef _PERFORM_TIMING_
323 mTimer.Stop();
324 mTimer.Reset();
325 mTimerMerge.Stop();
326 mTimerMerge.Reset();
327#endif
328}
329
330//__________________________________________________
332{
333 // init chip with the 1st unmasked pixel (entry "from" in the mChipData)
335 int chipId = curChipData->getChipID();
338 }
339
340 delete[] column1;
341 delete[] column2;
342 column1 = new int[size];
343 column2 = new int[size];
344 column1[0] = column1[size - 1] = -1;
345 column2[0] = column2[size - 1] = -1;
346 // init chip with the 1st unmasked pixel (entry "from" in the mChipData)
347 prev = column1 + 1;
348 curr = column2 + 1;
349 resetColumn(curr);
350
351 pixels.clear();
352 preClusterHeads.clear();
353 preClusterIndices.clear();
354 auto pix = curChipData->getData()[first];
355 currCol = pix.getCol();
356 curr[pix.getRowDirect()] = 0; // can use getRowDirect since the pixel is not masked
357 // start the first pre-cluster
358 preClusterHeads.push_back(0);
359 preClusterIndices.push_back(0);
360 pixels.emplace_back(-1, first); // id of current pixel
361 noLeftCol = true; // flag that there is no column on the left to check yet
362}
363
364//__________________________________________________
365void Clusterer::ClustererThread::updateChip(const ChipPixelData* curChipData, uint32_t ip)
366{
367 const auto pix = curChipData->getData()[ip];
368 uint16_t row = pix.getRowDirect(); // can use getRowDirect since the pixel is not masked
369 if (currCol != pix.getCol()) { // switch the buffers
370 swapColumnBuffers();
371 resetColumn(curr);
372 noLeftCol = false;
373 if (pix.getCol() > currCol + 1) {
374 // no connection with previous column, this pixel cannot belong to any of the
375 // existing preclusters, create a new precluster and flag to check only the row above for next pixels of this column
376 currCol = pix.getCol();
377 addNewPrecluster(ip, row);
378 noLeftCol = true;
379 return;
380 }
381 currCol = pix.getCol();
382 }
383
384 bool orphan = true;
385
386 if (noLeftCol) { // check only the row above
387 if (curr[row - 1] >= 0) {
388 expandPreCluster(ip, row, curr[row - 1]); // attach to the precluster of the previous row
389 return;
390 }
391 } else {
392#ifdef _ALLOW_DIAGONAL_ALPIDE_CLUSTERS_
393 int neighbours[]{curr[row - 1], prev[row], prev[row + 1], prev[row - 1]};
394#else
395 int neighbours[]{curr[row - 1], prev[row]};
396#endif
397 for (auto pci : neighbours) {
398 if (pci < 0) {
399 continue;
400 }
401 if (orphan) {
402 expandPreCluster(ip, row, pci); // attach to the adjascent precluster
403 orphan = false;
404 continue;
405 }
406 // reassign precluster index to smallest one
407 if (preClusterIndices[pci] < preClusterIndices[curr[row]]) {
408 preClusterIndices[curr[row]] = preClusterIndices[pci];
409 } else {
410 preClusterIndices[pci] = preClusterIndices[curr[row]];
411 }
412 }
413 }
414 if (orphan) {
415 addNewPrecluster(ip, row); // start new precluster
416 }
417}
418
419//__________________________________________________
420void Clusterer::ClustererThread::fetchMCLabels(int digID, const ConstMCTruth* labelsDig, int& nfilled)
421{
422 // transfer MC labels to cluster
423 if (nfilled >= MaxLabels) {
424 return;
425 }
426 const auto& lbls = labelsDig->getLabels(digID);
427 for (int i = lbls.size(); i--;) {
428 int ic = nfilled;
429 for (; ic--;) { // check if the label is already present
430 if (labelsBuff[ic] == lbls[i]) {
431 return; // label is found, do nothing
432 }
433 }
434 labelsBuff[nfilled++] = lbls[i];
435 if (nfilled >= MaxLabels) {
436 break;
437 }
438 }
439 //
440}
441
442//__________________________________________________
444{
445 // reset
446#ifdef _PERFORM_TIMING_
447 mTimer.Stop();
448 mTimer.Reset();
449 mTimerMerge.Stop();
450 mTimerMerge.Reset();
451#endif
452}
453
454//__________________________________________________
456{
457 LOGP(info, "Clusterizer has {} chips registered", mChips.size());
458 LOGP(info, "Clusterizer squashes overflow pixels separated by {} BC and <= {} in row/col seeking down to {} neighbour ROFs", mMaxBCSeparationToSquash, mMaxRowColDiffToMask, mSquashingDepth);
459 LOGP(info, "Clusterizer masks overflow pixels separated by < {} BC and <= {} in row/col", mMaxBCSeparationToMask, mMaxRowColDiffToMask);
460
461#ifdef _PERFORM_TIMING_
462 auto& tmr = const_cast<TStopwatch&>(mTimer); // ugly but this is what root does internally
463 auto& tmrm = const_cast<TStopwatch&>(mTimerMerge);
464 LOG(info) << "Inclusive clusterization timing (w/o disk IO): Cpu: " << tmr.CpuTime()
465 << " Real: " << tmr.RealTime() << " s in " << tmr.Counter() << " slots";
466 LOG(info) << "Threads output merging timing : Cpu: " << tmrm.CpuTime()
467 << " Real: " << tmrm.RealTime() << " s in " << tmrm.Counter() << " slots";
468
469#endif
470}
471} // namespace o2::its3
int32_t i
Definition of a container to keep Monte Carlo truth external to simulation objects.
uint32_t col
Definition RawData.h:4
Definition of the SegmentationSuperAlpide class.
Definition of the ITS cluster finder.
gsl::span< const TruthElement > getLabels(uint32_t dataindex) const
void mergeAtBack(MCTruthContainer< TruthElement > const &other)
void addElement(uint32_t dataindex, TruthElement const &element, bool noElement=false)
static constexpr int MaxLabels
Definition Clusterer.h:76
void process(int nThreads, PixelReader &r, CompClusCont *compClus, PatternCont *patterns, ROFRecCont *vecROFRec, MCTruth *labelsCl=nullptr)
Definition Clusterer.cxx:32
void clear()
load the dictionary of cluster topologies
void print() const
static constexpr int MaxHugeClusWarn
Definition Clusterer.h:77
uint32_t getOrderedPixId(int pos) const
Definition PixelData.h:275
uint32_t getStartID() const
Definition PixelData.h:110
uint16_t getChipID() const
Definition PixelData.h:108
const std::vector< PixelData > & getData() const
Definition PixelData.h:115
static constexpr int MaxPatternBytes
static constexpr uint8_t MaxRowSpan
static constexpr uint8_t MaxColSpan
static constexpr unsigned short InvalidPatternID
Definition CompCluster.h:46
PixelReader class for the ITSMFT.
Definition PixelReader.h:34
virtual int decodeNextTrigger()=0
bool getDecodeNextAuto() const
Definition PixelReader.h:67
void setDecodeNextAuto(bool v)
Definition PixelReader.h:68
const o2::InteractionRecord & getInteractionRecord() const
Definition PixelReader.h:58
virtual const o2::dataformats::ConstMCTruthContainerView< o2::MCCompLabel > * getDigitsMCTruth() const
Definition PixelReader.h:49
virtual bool getNextChipData(ChipPixelData &chipData)=0
const BCData & getBCData() const
Definition ROFRecord.h:58
GLsizeiptr size
Definition glcorearb.h:659
GLint first
Definition glcorearb.h:399
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLint GLint GLsizei GLint GLenum GLenum const void * pixels
Definition glcorearb.h:275
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
bool isDetITS3(T detID)
Definition SpecsV2.h:161
std::vector< unsigned char > PatternCont
Definition Clusterer.h:59
std::vector< itsmft::CompClusterExt > CompClusCont
Definition Clusterer.h:58
std::vector< itsmft::ROFRecord > ROFRecCont
Definition Clusterer.h:60
int64_t differenceInBC(const InteractionRecord &other) const
bool isInside(uint16_t row, uint16_t col) const
Definition Clusterer.h:86
std::vector< ThreadStat > stats
Definition Clusterer.h:147
void finishChipSingleHitFast(uint32_t hit, ChipPixelData *curChipData, CompClusCont *compClusPtr, PatternCont *patternsPtr, const ConstMCTruth *labelsDigPtr, MCTruth *labelsClusPTr)
void process(uint16_t chip, uint16_t nChips, CompClusCont *compClusPtr, PatternCont *patternsPtr, const ConstMCTruth *labelsDigPtr, MCTruth *labelsClPtr, const itsmft::ROFRecord &rofPtr)
void updateChip(const ChipPixelData *curChipData, uint32_t ip)
void fetchMCLabels(int digID, const ConstMCTruth *labelsDig, int &nfilled)
void initChip(const ChipPixelData *curChipData, uint32_t first)
void finishChip(ChipPixelData *curChipData, CompClusCont *compClus, PatternCont *patterns, const ConstMCTruth *labelsDig, MCTruth *labelsClus)
methods and transient data used within a thread
Definition Clusterer.h:114
std::vector< o2::mch::DsChannelId > chid
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< int > row