Project
Loading...
Searching...
No Matches
IDCFactorization.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#include "Framework/Logger.h"
18#include "TFile.h"
19#include "TPCBase/CalDet.h"
20#include <functional>
23#include "TKey.h"
24#include "TROOT.h"
25
26#if (defined(WITH_OPENMP) || defined(_OPENMP))
27#include <omp.h>
28#endif
29
30o2::tpc::IDCFactorization::IDCFactorization(const std::array<unsigned char, Mapper::NREGIONS>& groupPads, const std::array<unsigned char, Mapper::NREGIONS>& groupRows, const std::array<unsigned char, Mapper::NREGIONS>& groupLastRowsThreshold, const std::array<unsigned char, Mapper::NREGIONS>& groupLastPadsThreshold, const unsigned int groupNotnPadsSectorEdges, const unsigned int timeFrames, const unsigned int timeframesDeltaIDC, const std::vector<uint32_t>& crus)
31 : IDCGroupHelperSector{groupPads, groupRows, groupLastRowsThreshold, groupLastPadsThreshold, groupNotnPadsSectorEdges}, mTimeFrames{timeFrames}, mTimeFramesDeltaIDC{timeframesDeltaIDC}, mCRUs{crus}
32{
34 if (mSides.size() == 1) {
35 mSideIndex[1] = 0;
36 }
37
38 for (int i = 0; i < mSides.size(); ++i) {
39 mIDCDelta.emplace_back(std::vector<IDCDelta<float>>{timeFrames / timeframesDeltaIDC + (timeFrames % timeframesDeltaIDC != 0)});
40 }
41
42 mIDCZero.resize(mSides.size());
43 mIDCOne.resize(mSides.size());
44
45 for (auto& idc : mIDCs) {
46 idc.resize(mTimeFrames);
47 }
48 // check if the input IDCs are grouped
49 for (int region = 0; region < Mapper::NREGIONS; ++region) {
50 if (mNIDCsPerCRU[region] != Mapper::PADSPERREGION[region]) {
51 mInputGrouped = true;
52 break;
53 }
54 }
55}
56
58
59std::vector<o2::tpc::Side> o2::tpc::IDCFactorization::getSides(const std::vector<uint32_t>& crus)
60{
61 std::unordered_map<o2::tpc::Side, bool> map;
62 for (auto cru : crus) {
63 map[CRU(cru).side()] = true;
64 }
65
66 std::vector<o2::tpc::Side> sides;
67 if (!(map.find(o2::tpc::Side::A) == map.end())) {
68 sides.emplace_back(o2::tpc::Side::A);
69 }
70 if (!(map.find(o2::tpc::Side::C) == map.end())) {
71 sides.emplace_back(o2::tpc::Side::C);
72 }
73 return sides;
74}
75
76void o2::tpc::IDCFactorization::dumpToFile(const char* outFileName, const char* outName) const
77{
78 TFile fOut(outFileName, "RECREATE");
79 fOut.WriteObject(this, outName);
80 fOut.Close();
81}
82
83void o2::tpc::IDCFactorization::dumpLargeObjectToFile(const char* outFileName, const char* outName) const
84{
85 LOGP(info, "Writing object to file using {} threads", sNThreads);
86 ROOT::EnableThreadSafety();
87 std::vector<std::unique_ptr<TFile>> fOutPar;
88 fOutPar.reserve(sNThreads);
89 for (int i = 0; i < sNThreads; ++i) {
90 const std::string outFile = fmt::format("{}_{}", outFileName, i);
91 fOutPar.emplace_back(TFile::Open(outFile.data(), "RECREATE"));
92 }
93
94#pragma omp parallel for num_threads(sNThreads)
95 for (int iTF = 0; iTF < mTimeFrames; ++iTF) {
96 IDCFactorizeSplit idcTmp;
97 idcTmp.relTF = iTF;
98 for (int icru = 0; icru < CRU::MaxCRU; ++icru) {
99 idcTmp.idcs[icru] = mIDCs[icru][iTF];
100 }
101#ifdef WITH_OPENMP
102 const int ithread = omp_get_thread_num();
103#else
104 const int ithread = 0;
105#endif
106 fOutPar[ithread]->WriteObject(&idcTmp, fmt::format("IDCs_TF{}", iTF).data());
107 }
108
109 // write empty dummy object to file which can be filled with the previously written IDCs
110 TFile fOut(outFileName, "RECREATE");
111 IDCFactorization idcFacTmp(mTimeFrames, mTimeFramesDeltaIDC, mCRUs);
112 idcFacTmp.setRun(getRun());
113 idcFacTmp.setTimeStamp(getTimeStamp());
114 fOut.WriteObject(&idcFacTmp, outName);
115 std::vector<int> nFiles{sNThreads};
116 fOut.WriteObject(&nFiles, "out_files");
117 fOut.Close();
118}
119
120std::unique_ptr<o2::tpc::IDCFactorization> o2::tpc::IDCFactorization::getLargeObjectFromFile(const char* inpFileName, const char* inName)
121{
122 LOGP(info, "Reading object from file using {} threads", sNThreads);
123 TFile fMeta(inpFileName, "READ");
124 if (fMeta.IsZombie()) {
125 LOGP(warning, "Input file {} not found! Returning", inpFileName);
126 }
127
128 o2::tpc::IDCFactorization* idcTmp = (o2::tpc::IDCFactorization*)fMeta.Get(inName);
129 std::unique_ptr<IDCFactorization> idc;
130 idc = std::make_unique<IDCFactorization>(idcTmp->getNTimeframes(), idcTmp->getTimeFramesDeltaIDC(), idcTmp->getCRUs());
131 idc->setRun(idcTmp->getRun());
132 idc->setTimeStamp(idcTmp->getTimeStamp());
133 std::vector<int>* out_files = (std::vector<int>*)fMeta.Get("out_files");
134 const int nFiles = out_files->front();
135 delete idcTmp;
136 delete out_files;
137
138 std::vector<std::unique_ptr<TFile>> fInp;
139 fInp.reserve(nFiles);
140 for (int i = 0; i < nFiles; ++i) {
141 const std::string inpFile = fmt::format("{}_{}", inpFileName, i);
142 fInp.emplace_back(TFile::Open(inpFile.data()));
143 if (fInp.front() == nullptr) {
144 LOGP(warning, "Input file {} not found!", inpFile.data());
145 }
146 }
147
148 ROOT::EnableThreadSafety();
149
150// loop over files
151#pragma omp parallel for num_threads(sNThreads)
152 for (int iFile = 0; iFile < nFiles; ++iFile) {
153 const auto* keys = fInp[iFile]->GetListOfKeys();
154 const auto nKeys = keys->GetEntries();
155 for (int iKey = 0; iKey < nKeys; ++iKey) {
156 const auto key = dynamic_cast<TKey*>(keys->At(iKey));
157 IDCFactorizeSplit* idcTmp = nullptr;
158 fInp[iFile]->GetObject(key->GetName(), idcTmp);
159
160 const int relTF = idcTmp->relTF;
161 if (relTF >= idc->getNTimeframes()) {
162 LOGP(warning, "stored TF {} is larger than max TF {}", relTF, idc->getNTimeframes());
163 continue;
164 }
165 for (int icru = 0; icru < CRU::MaxCRU; ++icru) {
166 auto idcVec = idcTmp->idcs[icru];
167 idc->setIDCs(std::move(idcVec), icru, relTF);
168 }
169 delete idcTmp;
170 }
171 }
172 for (auto& file : fInp) {
173 file->Close();
174 }
175 return idc;
176}
177
178void o2::tpc::IDCFactorization::dumpIDCZeroToFile(const Side side, const char* outFileName, const char* outName) const
179{
180 TFile fOut(outFileName, "RECREATE");
181 fOut.WriteObject(&mIDCZero[mSideIndex[side]], outName);
182 fOut.Close();
183}
184
185void o2::tpc::IDCFactorization::dumpIDCOneToFile(const Side side, const char* outFileName, const char* outName) const
186{
187 TFile fOut(outFileName, "RECREATE");
188 fOut.WriteObject(&mIDCOne[mSideIndex[side]], outName);
189 fOut.Close();
190}
191
192void o2::tpc::IDCFactorization::dumpToTree(const Side side, const char* outFileName)
193{
194 if (mIDCZero[mSideIndex[side]].mIDCZero.empty() || mIDCOne[mSideIndex[side]].mIDCOne.empty()) {
195 LOGP(info, "IDC0, IDC1 is empty! Returning");
196 return;
197 }
199 if (!mPadFlagsMap) {
200 createStatusMap();
201 }
202 helper.setPadStatusMap(*mPadFlagsMap);
203 helper.setIDCZero(&mIDCZero[mSideIndex[side]], side);
204 helper.setIDCOne(&mIDCOne[mSideIndex[side]], side);
205 helper.dumpToTree(side, outFileName);
206}
207
208void o2::tpc::IDCFactorization::dumpToTreeIDC1(const float integrationTimeOrbits, const char* outFileName) const
209{
210 o2::utils::TreeStreamRedirector pcstream(outFileName, "RECREATE");
211 pcstream.GetFile()->cd();
212 std::vector<double> timestamp;
213
214 for (int i = 0; i < getIDCOneVec(Side::A).size(); ++i) {
215 timestamp.emplace_back((mTimeStamp + i * integrationTimeOrbits * o2::constants::lhc::LHCOrbitMUS * 0.001) / 1000);
216 }
217
218 pcstream << "tree"
219 << "IDC1A=" << getIDCOne(o2::tpc::Side::A)
220 << "IDC1C=" << getIDCOne(o2::tpc::Side::C)
221 << "timestamp=" << timestamp
222 << "\n";
223
224 pcstream.Close();
225}
226
227void o2::tpc::IDCFactorization::dumpIDCDeltaToTree(const Side side, const int chunk, const char* outFileName)
228{
229 if (chunk >= getTimeFramesDeltaIDC()) {
230 LOGP(info, "index {} is larger than maximum index {} for IDCDelta", chunk, getTimeFramesDeltaIDC() - 1);
231 }
233 IDCGroupHelperSector parameters = *(dynamic_cast<const IDCGroupHelperSector*>(this));
234 helper.setGroupingParameter(&parameters, side);
235 helper.setIDCDelta(&mIDCDelta[mSideIndex[side]][chunk], side);
236 helper.setIDCOne(&mIDCOne[mSideIndex[side]], side);
237 helper.dumpToTreeIDCDelta(side, outFileName);
238}
239
241{
242 const unsigned int nIDCsSide = mNIDCsPerSector * o2::tpc::SECTORSPERSIDE;
243 for (auto& idcZero : mIDCZero) {
244 idcZero.clear();
245 idcZero.resize(nIDCsSide);
246 }
247
248#pragma omp parallel for num_threads(sNThreads)
249 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
250 const unsigned int cru = mCRUs[cruInd];
251 const o2::tpc::CRU cruTmp(cru);
252 const auto side = cruTmp.side();
253 const unsigned int region = cruTmp.region();
254 const auto factorIndexGlob = mRegionOffs[region] + mNIDCsPerSector * (cruTmp.sector() % o2::tpc::SECTORSPERSIDE);
255 for (unsigned int timeframe = 0; timeframe < mTimeFrames; ++timeframe) {
256 for (unsigned int idcs = 0; idcs < mIDCs[cru][timeframe].size(); ++idcs) {
257 if ((mIDCs[cru][timeframe][idcs] == -1) || (mIDCs[cru][timeframe][idcs] == 0)) {
258 continue;
259 }
260 if (norm) {
261 mIDCs[cru][timeframe][idcs] *= Mapper::INVPADAREA[region];
262 }
263 const unsigned int indexGlob = (idcs % mNIDCsPerCRU[region]) + factorIndexGlob;
264 mIDCZero[mSideIndex[side]].fillValueIDCZero(mIDCs[cru][timeframe][idcs], indexGlob);
265 }
266 }
267 }
268
269// perform normalization per CRU (in case some CRUs lack data)
270#pragma omp parallel for num_threads(sNThreads)
271 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
272 const unsigned int cru = mCRUs[cruInd];
273 const o2::tpc::CRU cruTmp(cru);
274 const auto side = cruTmp.side();
275 const unsigned int region = cruTmp.region();
276 const auto indexStart = mRegionOffs[region] + mNIDCsPerSector * (cruTmp.sector() % o2::tpc::SECTORSPERSIDE);
277 const auto indexEnd = indexStart + mNIDCsPerCRU[region];
278
279 const auto normFac = getNIntegrationIntervals(cru);
280 if (normFac == 0) {
281 LOGP(info, "number of integration intervals is zero for CRU {}! Skipping normalization of IDC0", cru);
282 continue;
283 }
284
285 std::transform(mIDCZero[mSideIndex[side]].mIDCZero.begin() + indexStart, mIDCZero[mSideIndex[side]].mIDCZero.begin() + indexEnd, mIDCZero[mSideIndex[side]].mIDCZero.begin() + indexStart, [normVal = normFac](auto& val) { return val / normVal; });
286 }
287}
288
290{
291#pragma omp parallel for num_threads(sNThreads)
292 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
293 const unsigned int cru = mCRUs[cruInd];
294 const o2::tpc::CRU cruTmp(cru);
295 const Side side = cruTmp.side();
296 const int region = cruTmp.region();
297 const int sector = cruTmp.sector();
298 std::vector<unsigned int> padTmp;
299 padTmp.reserve(Mapper::PADSPERROW[region].back());
300 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
301 float idcRow = 0;
302 int count = 0;
303 const unsigned int integrationInterval = 0;
304 // loop over pad row and calculate mean of IDC0
305 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
306 const auto index = getIndexUngrouped(sector, region, lrow, pad, integrationInterval);
307 const auto idcZeroVal = mIDCZero[mSideIndex[side]].mIDCZero[index];
308 if (idcZeroVal > 0) {
309 const auto globalPad = Mapper::getGlobalPadNumber(lrow, pad, region);
310 const float gain = (mGainMap ? mGainMap->getValue(sector, globalPad) : 1);
311 idcRow += idcZeroVal / gain;
312 ++count;
313 } else {
314 padTmp.emplace_back(pad);
315 }
316 }
317 // loop over dead pads and set value
318 for (const auto pad : padTmp) {
319 const float meanIDC0 = idcRow / count;
320 const auto globalPad = Mapper::getGlobalPadNumber(lrow, pad, region);
321 const float gain = (mGainMap ? mGainMap->getValue(sector, globalPad) : 1);
322 const float idcZero = gain * meanIDC0;
323 const auto index = getIndexUngrouped(sector, region, lrow, pad, integrationInterval);
324 mIDCZero[mSideIndex[side]].setValueIDCZero(idcZero, index);
325 }
326 padTmp.clear();
327 }
328 }
329}
330
332{
333 const unsigned int integrationIntervals = getNIntegrationIntervals();
334 for (auto& idcOne : mIDCOne) {
335 idcOne.mIDCOne.clear();
336 idcOne.mIDCOneMedian.clear();
337 idcOne.mIDCOneRMS.clear();
338 idcOne.mIDCOne.resize(integrationIntervals);
339 idcOne.mIDCOneMedian.resize(integrationIntervals);
340 idcOne.mIDCOneRMS.resize(integrationIntervals);
341 std::fill(idcOne.mIDCOne.begin(), idcOne.mIDCOne.end(), 1);
342 }
343
344#pragma omp parallel for num_threads(sNThreads)
345 for (unsigned int timeframe = 0; timeframe < mTimeFrames; ++timeframe) {
346 const unsigned int nIDCsSide = mNIDCsPerSector * SECTORSPERSIDE;
347 std::vector<std::vector<std::vector<float>>> idcOneSafe(mSides.size()); // side -> intervals -> IDCs
348
349 // reserve memory
350 for (auto& interavalvec : idcOneSafe) {
351 interavalvec.resize(mIntegrationIntervalsPerTF[timeframe]);
352 for (auto& idcsPerInteraval : interavalvec) {
353 idcsPerInteraval.reserve(nIDCsSide);
354 }
355 }
356
357 // loop over CRUs and fill idcOneSafe vector with IDC/IDC0 values
358 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
359 const unsigned int cru = mCRUs[cruInd];
360 const o2::tpc::CRU cruTmp(cru);
361 const unsigned int region = cruTmp.region();
362 const auto side = cruTmp.side();
363 const unsigned int indexSide = mSideIndex[side];
364 const auto factorIndexGlob = mRegionOffs[region] + mNIDCsPerSector * (cruTmp.sector() % o2::tpc::SECTORSPERSIDE);
365 unsigned int integrationIntervalOffset = 0;
366 calcIDCOne(mIDCs[cru][timeframe], mNIDCsPerCRU[region], integrationIntervalOffset, factorIndexGlob, cru, idcOneSafe[indexSide], &mIDCZero[indexSide], mInputGrouped ? nullptr : mPadFlagsMap.get(), mUsePadStatusMap);
367 }
368
369 // calculate global offset for interval
370 int offsInterval = 0;
371 for (unsigned int tf = 0; tf < timeframe; ++tf) {
372 offsInterval += mIntegrationIntervalsPerTF[tf];
373 }
374
375 for (int side = 0; side < mSides.size(); ++side) {
376 const unsigned int indexSide = mSideIndex[side];
377 for (int interval = 0; interval < mIntegrationIntervalsPerTF[timeframe]; ++interval) {
378
379 // calculate robust average for each slice
380 const int intervalGlobal = interval + offsInterval;
381 if (!idcOneSafe[indexSide].empty()) {
382 RobustAverage average(std::move(idcOneSafe[indexSide][interval]));
383 const float mean = average.getTrunctedMean(0.05, 0.95);
384 const float median = average.getMedian();
385 const float rms = average.getStdDev();
386 if (mean != 0) {
387 mIDCOne[indexSide].mIDCOne[intervalGlobal] = mean;
388 }
389 mIDCOne[indexSide].mIDCOneMedian[intervalGlobal] = median;
390 mIDCOne[indexSide].mIDCOneRMS[intervalGlobal] = rms;
391 }
392 }
393 }
394 }
395}
396
397template <typename DataVec>
398void o2::tpc::IDCFactorization::calcIDCOne(const DataVec& idcsData, const int idcsPerCRU, const int integrationIntervalOffset, const unsigned int indexOffset, const CRU cru, std::vector<std::vector<float>>& idcOneTmp, const IDCZero* idcZero, const CalDet<PadFlags>* flagMap, const bool usePadStatusMap)
399{
400 for (unsigned int idcs = 0; idcs < idcsData.size(); ++idcs) {
401 if ((idcsData[idcs] == -1) || (idcsData[idcs] == 0)) {
402 continue;
403 }
404
405 const unsigned int integrationInterval = idcs / idcsPerCRU + integrationIntervalOffset;
406 const unsigned int localPad = (idcs % idcsPerCRU);
407 const unsigned int indexGlob = localPad + indexOffset;
408 const auto idcZeroVal = idcZero ? idcZero->mIDCZero[indexGlob] : 1;
409
410 // check pad in case of input is not grouped
411 if (usePadStatusMap && flagMap) {
412 const o2::tpc::PadFlags flag = flagMap->getCalArray(cru).getValue(localPad);
413 if ((flag & PadFlags::flagSkip) == PadFlags::flagSkip) {
414 continue;
415 }
416 }
417
418 if (idcOneTmp.size() <= integrationInterval) {
419 LOGP(error, "integrationInterval {} is larger than maximum: {}", integrationInterval, idcOneTmp.size());
420 break;
421 }
422 const double epsilon = 0.001;
423 if (std::abs(idcZeroVal) > epsilon) {
424 idcOneTmp[integrationInterval].emplace_back(idcsData[idcs] / idcZeroVal);
425 }
426 }
427}
428
430{
431 const unsigned int nIDCsSide = mNIDCsPerSector * SECTORSPERSIDE;
432 for (auto side : mSides) {
433 for (unsigned int i = 0; i < getNChunks(side); ++i) {
434 const auto idcsSide = nIDCsSide * getNIntegrationIntervalsInChunk(i);
435 mIDCDelta[mSideIndex[side]][i].getIDCDelta().clear();
436 mIDCDelta[mSideIndex[side]][i].getIDCDelta().resize(idcsSide);
437 }
438 }
439
440#pragma omp parallel for num_threads(sNThreads)
441 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
442 const unsigned int cru = mCRUs[cruInd];
443 const o2::tpc::CRU cruTmp(cru);
444 const unsigned int region = cruTmp.region();
445 const auto side = cruTmp.side();
446 const auto factorIndexGlob = mRegionOffs[region] + mNIDCsPerSector * (cruTmp.sector() % o2::tpc::SECTORSPERSIDE);
447 unsigned int integrationIntervallast = 0;
448 unsigned int integrationIntervallastLocal = 0;
449 unsigned int lastChunk = 0;
450
451 for (unsigned int timeframe = 0; timeframe < mTimeFrames; ++timeframe) {
452 const unsigned int chunk = getChunk(timeframe);
453 if (lastChunk != chunk) {
454 integrationIntervallastLocal = 0;
455 }
456
457 for (unsigned int idcs = 0; idcs < mIDCs[cru][timeframe].size(); ++idcs) {
458 if ((mIDCs[cru][timeframe][idcs] == -1) || (mIDCs[cru][timeframe][idcs] == 0)) {
459 continue;
460 }
461 const unsigned int intervallocal = idcs / mNIDCsPerCRU[region];
462 const unsigned int integrationIntervalGlobal = intervallocal + integrationIntervallast;
463 const unsigned int integrationIntervalLocal = intervallocal + integrationIntervallastLocal;
464 const unsigned int localPad = idcs % mNIDCsPerCRU[region];
465 const unsigned int indexGlob = localPad + factorIndexGlob;
466 const auto idcZero = mIDCZero[mSideIndex[side]].mIDCZero[indexGlob];
467 const auto idcOne = mIDCOne[mSideIndex[side]].mIDCOne[integrationIntervalGlobal];
468 const auto mult = idcZero * idcOne;
469 auto val = (mult != 0) ? mIDCs[cru][timeframe][idcs] / mult - 1 : 0;
470 if (mUsePadStatusMap && !mInputGrouped && mPadFlagsMap) {
471 const o2::tpc::PadFlags flag = mPadFlagsMap->getCalArray(cru).getValue(localPad);
472 if ((flag & PadFlags::flagSkip) == PadFlags::flagSkip) {
473 val = 0;
474 }
475 }
476 mIDCDelta[mSideIndex[side]][chunk].setIDCDelta(indexGlob + integrationIntervalLocal * nIDCsSide, val);
477 }
478
479 const unsigned int intervals = mIntegrationIntervalsPerTF[timeframe];
480 integrationIntervallast += intervals;
481 integrationIntervallastLocal += intervals;
482 lastChunk = chunk;
483 }
484 }
485}
486
488{
489 createStatusMapOutlier();
490 checkFECs();
491 checkNeighbourOutliers();
492}
493
495{
496 const static auto& paramIDCGroup = ParameterIDCGroup::Instance();
497 mPadFlagsMap = std::make_unique<CalDet<PadFlags>>(CalDet<PadFlags>("flags", PadSubset::Region));
498
499 std::array<std::vector<o2::tpc::CRU>, SIDES> crus;
500 const int crusPerSide = CRU::MaxCRU / 2;
501 crus[Side::A].reserve(crusPerSide);
502 crus[Side::C].reserve(crusPerSide);
503 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
504 o2::tpc::CRU cru(mCRUs[cruInd]);
505 crus[cru.side()].emplace_back(cru);
506 }
507
508 // normalize IDC per stack mean
509 const auto stacksMedian = getStackMedian();
510 const std::array<int, Mapper::NREGIONS> stackInSector{0, 0, 0, 0, 1, 1, 2, 2, 3, 3};
511
512 for (const auto side : mSides) {
513 std::vector<o2::tpc::RobustAverage> average; // [0]: all pads except edges and centre, centre (cross), edge pads
514 average.reserve(Mapper::PADROWS);
515 for (int row = 0; row < Mapper::PADROWS; ++row) {
516 average.emplace_back(Mapper::PADSPERROW[Mapper::REGION[row]].back());
517 }
518
519 for (int iter = 0; iter < 2; ++iter) {
520 std::vector<std::pair<float, float>> mean;
521 if (iter == 1) {
522 mean.reserve(Mapper::PADROWS);
523 for (int row = 0; row < Mapper::PADROWS; ++row) {
524 mean.emplace_back(average[row].getFilteredAverage());
525 }
526
527 if (debug) {
528 o2::utils::TreeStreamRedirector pcstream(fmt::format("rob_average_run_{}_side_{}.root", mRun, static_cast<int>(side)).data(), "RECREATE");
529 for (int row = 0; row < average.size(); ++row) {
530 std::vector<float> values = average[row].getValues();
531 int nValues = average[row].getValues().size();
532 float median = average[row].getMedian();
533 pcstream << "tree"
534 << "values=" << values
535 << "mean=" << mean[row]
536 << "row=" << row
537 << "nValues=" << nValues
538 << "median=" << median
539 << "\n";
540 }
541 }
542 }
543
544 for (const auto& cru : crus[side]) {
545 const int region = cru.region();
546 const int sector = cru.sector();
547
548 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
549 // loop over pads in row in the first iteration and calculate median at the end
550 // in the second iteration check if the IDC value is to far away from the median
551 const unsigned int integrationInterval = 0;
552 const int globalRow = Mapper::ROWOFFSET[region] + lrow;
553
554 // loop over pad row and calculate mean of IDC0
555 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
556 const auto index = getIndexUngrouped(sector, region, lrow, pad, integrationInterval);
557 auto idcZeroVal = mIDCZero[mSideIndex[side]].mIDCZero[index];
558
559 // correct by gain map: expected improvement at the cross regions with large gain variation
560 if (mGainMap) {
561 const auto globalPad = Mapper::getGlobalPadNumber(lrow, pad, region);
562 const float gain = mGainMap->getValue(sector, globalPad);
563 idcZeroVal /= gain;
564 }
565
566 const bool checkSimple = (idcZeroVal != -1) && (idcZeroVal != 0);
567 if (checkSimple) {
568 idcZeroVal /= stacksMedian[sector * GEMSTACKSPERSECTOR + stackInSector[region]];
569 }
570
571 if (iter == 0) {
572 // exclude obvious outliers
573 if (checkSimple) {
574 average[globalRow].addValue(idcZeroVal);
575 }
576 } else {
577 const int relPad = static_cast<int>(pad) - static_cast<int>(Mapper::PADSPERROW[region][lrow] / 2);
578 // nPads pads on the left and nPads pads on the right = 2 * nPads
579 const int nPadsCentre = 2;
580 const bool isSecCentreEdge = ((relPad >= -nPadsCentre) && (relPad < nPadsCentre)) || ((pad == 0) || (pad == Mapper::PADSPERROW[region][lrow] - 1));
581 const float idc0MedianMin = isSecCentreEdge ? paramIDCGroup.minmaxIDC0MedianCentreEdge : paramIDCGroup.minIDC0Median;
582 const float idc0MedianMax = isSecCentreEdge ? paramIDCGroup.minmaxIDC0MedianCentreEdge : paramIDCGroup.maxIDC0Median;
583
584 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][lrow] + pad;
586 if ((idcZeroVal == -1) || (idcZeroVal == 0)) {
587 flag = o2::tpc::PadFlags::flagDeadPad | o2::tpc::PadFlags::flagSkip | mPadFlagsMap->getCalArray(cru).getValue(padInRegion);
588 } else if (idcZeroVal > mean[globalRow].first + mean[globalRow].second * idc0MedianMax) {
589 flag = o2::tpc::PadFlags::flagHighPad | o2::tpc::PadFlags::flagSkip | mPadFlagsMap->getCalArray(cru).getValue(padInRegion);
590 } else if (idcZeroVal < mean[globalRow].first - mean[globalRow].second * idc0MedianMin) {
591 flag = o2::tpc::PadFlags::flagLowPad | o2::tpc::PadFlags::flagSkip | mPadFlagsMap->getCalArray(cru).getValue(padInRegion);
592 }
593 mPadFlagsMap->getCalArray(cru).setValue(padInRegion, flag);
594 }
595 } // loop over pads in row
596 } // end iteration
597 }
598 }
599 }
600}
601
602void o2::tpc::IDCFactorization::checkFECs(const float maxOutliersPerFEC)
603{
604 /*
605 1. loop over FECs
606 2. count number of outliers
607 3. compare number of outliers with pads in FEC
608 4. if too many pads are outliers, the whole FEC is masked
609 */
610
611 if (!mPadFlagsMap) {
612 createStatusMapOutlier();
613 }
614
615 std::array<int, FECInfo::FECsTotal> nOutliersPerFEC{0};
616 const Mapper& mapper = Mapper::instance();
617 for (int iter = 0; iter < 2; ++iter) {
618 for (unsigned int sector = 0; sector < Mapper::NSECTORS; ++sector) {
619 for (unsigned int region = 0; region < Mapper::NREGIONS; ++region) {
620 const CRU cru(Sector(sector), region);
621
622 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
623 // loop over pad row and calculate mean of IDC0
624 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
625 const FECInfo& fec = mapper.fecInfo(Mapper::getGlobalPadNumber(lrow, pad, region));
626 const size_t fecIndex = sector * FECInfo::FECsPerSector + fec.getIndex();
627 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][lrow] + pad;
628 o2::tpc::PadFlags flag = mPadFlagsMap->getCalArray(cru).getValue(padInRegion);
629
630 if (iter == 0) {
631 // count outliers per FEC
632 if ((flag & PadFlags::flagSkip) == PadFlags::flagSkip) {
633 ++nOutliersPerFEC[fecIndex];
634 }
635 } else {
636 // mark whole FEC
637 if (nOutliersPerFEC[fecIndex] / static_cast<float>(FECInfo::ChannelsPerFEC) > maxOutliersPerFEC) {
638 // turn bit off
640 // set bit
642 mPadFlagsMap->getCalArray(cru).setValue(padInRegion, flag);
643 }
644 }
645 } // loop over pads in row
646 } // end iteration
647 } // loop region
648 } // loop sector
649 } // loop iter
650}
651
652void o2::tpc::IDCFactorization::checkNeighbourOutliers(const int maxIter, const int nOutliersNeighbours)
653{
654 for (int iter = 0; iter < maxIter; ++iter) {
655 int countTotOutlier = 0;
656 for (unsigned int sector = 0; sector < Mapper::NSECTORS; ++sector) {
657 for (unsigned int region = 0; region < Mapper::NREGIONS; ++region) {
658 const CRU cru(Sector(sector), region);
659 const int maxRow = Mapper::ROWSPERREGION[region] - 1;
660 for (int lrow = 0; lrow <= maxRow; ++lrow) {
661 // find first row of the cluster
662 const int nPadsInRow = 1; // +- one pad in row direction
663 const int nPadsInPad = 2; // +- two pads in pad direction
664 const int rowStart = std::clamp(lrow - nPadsInRow, 0, maxRow);
665 const int rowEnd = std::clamp(lrow + nPadsInRow, 0, maxRow);
666 const int addPadsStart = Mapper::ADDITIONALPADSPERROW[region][lrow];
667
668 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
669 int countOutliers = 0;
670 // loop ove the rows from the cluster
671 for (int rowCl = rowStart; rowCl <= rowEnd; ++rowCl) {
672 // shift local pad in row in case current row from the cluster has more pads in the row
673 const int addPadsCl = Mapper::ADDITIONALPADSPERROW[region][rowCl];
674 const int diffAddPads = addPadsCl - addPadsStart;
675 const int padClCentre = pad + diffAddPads;
676
677 const int maxPad = Mapper::PADSPERROW[region][rowCl] - 1;
678 const int padStart = std::clamp(padClCentre - nPadsInPad, 0, maxPad);
679 const int padEnd = std::clamp(padClCentre + nPadsInPad, 0, maxPad);
680 for (int padCl = padStart; padCl <= padEnd; ++padCl) {
681 // skip for current cluster position as the charge there is not effected from the threshold
682 if ((padCl == pad) && (rowCl == lrow)) {
683 continue;
684 }
685
686 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][rowCl] + padCl;
687 o2::tpc::PadFlags flag = mPadFlagsMap->getCalArray(cru).getValue(padInRegion);
688 if ((flag & PadFlags::flagSkip) == PadFlags::flagSkip) {
689 ++countOutliers;
690 }
691 }
692 }
693
694 if (countOutliers >= nOutliersNeighbours) {
695 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][lrow] + pad;
696 o2::tpc::PadFlags flag = mPadFlagsMap->getCalArray(cru).getValue(padInRegion);
697
699 ++countTotOutlier;
700 }
701
702 // turn bit off
704 // set bit
706 mPadFlagsMap->getCalArray(cru).setValue(padInRegion, flag);
707 }
708 }
709 }
710 }
711 }
712 LOGP(debug, "iter: {} count total outlier: {}", iter, countTotOutlier);
713 if (countTotOutlier == 0) {
714 break;
715 }
716 }
717}
718
720{
721 if (!mPadFlagsMap) {
722 LOGP(info, "getMeanZ(): Status map not set returning");
723 return -1;
724 }
725 return IDCCCDBHelper<>::getMeanIDC0(side, mIDCZero[mSideIndex[side]], mPadFlagsMap.get());
726}
727
728float o2::tpc::IDCFactorization::getIDCValUngrouped(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const
729{
730 unsigned int timeFrame = 0;
731 unsigned int interval = 0;
732 getTF(integrationInterval, timeFrame, interval);
733 if (mIDCs[sector * Mapper::NREGIONS + region][timeFrame].empty()) {
734 return 0.f;
735 }
736
737 const int index = interval * mNIDCsPerCRU[region] + mOffsRow[region][getGroupedRow(region, urow)] + getGroupedPad(region, urow, upad) + getOffsetForEdgePad(upad, urow, region);
738 return mIDCs[sector * Mapper::NREGIONS + region][timeFrame][index];
739}
740
741void o2::tpc::IDCFactorization::getTF(unsigned int integrationInterval, unsigned int& timeFrame, unsigned int& interval) const
742{
743 const auto& intervalsPerTF = mIntegrationIntervalsPerTF.empty() ? getAllIntegrationIntervalsPerTF() : mIntegrationIntervalsPerTF;
744 unsigned int nintervals = 0;
745 unsigned int intervalTmp = 0;
746 for (unsigned int tf = 0; tf < mTimeFrames; ++tf) {
747 nintervals += intervalsPerTF[tf];
748 if (integrationInterval < nintervals) {
749 timeFrame = tf;
750 interval = integrationInterval - intervalTmp;
751 return;
752 }
753 intervalTmp = nintervals;
754 }
755}
756
758{
759 // find three consecutive time frames
760 const int nMaxTFs = 3;
761 const int nTFs = (mTimeFrames >= nMaxTFs) ? nMaxTFs : mTimeFrames;
762
763 // store the number of integration intervals like 10 -> 11 -> 11
764 std::vector<unsigned int> ordering;
765 ordering.reserve(nTFs);
766 unsigned int order = 0;
767 int firstTF = 0;
768
769 // find ordering
770 for (int iTF = 0; iTF < mTimeFrames; ++iTF) {
771 bool found = false;
772 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
773 const unsigned int cru = mCRUs[cruInd];
774 if (!mIDCs[cru][iTF].empty()) {
775 order = mIDCs[cru][iTF].size() / mNIDCsPerCRU[CRU(cru).region()];
776 ordering.emplace_back(order);
777 found = true;
778 break;
779 }
780 }
781 if (ordering.size() == nTFs) {
782 firstTF = iTF - nTFs + 1;
783 break;
784 }
785 if (!found) {
786 ordering.clear();
787 }
788 }
789
790 // in case TFs < 3
791 if ((nTFs == mTimeFrames) && (ordering.size() == nTFs)) {
792 return ordering;
793 }
794
795 if (ordering.size() != nTFs) {
796 LOGP(warning, "Couldnt find {} consecutive TFs with data", nTFs);
797 ordering = std::vector<unsigned int>(nTFs, order);
798 firstTF = 0;
799 } else {
800 std::sort(ordering.begin(), ordering.end());
801 }
802
803 std::vector<unsigned int> integrationIntervalsPerTF(mTimeFrames);
804 for (int iter = 0; iter < 2; ++iter) {
805 const int end = (iter == 0) ? (mTimeFrames - firstTF) : firstTF;
806 for (int iTFLoop = 0; iTFLoop < end; ++iTFLoop) {
807 const int iTF = (iter == 0) ? (firstTF + iTFLoop) : (firstTF - iTFLoop - 1);
808 int intervalsInTF = -1;
809 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
810 const unsigned int cru = mCRUs[cruInd];
811 if (!mIDCs[cru][iTF].empty()) {
812 intervalsInTF = mIDCs[cru][iTF].size() / mNIDCsPerCRU[CRU(cru).region()];
813 break;
814 }
815 }
816
817 if (intervalsInTF == -1) {
818 if ((ordering.size() != nTFs)) {
819 if (iTF == 0) {
820 intervalsInTF = ordering.front();
821 } else if (iTF == 1) {
822 intervalsInTF = ordering[1];
823 }
824 } else {
825 const int m1Val = (iter == 0) ? -1 : +1;
826 const int intervalsLastm1 = integrationIntervalsPerTF[iTF + m1Val];
827 const int intervalsLastm2 = integrationIntervalsPerTF[iTF + 2 * m1Val];
828 if ((intervalsLastm1 == ordering[1]) && (intervalsLastm2 == ordering[2])) {
829 intervalsInTF = ordering[0];
830 } else {
831 intervalsInTF = ordering[1];
832 }
833 }
834 }
835
836 if (intervalsInTF < -1) {
837 LOGP(warning, "interval is smaller than 0!");
838 continue;
839 }
840
841 integrationIntervalsPerTF[iTF] = intervalsInTF;
842 }
843 }
844 return integrationIntervalsPerTF;
845}
846
847void o2::tpc::IDCFactorization::factorizeIDCs(const bool norm, const bool calcDeltas)
848{
849 using timer = std::chrono::high_resolution_clock;
850
851 LOGP(info, "Using {} threads for factorization of IDCs", sNThreads);
852 LOGP(info, "Checking received IDCs for consistency");
853 checkReceivedIDCs();
854
855 LOGP(info, "Calculating IDC0");
856
857 auto start = timer::now();
858 calcIDCZero(norm);
859 auto stop = timer::now();
860 std::chrono::duration<float> time = stop - start;
861 float totalTime = time.count();
862 LOGP(info, "IDCZero time: {}", time.count());
863
864 if (!mInputGrouped) {
865 LOGP(info, "Creating pad status map");
866 start = timer::now();
867 createStatusMap();
868 stop = timer::now();
869 time = stop - start;
870 LOGP(info, "Pad-by-pad status map time: {}", time.count());
871 totalTime += time.count();
872 }
873
874 start = timer::now();
875 mIntegrationIntervalsPerTF = getAllIntegrationIntervalsPerTF();
876 stop = timer::now();
877 time = stop - start;
878 totalTime = time.count();
879 LOGP(info, "Getting integration intervals for all TFs time: {}", time.count());
880
881 LOGP(info, "Calculating IDC1");
882 start = timer::now();
883 calcIDCOne();
884 stop = timer::now();
885 time = stop - start;
886 LOGP(info, "IDC1 time: {}", time.count());
887 totalTime += time.count();
888
889 if (calcDeltas) {
890 LOGP(info, "Calculating IDCDelta");
891 start = timer::now();
892 calcIDCDelta();
893 stop = timer::now();
894 time = stop - start;
895 LOGP(info, "IDCDelta time: {}", time.count());
896 totalTime += time.count();
897 }
898
899 LOGP(info, "Factorization done. Total time: {}", totalTime);
900}
901
902unsigned long o2::tpc::IDCFactorization::getNIntegrationIntervalsInChunk(const unsigned int chunk) const
903{
904 const auto firstTF = chunk * getNTFsPerChunk(0);
905 const auto& intervalsPerTF = mIntegrationIntervalsPerTF.empty() ? getAllIntegrationIntervalsPerTF() : mIntegrationIntervalsPerTF;
906 const auto intervals = std::reduce(intervalsPerTF.begin() + firstTF, intervalsPerTF.begin() + firstTF + getNTFsPerChunk(chunk));
907 return intervals;
908}
909
910unsigned long o2::tpc::IDCFactorization::getNIntegrationIntervalsToChunk(const unsigned int chunk) const
911{
912 const auto chunkTF = chunk * getNTFsPerChunk(0);
913 const auto& intervalsPerTF = mIntegrationIntervalsPerTF.empty() ? getAllIntegrationIntervalsPerTF() : mIntegrationIntervalsPerTF;
914 const auto intervals = std::reduce(intervalsPerTF.begin(), intervalsPerTF.begin() + chunkTF);
915 return intervals;
916}
917
919{
920 std::size_t sum = 0;
921 for (auto&& idcsTF : mIDCs[cru]) {
922 sum += idcsTF.size();
923 }
924 return sum / mNIDCsPerCRU[CRU(cru).region()];
925}
926
928{
929 const auto& intervalsPerTF = mIntegrationIntervalsPerTF.empty() ? getAllIntegrationIntervalsPerTF() : mIntegrationIntervalsPerTF;
930 const auto intervals = std::reduce(intervalsPerTF.begin(), intervalsPerTF.end());
931 return intervals;
932}
933
934void o2::tpc::IDCFactorization::getLocalIntegrationInterval(const unsigned int integrationInterval, unsigned int& chunk, unsigned int& localintegrationInterval) const
935{
936 const auto& intervalsPerTF = mIntegrationIntervalsPerTF.empty() ? getAllIntegrationIntervalsPerTF() : mIntegrationIntervalsPerTF;
937 unsigned int nintervals = 0;
938 unsigned int nitervalsChunk = 0;
939 unsigned int globalTF = 0;
940 for (unsigned int ichunk = 0; ichunk < getNChunks(mSides.front()); ++ichunk) {
941 const auto nTFsPerChunk = getNTFsPerChunk(ichunk);
942 for (unsigned int tf = 0; tf < nTFsPerChunk; ++tf) {
943 nintervals += intervalsPerTF[globalTF];
944 if (integrationInterval < nintervals) {
945 chunk = getChunk(globalTF);
946 localintegrationInterval = integrationInterval - nitervalsChunk;
947 return;
948 }
949 ++globalTF;
950 }
951 nitervalsChunk = nintervals;
952 }
953}
954
955unsigned int o2::tpc::IDCFactorization::getChunk(const unsigned int timeframe) const
956{
957 return timeframe / mTimeFramesDeltaIDC;
958}
959
960unsigned int o2::tpc::IDCFactorization::getNTFsPerChunk(const unsigned int chunk) const
961{
962 const unsigned int remain = mTimeFrames % mTimeFramesDeltaIDC;
963 return ((chunk == getNChunks(mSides.front()) - 1) && remain) ? remain : mTimeFramesDeltaIDC;
964}
965
966std::vector<unsigned int> o2::tpc::IDCFactorization::getIntegrationIntervalsPerTF(const int cru) const
967{
968 const uint32_t cruTmp = (cru < 0) ? mCRUs.front() : cru;
969 const auto region = CRU(cruTmp).region();
970 std::vector<unsigned int> integrationIntervalsPerTF;
971 integrationIntervalsPerTF.reserve(mTimeFrames);
972 for (unsigned int tf = 0; tf < mTimeFrames; ++tf) {
973 integrationIntervalsPerTF.emplace_back(mIDCs[cruTmp][tf].size() / mNIDCsPerCRU[region]);
974 }
975 return integrationIntervalsPerTF;
976}
977
979{
980 for (auto& tf : mIDCs) {
981 for (auto& idcs : tf) {
982 idcs.clear();
983 }
984 }
985}
986
987void o2::tpc::IDCFactorization::drawIDCDeltaHelper(const bool type, const Sector sector, const unsigned int integrationInterval, const IDCDeltaCompression compression, const std::string filename, const float minZ, const float maxZ) const
988{
989 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc;
990
991 unsigned int chunk = 0;
992 unsigned int localintegrationInterval = 0;
993 getLocalIntegrationInterval(integrationInterval, chunk, localintegrationInterval);
994 const std::string zAxisTitle = IDCDrawHelper::getZAxisTitle(IDCType::IDCDelta, compression);
995
997 switch (compression) {
999 default: {
1000 idcFunc = [this, chunk, localintegrationInterval](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
1001 return this->getIDCDeltaVal(sector, region, irow, pad, chunk, localintegrationInterval);
1002 };
1003 drawFun.mIDCFunc = idcFunc;
1004 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZ) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZ);
1005 break;
1006 }
1008 const auto idcDeltaMedium = this->getIDCDeltaMediumCompressed(chunk, sector.side());
1009 idcFunc = [this, &idcDeltaMedium, localintegrationInterval = localintegrationInterval](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
1010 return idcDeltaMedium.getValue(this->getIndexUngrouped(sector, region, irow, pad, localintegrationInterval));
1011 };
1012 drawFun.mIDCFunc = idcFunc;
1013 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZ) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZ);
1014 break;
1015 }
1017 const auto idcDeltaHigh = this->getIDCDeltaHighCompressed(chunk, sector.side());
1018 idcFunc = [this, &idcDeltaHigh, localintegrationInterval](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
1019 return idcDeltaHigh.getValue(this->getIndexUngrouped(sector, region, irow, pad, localintegrationInterval));
1020 };
1021 drawFun.mIDCFunc = idcFunc;
1022 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZ) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZ);
1023 break;
1024 }
1025 }
1026}
1027
1028void o2::tpc::IDCFactorization::drawIDCZeroHelper(const bool type, const Sector sector, const std::string filename, const float minZ, const float maxZ) const
1029{
1030 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc = [this](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
1031 return this->getIDCZeroVal(sector, region, irow, pad);
1032 };
1033
1034 IDCDrawHelper::IDCDraw drawFun;
1035 drawFun.mIDCFunc = idcFunc;
1036 const std::string zAxisTitle = IDCDrawHelper::getZAxisTitle(IDCType::IDCZero);
1037
1038 const float maxZTmp = (maxZ == -1) ? (3 * getMeanZ(sector.side())) : maxZ;
1039 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZTmp) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZTmp);
1040}
1041
1042void o2::tpc::IDCFactorization::drawIDCHelper(const bool type, const Sector sector, const unsigned int integrationInterval, const std::string filename, const float minZ, const float maxZ, const bool drawGIF, const int run) const
1043{
1044 const float maxZTmp = (maxZ == -1) ? (3 * getMeanZ(sector.side())) : maxZ;
1045 const std::string zAxisTitleDraw = IDCDrawHelper::getZAxisTitle(IDCType::IDC);
1046 if (!drawGIF) {
1047 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc = [this, integrationInterval](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
1048 return this->getIDCValUngrouped(sector, region, irow, pad, integrationInterval);
1049 };
1050
1051 IDCDrawHelper::IDCDraw drawFun;
1052 drawFun.mIDCFunc = idcFunc;
1053
1054 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitleDraw, filename, minZ, maxZTmp) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitleDraw, filename, minZ, maxZTmp);
1055 } else {
1056 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc = [this](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad, const unsigned int slice) {
1057 return this->getIDCValUngrouped(sector, region, irow, pad, slice);
1058 };
1059
1060 if (mIDCOne.front().mIDCOne.empty()) {
1061 LOGP(warning, "Factorised IDCs are missing!");
1062 return;
1063 }
1064
1065 std::function<float(const o2::tpc::Side, const unsigned int)> idcOneFunc = [this](const o2::tpc::Side side, const unsigned int slice) {
1066 return this->getIDCOneVec(side)[slice];
1067 };
1068
1069 IDCDrawHelper::IDCDrawGIF drawFun;
1070 drawFun.mIDCFunc = idcFunc;
1071 drawFun.mIDCOneFunc = idcOneFunc;
1072 const int nSlices = (integrationInterval == 0) ? getNIntegrationIntervals() : integrationInterval;
1073 IDCDrawHelper::drawSideGIF(drawFun, nSlices, zAxisTitleDraw, filename, minZ, maxZ, run);
1074 }
1075}
1076
1077void o2::tpc::IDCFactorization::setGainMap(const char* inpFile, const char* mapName)
1078{
1079 TFile f(inpFile, "READ");
1080 o2::tpc::CalDet<float>* gainMap = nullptr;
1081 f.GetObject(mapName, gainMap);
1082
1083 if (!gainMap) {
1084 LOGP(info, "GainMap {} not found returning", mapName);
1085 return;
1086 }
1087 setGainMap(*gainMap);
1088 delete gainMap;
1089}
1090
1092{
1093 mGainMap = std::make_unique<CalDet<float>>(gainmap);
1094}
1095
1096void o2::tpc::IDCFactorization::setPadFlagMap(const char* inpFile, const char* mapName)
1097{
1098 TFile f(inpFile, "READ");
1099 o2::tpc::CalDet<PadFlags>* statusmap = nullptr;
1100 f.GetObject(mapName, statusmap);
1101
1102 if (!statusmap) {
1103 LOGP(info, "Pad flag map {} not found returning", mapName);
1104 return;
1105 }
1106 setPadFlagMap(*statusmap);
1107 delete statusmap;
1108}
1109
1111{
1112 mPadFlagsMap = std::make_unique<CalDet<PadFlags>>(flagmap);
1113}
1114
1115void o2::tpc::IDCFactorization::drawPadFlagMap(const bool type, const Sector sector, const std::string filename, const PadFlags flag) const
1116{
1117 if (!mPadFlagsMap) {
1118 LOGP(info, "Status map not set returning");
1119 return;
1120 }
1121
1122 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc = [this, flag](const unsigned int sector, const unsigned int region, const unsigned int row, const unsigned int pad) {
1123 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][row] + pad;
1124 const auto flagDraw = mPadFlagsMap->getCalArray(region + sector * Mapper::NREGIONS).getValue(padInRegion);
1125 if ((flagDraw & flag) == flag) {
1126 return 1;
1127 } else {
1128 return 0;
1129 }
1130 };
1131
1132 IDCDrawHelper::IDCDraw drawFun;
1133 drawFun.mIDCFunc = idcFunc;
1134 const std::string zAxisTitle = "status flag";
1135 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename);
1136}
1137
1138void o2::tpc::IDCFactorization::dumpPadFlagMap(const char* outFile, const char* mapName)
1139{
1140 if (!mPadFlagsMap) {
1141 LOGP(info, "Status map not set returning");
1142 }
1143 TFile fOut(outFile, "RECREATE");
1144 fOut.WriteObject(mPadFlagsMap.get(), mapName);
1145 fOut.Close();
1146}
1147
1148void o2::tpc::IDCFactorization::normIDCZero(const int type)
1149{
1150 if ((type == 0) && !mGainMap) {
1151 LOGP(info, "Gain map not set returning");
1152 return;
1153 }
1154
1155 const std::array<int, Mapper::NREGIONS> stackInSector{0, 0, 0, 0, 1, 1, 2, 2, 3, 3};
1156 std::array<float, o2::tpc::GEMSTACKS> stacksMedian;
1157 if (type == 1) {
1158 stacksMedian = getStackMedian();
1159 }
1160
1161#pragma omp parallel for num_threads(sNThreads)
1162 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
1163 const unsigned int cru = mCRUs[cruInd];
1164 const o2::tpc::CRU cruTmp(cru);
1165 const Side side = cruTmp.side();
1166 const int region = cruTmp.region();
1167 const int sector = cruTmp.sector();
1168 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
1169 const unsigned int integrationInterval = 0;
1170 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
1171 const auto index = getIndexUngrouped(sector, region, lrow, pad, integrationInterval);
1172 const auto globalPad = Mapper::getGlobalPadNumber(lrow, pad, region);
1173 float normVal = 1;
1174 if (type == 0) {
1175 normVal = mGainMap->getValue(sector, globalPad);
1176 } else if (type == 1) {
1177 normVal = stacksMedian[sector * GEMSTACKSPERSECTOR + stackInSector[region]];
1178 }
1179 mIDCZero[mSideIndex[side]].mIDCZero[index] /= normVal;
1180 }
1181 }
1182 }
1183}
1184
1185std::array<float, o2::tpc::GEMSTACKS> o2::tpc::IDCFactorization::getStackMedian() const
1186{
1187 std::array<float, GEMSTACKS> median;
1188 for (const auto side : mSides) {
1189 const auto tmpMedian = IDCCCDBHelper<>::getStackMedian(mIDCZero[mSideIndex[side]], side);
1191 std::copy(tmpMedian.begin(), tmpMedian.end(), median.begin() + offs);
1192 }
1193 return median;
1194}
1195
1196bool o2::tpc::IDCFactorization::checkReceivedIDCs()
1197{
1198 bool idcsGood = true;
1199 for (unsigned int timeframe = 0; timeframe < mTimeFrames; ++timeframe) {
1200 // check number of received slices for each CRU - this should be the same value or if no data has been received empty -
1201 std::unordered_map<int, std::vector<int>> receivedIDCsSize; // number of received slices, CRUs
1202 for (unsigned int cruInd = 0; cruInd < mCRUs.size(); ++cruInd) {
1203 const unsigned int cru = mCRUs[cruInd];
1204 const o2::tpc::CRU cruTmp(cru);
1205 const unsigned int region = cruTmp.region();
1206 const int nSlices = mIDCs[cru][timeframe].size() / mNIDCsPerCRU[region];
1207 if (nSlices > 0) {
1208 if (receivedIDCsSize[nSlices].empty()) {
1209 receivedIDCsSize[nSlices].reserve(mCRUs.size());
1210 }
1211 receivedIDCsSize[nSlices].emplace_back(cru);
1212 }
1213 }
1214 // if more than one slice has been found use the IDCs from CRUs which have the most equal slices
1215 if (receivedIDCsSize.size() > 1) {
1216 LOGP(warning, "Received inconsistent IDCs!");
1217 idcsGood = false;
1218 // check which slice has the most CRUs - reset all other -
1219 std::pair<int, int> nSlicesMostCRUs = {-1, -1}; // nSlices, nCRUs
1220 for (auto nSlices : receivedIDCsSize) {
1221 const int nCRUs = nSlices.second.size();
1222 // store number of slices with most CRUs
1223 if (nCRUs > nSlicesMostCRUs.second) {
1224 nSlicesMostCRUs = {nSlices.first, nCRUs};
1225 }
1226 for (auto cru : nSlices.second) {
1227 LOGP(info, "Received {} slices for CRU {}", nSlices.first, cru);
1228 }
1229 }
1230 LOGP(info, "Setting {} slices for {} valid CRUs", nSlicesMostCRUs.first, nSlicesMostCRUs.second);
1231 for (auto nSlices : receivedIDCsSize) {
1232 if (nSlices.first != nSlicesMostCRUs.first) {
1233 for (auto cru : nSlices.second) {
1234 LOGP(info, "Clearing IDCs for CRU {}", cru);
1235 mIDCs[cru][timeframe].clear();
1236 }
1237 }
1238 }
1239 }
1240 }
1241 return idcsGood;
1242}
1243
1244template void o2::tpc::IDCFactorization::calcIDCOne(const o2::pmr::vector<float>&, const int, const int, const unsigned int, const CRU, std::vector<std::vector<float>>&, const IDCZero*, const CalDet<PadFlags>*, const bool);
int16_t time
Definition RawEventData.h:4
int32_t i
helper class for accessing IDCs from CCDB
helper class for drawing IDCs per region/side
class for aggregating IDCs for the full TPC (all sectors) and factorization of aggregated IDCs
Header to collect LHC related constants.
uint32_t side
Definition RawData.h:0
class for performing robust averaging and outlier filtering
std::ostringstream debug
StringRef key
unsigned char region() const
Definition CRU.h:64
Side side() const
Definition CRU.h:62
@ MaxCRU
Definition CRU.h:31
const Sector sector() const
Definition CRU.h:65
const T getValue(const size_t channel) const
Definition CalArray.h:96
const CalType & getCalArray(size_t position) const
Definition CalDet.h:63
static constexpr int FECsPerSector
Definition FECInfo.h:28
static constexpr int ChannelsPerFEC
Definition FECInfo.h:27
unsigned char getIndex() const
Definition FECInfo.h:41
void setIDCZero(IDCZero *idcZero, const Side side=Side::A)
setting the 0D-IDCs
void setIDCDelta(IDCDelta< DataT > *idcDelta, const Side side=Side::A)
setting the IDCDelta class member
void setIDCOne(IDCOne *idcOne, const Side side=Side::A)
setting the 1D-IDCs
static float getMeanIDC0(const Side side, const IDCZero &idcZero, const CalDet< PadFlags > *outlierMap)
void dumpToTree(const Side side, const char *outFileName="IDCCCDBTree.root") const
void setPadStatusMap(const CalDet< PadFlags > &outliermap)
static float getStackMedian(const IDCZero &idczero, const Sector sector, const GEMstack stack)
void setGroupingParameter(IDCGroupHelperSector *helperSector, const Side side=Side::A)
setting the grouping parameters
void dumpToTreeIDCDelta(const Side side, const char *outFileName="IDCCCDBTreeDeltaIDC.root") const
static void drawSideGIF(const IDCDrawGIF &idcs, const unsigned int slices, const std::string zAxisTitle, const std::string filename="IDCs", const float minZ=0, const float maxZ=-1, const int run=-1)
make a GIF of IDCs for A and C side and the 1D IDCs
static std::string getZAxisTitle(const IDCType type, const IDCDeltaCompression compression=IDCDeltaCompression::NO)
static void drawSide(const IDCDraw &idc, const o2::tpc::Side side, const std::string zAxisTitle, const std::string filename, const float minZ=0, const float maxZ=-1)
std::array< float, o2::tpc::GEMSTACKS > getStackMedian() const
void setGainMap(const char *inpFile, const char *mapName)
const std::vector< Side > & getSides() const
void calcIDCOne()
calculate I_1(t) = <I(r,\phi,t) / I_0(r,\phi)>_{r,\phi}
void dumpIDCDeltaToTree(const Side side, const int chunk=0, const char *outFileName="IDCDeltaTree.root")
dump IDCDelta to a TTree
unsigned int getNTimeframes() const
void setRun(const int run)
void reset()
resetting aggregated IDCs
void fillIDCZeroDeadPads()
fill I_0 values in case of dead pads,FECs etc.
unsigned long getNIntegrationIntervalsInChunk(const unsigned int chunk) const
void dumpToTree(const Side side, const char *outFileName="IDCTree.root")
unsigned int getTimeFramesDeltaIDC() const
IDCFactorization()=default
default constructor for ROOT I/O
unsigned long getNIntegrationIntervalsToChunk(const unsigned int chunk) const
void factorizeIDCs(const bool norm, const bool calcDeltas)
void dumpPadFlagMap(const char *outFile, const char *mapName)
void dumpIDCZeroToFile(const Side side, const char *outFileName="IDCZero.root", const char *outName="IDC0") const
dump the IDC0 to file
float getIDCValUngrouped(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const
void checkNeighbourOutliers(const int maxIter=10, const int nOutliersNeighbours=8)
std::vector< unsigned int > getIntegrationIntervalsPerTF(const int cru=-1) const
std::vector< unsigned int > getAllIntegrationIntervalsPerTF() const
unsigned long getNIntegrationIntervals() const
void dumpToFile(const char *outFileName="IDCFactorized.root", const char *outName="IDCFactorized") const
void calcIDCZero(const bool norm)
static std::unique_ptr< IDCFactorization > getLargeObjectFromFile(const char *inpFileName="IDCFactorized.root", const char *inName="IDCFactorized")
void dumpLargeObjectToFile(const char *outFileName="IDCFactorized.root", const char *outName="IDCFactorized") const
void calcIDCDelta()
calculate \Delta I(r,\phi,t) = I(r,\phi,t) / ( I_0(r,\phi) * I_1(t) )
void checkFECs(const float maxOutliersPerFEC=0.7)
void getLocalIntegrationInterval(const unsigned int integrationInterval, unsigned int &chunk, unsigned int &localintegrationInterval) const
float getMeanZ(const Side side) const
void dumpToTreeIDC1(const float integrationTimeOrbits=12, const char *outFileName="IDC1Tree.root") const
void dumpIDCOneToFile(const Side side, const char *outFileName="IDCOne.root", const char *outName="IDC1") const
dump the IDC1 to file
void createStatusMapOutlier(const bool debug=false)
void setTimeStamp(const long timeStamp)
void setPadFlagMap(const char *inpFile, const char *mapName)
Helper class for accessing grouped pads for one sector.
std::array< unsigned int, Mapper::NREGIONS > mNIDCsPerCRU
total number of IDCs per region per integration interval
static constexpr unsigned int ROWOFFSET[NREGIONS]
offset to calculate local row from global row
Definition Mapper.h:533
static GlobalPadNumber getGlobalPadNumber(const unsigned int lrow, const unsigned int pad, const unsigned int region)
Definition Mapper.h:64
static const std::vector< unsigned int > ADDITIONALPADSPERROW[NREGIONS]
additional pads per row compared to first row
Definition Mapper.h:543
static const std::vector< unsigned int > PADSPERROW[NREGIONS]
number of pads per row in region
Definition Mapper.h:567
static const std::vector< unsigned int > OFFSETCRULOCAL[NREGIONS]
row offset in cru for given local pad row
Definition Mapper.h:555
const FECInfo & fecInfo(GlobalPadNumber padNumber) const
Definition Mapper.h:52
static constexpr float INVPADAREA[NREGIONS]
inverse size of the pad area padwidth*padLength
Definition Mapper.h:536
static Mapper & instance(const std::string mappingDir="")
Definition Mapper.h:44
static constexpr unsigned int ROWSPERREGION[NREGIONS]
number of pad rows for region
Definition Mapper.h:532
static constexpr unsigned int NSECTORS
total number of sectors in the TPC
Definition Mapper.h:526
static constexpr unsigned int NREGIONS
total number of regions in one sector
Definition Mapper.h:527
static constexpr unsigned int PADROWS
total number of pad rows
Definition Mapper.h:528
static constexpr unsigned int PADSPERREGION[NREGIONS]
number of pads per CRU
Definition Mapper.h:530
static constexpr unsigned REGION[PADROWS]
region for global pad row
Definition Mapper.h:537
float getTrunctedMean(float low, float high)
GLint GLsizei count
Definition glcorearb.h:399
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLuint end
Definition glcorearb.h:469
GLuint index
Definition glcorearb.h:781
GLdouble f
Definition glcorearb.h:310
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLenum GLsizei GLsizei GLint * values
Definition glcorearb.h:1576
GLboolean * data
Definition glcorearb.h:298
GLuint GLfloat * val
Definition glcorearb.h:1582
GLuint start
Definition glcorearb.h:469
GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat maxZ
Definition glcorearb.h:2910
GLfloat GLfloat minZ
Definition glcorearb.h:2910
constexpr double LHCOrbitMUS
std::vector< T, o2::pmr::polymorphic_allocator< T > > vector
constexpr unsigned char SECTORSPERSIDE
Definition Defs.h:40
@ IDC
integrated and grouped IDCs
@ IDCZero
IDC0: I_0(r,\phi) = <I(r,\phi,t)>_t.
@ IDCDelta
IDCDelta: \Delta I(r,\phi,t) = I(r,\phi,t) / ( I_0(r,\phi) * I_1(t) )
constexpr unsigned char SIDES
Definition Defs.h:41
@ Region
Regions (up to 36*10)
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
@ C
Definition Defs.h:36
IDCDeltaCompression
IDC Delta IDC Compression types.
@ HIGH
high compression using char (data compression ratio ~5.5 when stored in CCDB)
@ NO
no compression using floats
@ MEDIUM
medium compression using short (data compression ratio 2 when stored in CCDB)
PadFlags
Definition Defs.h:100
@ flagLowPad
flag for pad with extremly low IDC value
@ flagNeighbour
flag if n neighbouring pads are outlier
@ flagGoodPad
flag for a good pad binary 0001
@ flagSkip
flag for defining a pad which is just ignored during the calculation of I1 and IDCDelta
@ flagDeadPad
flag for a dead pad binary 0010
@ flagHighPad
flag for pad with extremly high IDC value
@ flagFEC
flag for a whole masked FEC
constexpr unsigned short GEMSTACKSPERSECTOR
Definition Defs.h:57
DataT sum(const Vector< DataT, N > &a)
Definition Vector.h:107
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string filename()
void empty(int)
std::unique_ptr< GPUReconstructionTimeframe > tf
struct to access and set Delta IDCs
std::function< float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> mIDCFunc
function returning the value which will be drawn for sector, region, row, pad
std::array< std::vector< float >, CRU::MaxCRU > idcs
struct containing the ITPC0 values (integrated TPC clusters)
std::vector< float > mIDCZero
I_0(r,\phi) = <I(r,\phi,t)>_t.
std::vector< int > row