Project
Loading...
Searching...
No Matches
IDCCCDBHelper.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
17
18#include "TPCBase/CalDet.h"
19#include "TPCBase/Mapper.h"
21#include "TPCBase/Painter.h"
22
23#include "TStyle.h"
24#include "TLine.h"
25#include "TCanvas.h"
26#include "TH2Poly.h"
27#include "TGraph.h"
28#include "TText.h"
30#include <map>
31
32#include "TProfile.h"
33
34template <typename DataT>
36{
37 return (mIDCDelta[side] && mHelperSector[side]) ? mIDCDelta[side]->getNIDCs() / (mHelperSector[side]->getNIDCsPerSector() * SECTORSPERSIDE) : 0;
38}
39
40template <typename DataT>
42{
43 return mIDCOne[side] ? mIDCOne[side]->getNIDCs() : 0;
44}
45
46template <typename DataT>
47float o2::tpc::IDCCCDBHelper<DataT>::getIDCZeroVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad) const
48{
50 const auto side = Sector(sector).side();
51 return !mIDCZero[side] ? -1 : (mIDCZero[side]->getNIDC0() == Mapper::getNumberOfPadsPerSide()) ? mIDCZero[side]->getValueIDCZero(getUngroupedIndexGlobal(sector, region, urow, upad, 0))
52 : mIDCZero[side]->getValueIDCZero(mHelperSector[side]->getIndexUngrouped(sector, region, urow, upad, 0));
53}
54
55template <typename DataT>
56float o2::tpc::IDCCCDBHelper<DataT>::getIDCDeltaVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const
57{
58 const auto side = Sector(sector).side();
59 return (!mIDCDelta[side] || !mHelperSector[side]) ? -1 : mIDCDelta[side]->getValue(mHelperSector[side]->getIndexUngrouped(sector, region, urow, upad, integrationInterval));
60}
61
62template <typename DataT>
63float o2::tpc::IDCCCDBHelper<DataT>::getIDCOneVal(const o2::tpc::Side side, const unsigned int integrationInterval) const
64{
65 return !mIDCOne[side] ? -1 : mIDCOne[side]->getValueIDCOne(integrationInterval);
66}
67
68template <typename DataT>
69float o2::tpc::IDCCCDBHelper<DataT>::getIDCVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const
70{
71 return (getIDCDeltaVal(sector, region, urow, upad, integrationInterval) + 1.f) * getIDCZeroVal(sector, region, urow, upad) * getIDCOneVal(Sector(sector).side(), integrationInterval);
72}
73
74template <typename DataT>
75void o2::tpc::IDCCCDBHelper<DataT>::drawIDCZeroHelper(const bool type, const o2::tpc::Sector sector, const std::string filename, const float minZ, const float maxZ) const
76{
77 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) {
78 return this->getIDCZeroVal(sector, region, irow, pad);
79 };
82 drawFun.mIDCFunc = idcFunc;
83 const std::string zAxisTitle = IDCDrawHelper::getZAxisTitle(IDCType::IDCZero);
84 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZ) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZ);
85}
87template <typename DataT>
88void o2::tpc::IDCCCDBHelper<DataT>::drawIDCDeltaHelper(const bool type, const Sector sector, const unsigned int integrationInterval, const std::string filename, const float minZ, const float maxZ) const
90 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) {
91 return this->getIDCDeltaVal(sector, region, irow, pad, integrationInterval);
92 };
93
95 drawFun.mIDCFunc = idcFunc;
96 const std::string zAxisTitle = IDCDrawHelper::getZAxisTitle(IDCType::IDCDelta);
97 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZ) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZ);
98}
99
100template <typename DataT>
101void o2::tpc::IDCCCDBHelper<DataT>::drawIDCHelper(const bool type, const Sector sector, const unsigned int integrationInterval, const std::string filename, const float minZ, const float maxZ) const
102{
103 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) {
104 return this->getIDCVal(sector, region, irow, pad, integrationInterval);
105 };
106
108 drawFun.mIDCFunc = idcFunc;
109 const std::string zAxisTitle = IDCDrawHelper::getZAxisTitle(IDCType::IDC);
110 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename, minZ, maxZ) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename, minZ, maxZ);
111}
113template <typename DataT>
114void o2::tpc::IDCCCDBHelper<DataT>::drawPadFlagMap(const bool type, const Sector sector, const std::string filename, const PadFlags flag) const
115{
116 if (!mPadFlagsMap) {
117 LOGP(info, "Status map not set returning");
118 return;
119 }
121 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) {
122 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][row] + pad;
123 const auto flagDraw = mPadFlagsMap->getCalArray(region + sector * Mapper::NREGIONS).getValue(padInRegion);
124 if ((flagDraw & flag) == flag) {
125 return 1;
126 } else {
127 return 0;
128 }
129 };
130
132 drawFun.mIDCFunc = idcFunc;
133 const std::string zAxisTitle = "status flag";
134 type ? IDCDrawHelper::drawSide(drawFun, sector.side(), zAxisTitle, filename) : IDCDrawHelper::drawSector(drawFun, 0, Mapper::NREGIONS, sector, zAxisTitle, filename);
135}
136
137template <typename DataT>
138unsigned int o2::tpc::IDCCCDBHelper<DataT>::getUngroupedIndexGlobal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval)
139{
140 return IDCGroupHelperSector::getUngroupedIndexGlobal(sector, region, urow, upad, integrationInterval);
141}
142
143template <typename DataT>
144TCanvas* o2::tpc::IDCCCDBHelper<DataT>::drawIDCZeroCanvas(TCanvas* outputCanvas, std::string_view type, int nbins1D, float xMin1D, float xMax1D, int integrationInterval) const
145{
146 TCanvas* canv = nullptr;
147
148 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc;
149
150 if (type == "IDC0") {
151 idcFunc = [this](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
152 return this->getIDCZeroVal(sector, region, irow, pad);
153 };
154 } else if (type == "IDCDelta") {
155 idcFunc = [this, integrationInterval](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
156 return this->getIDCDeltaVal(sector, region, irow, pad, integrationInterval);
157 };
158 if (integrationInterval <= 0) {
159 integrationInterval = std::min(getNIntegrationIntervalsIDCDelta(Side::A), getNIntegrationIntervalsIDCDelta(Side::C));
160 }
161 } else if (type == "IDC") {
162 idcFunc = [this, integrationInterval](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
163 return this->getIDCVal(sector, region, irow, pad, integrationInterval);
164 };
165 if (integrationInterval <= 0) {
166 integrationInterval = std::min(getNIntegrationIntervalsIDCDelta(Side::A), getNIntegrationIntervalsIDCDelta(Side::C));
167 }
168 } else {
169 LOGP(error, "Please provide a valid IDC data type. 'IDC0', 'IDCDelta', or 'IDC'.");
170 return canv;
171 }
173 if (outputCanvas) {
174 canv = outputCanvas;
175 } else {
176 canv = new TCanvas(fmt::format("c_sides_{}", type.data()).data(), fmt::format("{}", type.data()).data(), 1000, 1000);
177 }
180 drawFun.mIDCFunc = idcFunc;
181 const std::string zAxisTitle = IDCDrawHelper::getZAxisTitle(IDCType::IDCZero);
182
183 auto hAside2D = IDCDrawHelper::drawSide(drawFun, Side::A, zAxisTitle);
184 auto hCside2D = IDCDrawHelper::drawSide(drawFun, Side::C, zAxisTitle);
185 hAside2D->SetTitle(fmt::format("{} (A-Side)", type.data()).data());
186 hCside2D->SetTitle(fmt::format("{} (C-Side)", type.data()).data());
187 hAside2D->SetMinimum(xMin1D);
188 hAside2D->SetMaximum(xMax1D);
189 hCside2D->SetMinimum(xMin1D);
190 hCside2D->SetMaximum(xMax1D);
191
192 auto hAside1D = IDCDrawHelper::drawSide(drawFun, fmt::format("{}", type.data()).data(), Side::A, nbins1D, xMin1D, xMax1D);
193 auto hCside1D = IDCDrawHelper::drawSide(drawFun, fmt::format("{}", type.data()).data(), Side::C, nbins1D, xMin1D, xMax1D);
194
195 canv->Divide(2, 2);
196 canv->cd(1);
197 hAside2D->Draw("colz");
198 hAside2D->SetStats(0);
199 hAside2D->SetTitleOffset(1.05, "XY");
200 hAside2D->SetTitleSize(0.05, "XY");
202
203 canv->cd(2);
204 hCside2D->Draw("colz");
205 hCside2D->SetStats(0);
206 hCside2D->SetTitleOffset(1.05, "XY");
207 hCside2D->SetTitleSize(0.05, "XY");
209
210 canv->cd(3);
211 hAside1D->Draw();
212
213 canv->cd(4);
214 hCside1D->Draw();
215
216 hAside2D->SetBit(TObject::kCanDelete);
217 hCside2D->SetBit(TObject::kCanDelete);
218 hAside1D->SetBit(TObject::kCanDelete);
219 hCside1D->SetBit(TObject::kCanDelete);
221 return canv;
222}
223
224template <typename DataT>
225TCanvas* o2::tpc::IDCCCDBHelper<DataT>::drawIDCZeroScale(TCanvas* outputCanvas, bool rejectOutlier) const
227 TCanvas* canv = nullptr;
228
229 if (outputCanvas) {
230 canv = outputCanvas;
231 } else {
232 canv = new TCanvas("c_sides_IDC0_scale", "IDC0", 1000, 1000);
233 }
234 canv->cd();
235 double xsides[] = {0, 1};
236 double yTotIDC0[] = {mScaleIDC0Aside, mScaleIDC0Cside};
237 auto gIDC0Scale = new TGraph(2, xsides, yTotIDC0);
238 gIDC0Scale->SetName("g_IDC0ScaleFactor");
239 gIDC0Scale->SetTitle("Scaling Factor (IDC0);Sides;IDC0 Total (arb. unit)");
240 gIDC0Scale->SetMarkerColor(kBlue);
241 gIDC0Scale->SetMarkerStyle(21);
242 gIDC0Scale->SetMarkerSize(1);
243 gIDC0Scale->GetXaxis()->SetLabelColor(0);
244 gIDC0Scale->GetXaxis()->CenterTitle();
245
246 gIDC0Scale->Draw("ap");
247 // Draw labels on the x axis
248 TText* t = new TText();
249 t->SetTextSize(0.035);
250 const char* labels[2] = {"A", "C"};
251 for (Int_t i = 0; i < 2; i++) {
252 t->DrawText(xsides[i], 0.5, labels[i]);
253 }
254 gIDC0Scale->GetXaxis()->SetLimits(-0.5, 1.5);
255 canv->Update();
256 canv->Modified();
257 gIDC0Scale->SetBit(TObject::kCanDelete);
258 t->SetBit(TObject::kCanDelete);
259 return canv;
260}
261
262template <typename DataT>
263TCanvas* o2::tpc::IDCCCDBHelper<DataT>::drawIDCZeroRadialProfile(TCanvas* outputCanvas, int nbinsY, float yMin, float yMax) const
264{
265 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) {
266 return this->getIDCZeroVal(sector, region, irow, pad);
267 };
268
269 TCanvas* canv = nullptr;
270
271 if (outputCanvas) {
272 canv = outputCanvas;
273 } else {
274 canv = new TCanvas("c_sides_IDC0_radialProfile", "IDC0", 1000, 1000);
275 }
276
278 drawFun.mIDCFunc = idcFunc;
279
280 const auto radialBinning = o2::tpc::painter::getRowBinningCM();
281
282 auto hAside2D = new TH2F("h_IDC0_radialProfile_Aside", "IDC0: Radial profile (A-Side)", radialBinning.size() - 1, radialBinning.data(), nbinsY, yMin, yMax);
283 hAside2D->GetXaxis()->SetTitle("x (cm)");
284 hAside2D->GetYaxis()->SetTitle("IDC0 (ADC)");
285 hAside2D->SetTitleOffset(1.05, "XY");
286 hAside2D->SetTitleSize(0.05, "XY");
287 hAside2D->SetStats(0);
288
289 auto hCside2D = new TH2F("h_IDC0_radialProfile_Cside", "IDC0: Radial profile (C-Side)", radialBinning.size() - 1, radialBinning.data(), nbinsY, yMin, yMax);
290 hCside2D->GetXaxis()->SetTitle("x (cm)");
291 hCside2D->GetYaxis()->SetTitle("IDC0 (ADC)");
292 hCside2D->SetTitleOffset(1.05, "XY");
293 hCside2D->SetTitleSize(0.05, "XY");
294 hCside2D->SetStats(0);
295
296 IDCDrawHelper::drawRadialProfile(drawFun, *hAside2D, Side::A);
297 IDCDrawHelper::drawRadialProfile(drawFun, *hCside2D, Side::C);
298
299 canv->Divide(1, 2);
300 canv->cd(1);
301 hAside2D->Draw("colz");
302 hAside2D->ProfileX("profile_ASide", 1, -1, "d,same");
303 // profA->Draw("d,same");
304 hAside2D->SetStats(0);
305
306 canv->cd(2);
307 hCside2D->Draw("colz");
308 hCside2D->ProfileX("profile_CSide", 1, -1, "d,same");
309 // profC->Draw("d,same");
310 hAside2D->SetStats(0);
311
312 hAside2D->SetBit(TObject::kCanDelete);
313 hCside2D->SetBit(TObject::kCanDelete);
314
315 return canv;
316}
317
318template <typename DataT>
319TCanvas* o2::tpc::IDCCCDBHelper<DataT>::drawIDCZeroStackCanvas(TCanvas* outputCanvas, Side side, std::string_view type, int nbins1D, float xMin1D, float xMax1D, int integrationInterval) const
320{
321 TCanvas* canv = nullptr;
322
323 std::function<float(const unsigned int, const unsigned int, const unsigned int, const unsigned int)> idcFunc;
324
325 if (type == "IDC0") {
326 idcFunc = [this](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
327 return this->getIDCZeroVal(sector, region, irow, pad);
328 };
329 } else if (type == "IDCDelta") {
330 idcFunc = [this, integrationInterval](const unsigned int sector, const unsigned int region, const unsigned int irow, const unsigned int pad) {
331 return this->getIDCDeltaVal(sector, region, irow, pad, integrationInterval);
332 };
333 if (integrationInterval <= 0) {
334 integrationInterval = std::min(getNIntegrationIntervalsIDCDelta(Side::A), getNIntegrationIntervalsIDCDelta(Side::C));
335 }
336 } else {
337 LOGP(error, "Please provide a valid IDC data type. 'IDC0' or 'IDCDelta'.");
338 return canv;
339 }
340
341 if (outputCanvas) {
342 canv = outputCanvas;
343 } else {
344 canv = new TCanvas(fmt::format("c_GEMStacks_{}_1D_{}Side", type.data(), (side == Side::A) ? "A" : "C").data(), fmt::format("{} value 1D distribution for each GEM stack {}-Side", type.data(), (side == Side::A) ? "A" : "C").data(), 1000, 1000);
345 }
346
348 drawFun.mIDCFunc = idcFunc;
349 IDCDrawHelper::drawIDCZeroStackCanvas(drawFun, side, type, nbins1D, xMin1D, xMax1D, *canv, integrationInterval);
350
351 return canv;
352}
353
354template <typename DataT>
355TCanvas* o2::tpc::IDCCCDBHelper<DataT>::drawIDCOneCanvas(TCanvas* outputCanvas, int nbins1D, float xMin1D, float xMax1D, int integrationIntervals) const
356{
357 TCanvas* canv = nullptr;
358
359 if (outputCanvas) {
360 canv = outputCanvas;
361 } else {
362 canv = new TCanvas("c_sides_IDC1_1D", "IDC1 1D distribution for each side", 1000, 1000);
363 }
364
365 auto hAside1D = new TH1F("h_IDC1_1D_ASide", "IDC1 distribution over integration intervals A-Side", nbins1D, xMin1D, xMax1D);
366 auto hCside1D = new TH1F("h_IDC1_1D_CSide", "IDC1 distribution over integration intervals C-Side", nbins1D, xMin1D, xMax1D);
367
368 hAside1D->GetXaxis()->SetTitle("IDC1");
369 hAside1D->SetTitleOffset(1.05, "XY");
370 hAside1D->SetTitleSize(0.05, "XY");
371 hCside1D->GetXaxis()->SetTitle("IDC1");
372 hCside1D->SetTitleOffset(1.05, "XY");
373 hCside1D->SetTitleSize(0.05, "XY");
374
375 if (integrationIntervals <= 0) {
376 integrationIntervals = std::min(mIDCOne[Side::A]->getNIDCs(), mIDCOne[Side::C]->getNIDCs());
377 }
378
379 for (unsigned int integrationInterval = 0; integrationInterval < integrationIntervals; ++integrationInterval) {
380 hAside1D->Fill(getIDCOneVal(Side::A, integrationInterval));
381 hCside1D->Fill(getIDCOneVal(Side::C, integrationInterval));
382 }
383
384 canv->Divide(1, 2);
385 canv->cd(1);
386 hAside1D->Draw();
387 canv->cd(2);
388 hCside1D->Draw();
389
390 hAside1D->SetBit(TObject::kCanDelete);
391 hCside1D->SetBit(TObject::kCanDelete);
392
393 return canv;
394}
395
396template <typename DataT>
397TCanvas* o2::tpc::IDCCCDBHelper<DataT>::drawFourierCoeff(TCanvas* outputCanvas, Side side, int nbins1D, float xMin1D, float xMax1D) const
398{
399 TCanvas* canv = nullptr;
400
401 if (outputCanvas) {
402 canv = outputCanvas;
403 } else {
404 canv = new TCanvas(fmt::format("c_FourierCoefficients_1D_{}Side", (side == Side::A) ? "A" : "C").data(), fmt::format("1D distributions of Fourier Coefficients ({}-Side)", (side == Side::A) ? "A" : "C").data(), 1000, 1000);
405 }
406
407 std::vector<TH1F*> histos;
408
409 for (int i = 0; i < mFourierCoeff[side]->getNCoefficientsPerTF(); i++) {
410 histos.emplace_back(new TH1F(fmt::format("h_FourierCoeff{}_{}Side", i, (side == Side::A) ? "A" : "C").data(), fmt::format("1D distribution of Fourier Coefficient {} ({}-Side)", i, (side == Side::A) ? "A" : "C").data(), nbins1D, xMin1D, xMax1D));
411 histos.back()->GetXaxis()->SetTitle(fmt::format("Fourier Coefficient {}", i).data());
412 histos.back()->SetBit(TObject::kCanDelete);
413 }
414
415 const auto& coeffs = mFourierCoeff[side]->getFourierCoefficients();
416 const auto nCoeffPerTF = mFourierCoeff[side]->getNCoefficientsPerTF();
417
418 for (int i = 0; i < mFourierCoeff[side]->getNCoefficients(); i++) {
419 histos.at(i % nCoeffPerTF)->Fill(coeffs.at(i));
420 }
421
422 canv->DivideSquare(mFourierCoeff[side]->getNCoefficientsPerTF());
423
424 size_t pad = 1;
425
426 for (const auto& hist : histos) {
427 canv->cd(pad);
428 hist->SetTitleOffset(1.05, "XY");
429 hist->SetTitleSize(0.05, "XY");
430 hist->Draw();
431 pad++;
432 }
433
434 return canv;
435}
436
437template <typename DataT>
438void o2::tpc::IDCCCDBHelper<DataT>::dumpToTree(const Side side, const char* outFileName) const
439{
440 const Mapper& mapper = Mapper::instance();
441 o2::utils::TreeStreamRedirector pcstream(outFileName, "RECREATE");
442 pcstream.GetFile()->cd();
443
444 const unsigned int nIDCsSector = Mapper::getPadsInSector() * Mapper::NSECTORS;
445 std::vector<int> vRow(nIDCsSector);
446 std::vector<int> vPad(nIDCsSector);
447 std::vector<float> vXPos(nIDCsSector);
448 std::vector<float> vYPos(nIDCsSector);
449 std::vector<float> vGlobalXPos(nIDCsSector);
450 std::vector<float> vGlobalYPos(nIDCsSector);
451 std::vector<float> idcsZero(nIDCsSector);
452 std::vector<unsigned int> sectorv(nIDCsSector);
453 std::vector<int> outlier(nIDCsSector);
454 std::vector<float> stackMedian(nIDCsSector);
455
456 if (!mIDCZero[side]) {
457 LOGP(info, "IDC0 not set. returning");
458 return;
459 }
460 auto stacksMedianTmp = getStackMedian(*mIDCZero[side], side);
461 const std::array<int, Mapper::NREGIONS> stackInSector{0, 0, 0, 0, 1, 1, 2, 2, 3, 3};
462 unsigned int index = 0;
463 const unsigned int secStart = (side == Side::A) ? 0 : SECTORSPERSIDE;
464 const unsigned int secEnd = (side == Side::A) ? SECTORSPERSIDE : Mapper::NSECTORS;
465 for (unsigned int sector = secStart; sector < secEnd; ++sector) {
466 for (unsigned int region = 0; region < Mapper::NREGIONS; ++region) {
467 for (unsigned int irow = 0; irow < Mapper::ROWSPERREGION[region]; ++irow) {
468 for (unsigned int ipad = 0; ipad < Mapper::PADSPERROW[region][irow]; ++ipad) {
469 const auto padNum = Mapper::getGlobalPadNumber(irow, ipad, region);
470 const auto padTmp = (sector < SECTORSPERSIDE) ? ipad : (Mapper::PADSPERROW[region][irow] - ipad - 1); // C-Side is mirrored
471 const auto& padPosLocal = mapper.padPos(padNum);
472 vRow[index] = padPosLocal.getRow();
473 vPad[index] = padPosLocal.getPad();
474 vXPos[index] = mapper.getPadCentre(padPosLocal).X();
475 vYPos[index] = mapper.getPadCentre(padPosLocal).Y();
476 const GlobalPosition2D globalPos = mapper.LocalToGlobal(LocalPosition2D(vXPos[index], -vYPos[index]), Sector(sector));
477 vGlobalXPos[index] = globalPos.X();
478 vGlobalYPos[index] = globalPos.Y();
479 idcsZero[index] = getIDCZeroVal(sector, region, irow, padTmp);
480 if (mPadFlagsMap) {
481 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][irow] + padTmp;
482 outlier[index] = static_cast<int>(mPadFlagsMap->getCalArray(region + sector * Mapper::NREGIONS).getValue(padInRegion));
483 }
484 stackMedian[index] = stacksMedianTmp[(sector - secStart) * GEMSTACKSPERSECTOR + stackInSector[region]];
485 sectorv[index++] = sector;
486 }
487 }
488 }
489 }
490 std::vector<float> idcOne = !mIDCOne[side] ? std::vector<float>() : mIDCOne[side]->mIDCOne;
491
492 pcstream << "tree"
493 << "IDC0.=" << idcsZero
494 << "IDC1=" << idcOne
495 << "pad.=" << vPad
496 << "row.=" << vRow
497 << "lx.=" << vXPos
498 << "ly.=" << vYPos
499 << "gx.=" << vGlobalXPos
500 << "gy.=" << vGlobalYPos
501 << "outlier.=" << outlier
502 << "sector.=" << sectorv
503 << "stacksMedian=" << stackMedian
504 << "\n";
505 pcstream.Close();
506}
507
508template <typename DataT>
509void o2::tpc::IDCCCDBHelper<DataT>::dumpToTreeIDCDelta(const Side side, const char* outFileName) const
510{
511 o2::utils::TreeStreamRedirector pcstream(outFileName, "RECREATE");
512 pcstream.GetFile()->cd();
513 const unsigned int nIDCsSector = Mapper::getPadsInSector() * SECTORSPERSIDE;
514 std::vector<float> idcsDelta(nIDCsSector);
515 for (int integrationInterval = 0; integrationInterval < getNIntegrationIntervalsIDCDelta(side); ++integrationInterval) {
516 unsigned int index = 0;
517 const unsigned int secStart = (side == Side::A) ? 0 : SECTORSPERSIDE;
518 const unsigned int secEnd = (side == Side::A) ? SECTORSPERSIDE : Mapper::NSECTORS;
519 for (unsigned int sector = secStart; sector < secEnd; ++sector) {
520 for (unsigned int region = 0; region < Mapper::NREGIONS; ++region) {
521 for (unsigned int irow = 0; irow < Mapper::ROWSPERREGION[region]; ++irow) {
522 for (unsigned int ipad = 0; ipad < Mapper::PADSPERROW[region][irow]; ++ipad) {
523 const auto padTmp = (sector < SECTORSPERSIDE) ? ipad : (Mapper::PADSPERROW[region][irow] - ipad - 1); // C-Side is mirrored
524 idcsDelta[index++] = getIDCDeltaVal(sector, region, irow, padTmp, integrationInterval);
525 }
526 }
527 }
528 }
529 float idcOne = getIDCOneVal(side, integrationInterval);
530
531 pcstream << "tree"
532 << "IDC1=" << idcOne
533 << "IDCDelta.=" << idcsDelta
534 << "\n";
535 }
536 pcstream.Close();
537}
538
539template <typename DataT>
541{
542 o2::utils::TreeStreamRedirector pcstream(outFileName, "RECREATE");
543 pcstream.GetFile()->cd();
544
545 for (int iside = 0; iside < SIDES; ++iside) {
546 const Side side = (iside == 0) ? Side::A : Side::C;
547
548 if (!mFourierCoeff[side]) {
549 continue;
550 }
551
552 const int nTFs = mFourierCoeff[side]->getNCoefficients() / mFourierCoeff[side]->getNCoefficientsPerTF();
553 for (int iTF = 0; iTF < nTFs; ++iTF) {
554 std::vector<float> coeff;
555 std::vector<int> ind;
556 int coeffPerTF = mFourierCoeff[side]->getNCoefficientsPerTF();
557 for (int i = 0; i < coeffPerTF; ++i) {
558 const int index = mFourierCoeff[side]->getIndex(iTF, i);
559 coeff.emplace_back((*mFourierCoeff[side])(index));
560 ind.emplace_back(i);
561 }
562
563 pcstream << "tree"
564 << "iTF=" << iTF
565 << "index=" << ind
566 << "coeffPerTF=" << coeffPerTF
567 << "coeff.=" << coeff
568 << "side=" << iside
569 << "\n";
570 }
571 }
572 pcstream.Close();
573}
574
575template <typename DataT>
577{
578 CalDet<float> calIDC0("IDC0");
579 for (unsigned int cru = 0; cru < CRU::MaxCRU; ++cru) {
580 const o2::tpc::CRU cruTmp(cru);
581 const Side side = cruTmp.side();
582 const int region = cruTmp.region();
583 const int sector = cruTmp.sector();
584 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
585 const unsigned int integrationInterval = 0;
586 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
587 const auto idcZero = getIDCZeroVal(sector, region, lrow, pad);
588 calIDC0.setValue(sector, Mapper::ROWOFFSET[region] + lrow, pad, idcZero);
589 }
590 }
591 }
592 return calIDC0;
593}
594
595template <typename DataT>
596std::vector<o2::tpc::CalDet<float>> o2::tpc::IDCCCDBHelper<DataT>::getIDCDeltaCalDet() const
597{
598 const unsigned int nIntervalsA = getNIntegrationIntervalsIDCDelta(Side::A);
599 const unsigned int nIntervalsC = getNIntegrationIntervalsIDCDelta(Side::C);
600 if (nIntervalsA != nIntervalsC) {
601 LOGP(info, "Number of integration interval for A Side {} unequal for C side {}", nIntervalsA, nIntervalsC);
602 }
603
604 std::vector<CalPad> calIDCDelta(std::max(nIntervalsA, nIntervalsC));
605 for (auto& calpadIDC : calIDCDelta) {
606 calpadIDC = CalPad("IDCDelta", PadSubset::ROC);
607 }
608
609 for (unsigned int cru = 0; cru < CRU::MaxCRU; ++cru) {
610 const o2::tpc::CRU cruTmp(cru);
611 const Side side = cruTmp.side();
612 const int region = cruTmp.region();
613 const int sector = cruTmp.sector();
614
615 const unsigned int nIntervals = getNIntegrationIntervalsIDCDelta(side);
616 for (unsigned int integrationInterval = 0; integrationInterval < nIntervals; ++integrationInterval) {
617 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
618 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
619 const auto idcdelta = getIDCDeltaVal(sector, region, lrow, pad, integrationInterval);
620 calIDCDelta[integrationInterval].setValue(sector, Mapper::ROWOFFSET[region] + lrow, pad, idcdelta);
621 }
622 }
623 }
624 }
625 return calIDCDelta;
626}
627
628template <typename DataT>
630{
631 if (!mIDCZero[Side::A] && !mIDCZero[Side::C]) {
632 LOGP(info, "IDC0 not set. returning");
633 }
634
635 int nCRUS = 0;
636 if (mIDCZero[Side::A]) {
637 nCRUS += CRU::MaxCRU / 2;
638 }
639
640 if (mIDCZero[Side::C]) {
641 nCRUS += CRU::MaxCRU / 2;
642 }
643
644 std::vector<uint32_t> crus(nCRUS);
645 std::iota(crus.begin(), crus.end(), mIDCZero[Side::A] ? 0 : (CRU::MaxCRU / 2));
646
647 IDCFactorization idc(1, 1, crus);
648 if (mIDCZero[Side::A]) {
649 idc.setIDCZero(Side::A, *mIDCZero[Side::A]);
650 }
651 if (mIDCZero[Side::C]) {
652 idc.setIDCZero(Side::C, *mIDCZero[Side::C]);
653 }
654 idc.createStatusMap();
655 mPadFlagsMap = idc.getPadStatusMap();
656}
657
658template <typename DataT>
660{
661 if (!mIDCZero[side]) {
662 LOGP(info, "IDC0 not set. returning");
663 return;
664 }
665
666 *mIDCZero[side] /= factor;
667}
668
669template <typename DataT>
671{
672 if (rejectOutlier) {
673 createOutlierMap();
674 }
675
676 if (mIDCZero[Side::A]) {
677 mScaleIDC0Aside = getMeanIDC0(Side::A, *mIDCZero[Side::A], rejectOutlier ? mPadFlagsMap.get() : nullptr);
678 scaleIDC0(mScaleIDC0Aside, Side::A);
679 } else {
680 mScaleIDC0Aside = 0;
681 }
682
683 if (mIDCZero[Side::C]) {
684 mScaleIDC0Cside = getMeanIDC0(Side::C, *mIDCZero[Side::C], rejectOutlier ? mPadFlagsMap.get() : nullptr);
685 scaleIDC0(mScaleIDC0Cside, Side::C);
686 } else {
687 mScaleIDC0Cside = 0;
688 }
689
691 if (mScaleIDC0Aside == 0.0 || mScaleIDC0Cside == 0.0) {
692 LOGP(error, "Please check the IDC0 total A side {} and C side {}, is zero, therefore no scaling applied!", mScaleIDC0Aside, mScaleIDC0Cside);
693 mScaleIDC0Aside = 1.0;
694 mScaleIDC0Cside = 1.0;
695 }
696}
697
698template <typename DataT>
699float o2::tpc::IDCCCDBHelper<DataT>::getMeanIDC0(const Side side, const IDCZero& idcZero, const CalDet<PadFlags>* outlierMap)
700{
701 // 1. calculate mean per row
702 std::vector<double> idc0Sum(Mapper::PADROWS);
703 std::vector<unsigned int> idcChannelsCount(Mapper::PADROWS);
704 const unsigned int firstCRU = (side == Side::A) ? 0 : CRU::MaxCRU / 2;
705 for (unsigned int cru = firstCRU; cru < firstCRU + CRU::MaxCRU / 2; ++cru) {
706 const o2::tpc::CRU cruTmp(cru);
707 const int region = cruTmp.region();
708 const int sector = cruTmp.sector();
709 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
710 const unsigned int glRow = lrow + Mapper::ROWOFFSET[region];
711 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
712 const unsigned int padInRegion = Mapper::OFFSETCRULOCAL[region][lrow] + pad;
713 if (outlierMap) {
714 const o2::tpc::PadFlags flag = outlierMap->getCalArray(cru).getValue(padInRegion);
715 if ((flag & PadFlags::flagSkip) == PadFlags::flagSkip) {
716 continue;
717 }
718 }
719 const auto index = getUngroupedIndexGlobal(sector, region, lrow, pad, 0);
720 idc0Sum[glRow] += idcZero.getValueIDCZero(index);
721 ++idcChannelsCount[glRow];
722 }
723 }
724 }
725
726 // 2. calculate mean of means
727 double idc0SumGlobal = 0;
728 unsigned int idc0CountGlobal = 0;
729 for (unsigned int i = 0; i < Mapper::PADROWS; ++i) {
730 if (idcChannelsCount[i]) {
731 idc0SumGlobal += idc0Sum[i] / idcChannelsCount[i];
732 ++idc0CountGlobal;
733 }
734 }
735
736 if (idc0CountGlobal == 0) {
737 return -1;
738 }
739
740 return idc0SumGlobal /= idc0CountGlobal;
741}
742
743template <typename DataT>
745{
746 int cruStart = sector * CRU::CRUperSector;
747 if (stack == GEMstack::OROC1gem) {
748 cruStart += CRU::CRUperIROC;
749 } else if (stack == GEMstack::OROC2gem) {
750 cruStart += CRU::CRUperIROC + CRU::CRUperPartition;
751 } else if (stack == GEMstack::OROC3gem) {
752 cruStart += CRU::CRUperIROC + 2 * CRU::CRUperPartition;
753 }
754 int cruEnd = (stack == GEMstack::IROCgem) ? (cruStart + CRU::CRUperIROC) : (cruStart + CRU::CRUperPartition);
755
756 RobustAverage average;
757 average.reserve(Mapper::getNumberOfPads(stack));
758 for (int cru = cruStart; cru < cruEnd; ++cru) {
759 const o2::tpc::CRU cruTmp(cru);
760 const int region = cruTmp.region();
761 for (unsigned int lrow = 0; lrow < Mapper::ROWSPERREGION[region]; ++lrow) {
762 const unsigned int integrationInterval = 0;
763 for (unsigned int pad = 0; pad < Mapper::PADSPERROW[region][lrow]; ++pad) {
764 const auto index = getUngroupedIndexGlobal(sector, region, lrow, pad, integrationInterval);
765 const float idcZero = idczero.mIDCZero[index];
766 if ((idcZero != -1) && (idcZero != 0)) {
767 average.addValue(idcZero);
768 }
769 }
770 }
771 }
772 return average.getMedian();
773}
774
775template <typename DataT>
776std::array<float, o2::tpc::GEMSTACKSPERSECTOR * o2::tpc::SECTORSPERSIDE> o2::tpc::IDCCCDBHelper<DataT>::getStackMedian(const IDCZero& idczero, const Side side)
777{
778 const int stacksPerSide = o2::tpc::GEMSTACKSPERSECTOR * o2::tpc::SECTORSPERSIDE;
779 std::array<float, stacksPerSide> median;
780 unsigned int sectorStart = (side == Side::A) ? 0 : o2::tpc::SECTORSPERSIDE;
781 unsigned int sectorEnd = (side == Side::A) ? o2::tpc::SECTORSPERSIDE : Mapper::NSECTORS;
782 for (unsigned int sector = sectorStart; sector < sectorEnd; ++sector) {
783 for (int stackInSector = 0; stackInSector < GEMSTACKSPERSECTOR; ++stackInSector) {
784 median[(sector * GEMSTACKSPERSECTOR + stackInSector) % stacksPerSide] = IDCCCDBHelper<>::getStackMedian(idczero, Sector(sector), static_cast<GEMstack>(stackInSector));
785 }
786 }
787 return median;
788}
789
790template <typename DataT>
792{
793 std::pair<int, int> nOutliers{};
794 if (!mPadFlagsMap) {
795 LOGP(info, "Status map not set returning");
796 return nOutliers;
797 }
798 for (CRU cru; !cru.looped(); ++cru) {
799 const auto& calArray = mPadFlagsMap->getCalArray(cru);
800 for (auto flag : calArray.getData()) {
801 if ((flag & PadFlags::flagSkip) == PadFlags::flagSkip) {
802 (cru.side() == Side::A) ? ++nOutliers.first : ++nOutliers.second;
803 }
804 }
805 }
806 return nOutliers;
807}
808
809template class o2::tpc::IDCCCDBHelper<float>;
int32_t i
float & yMax
helper class for accessing IDCs from CCDB
This file provides the structs for storing the factorized IDC values and fourier coefficients to be s...
helper class for drawing IDCs per region/side
class for aggregating IDCs for the full TPC (all sectors) and factorization of aggregated IDCs
helper class for grouping of pads and rows for one sector
uint32_t side
Definition RawData.h:0
uint32_t stack
Definition RawData.h:1
class for performing robust averaging and outlier filtering
unsigned char region() const
Definition CRU.h:64
Side side() const
Definition CRU.h:62
bool looped() const
Definition CRU.h:75
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
void setValue(const int sec, const int globalPadInSector, const T &value)
Definition CalDet.h:257
CalDet< float > getIDCZeroCalDet() const
TCanvas * drawIDCOneCanvas(TCanvas *outputCanvas, int nbins1D, float xMin1D, float xMax1D, int integrationIntervals=-1) const
void scaleIDC0(const float factor, const Side side)
void dumpToFourierCoeffToTree(const char *outFileName="FourierCCDBTree.root") const
float getIDCZeroVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad) const
void setIDCZeroScale(const bool rejectOutlier=true)
TCanvas * drawIDCZeroScale(TCanvas *outputCanvas, const bool rejectOutlier=true) const
void createOutlierMap()
create the outlier map with the set unscaled IDC0 map
std::vector< CalDet< float > > getIDCDeltaCalDet() const
TCanvas * drawIDCZeroCanvas(TCanvas *outputCanvas, std::string_view type, int nbins1D, float xMin1D, float xMax1D, int integrationInterval=-1) const
TCanvas * drawFourierCoeff(TCanvas *outputCanvas, Side side, int nbins1D, float xMin1D, float xMax1D) const
std::pair< int, int > getNOutliers() const
float getIDCVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const
static float getMeanIDC0(const Side side, const IDCZero &idcZero, const CalDet< PadFlags > *outlierMap)
unsigned int getNIntegrationIntervalsIDCDelta(const o2::tpc::Side side) const
float getIDCOneVal(const o2::tpc::Side side, const unsigned int integrationInterval) const
void dumpToTree(const Side side, const char *outFileName="IDCCCDBTree.root") const
static float getStackMedian(const IDCZero &idczero, const Sector sector, const GEMstack stack)
unsigned int getNIntegrationIntervalsIDCOne(const o2::tpc::Side side) const
TCanvas * drawIDCZeroStackCanvas(TCanvas *outputCanvas, Side side, std::string_view type, int nbins1D, float xMin1D, float xMax1D, int integrationInterval=-1) const
TCanvas * drawIDCZeroRadialProfile(TCanvas *outputCanvas, int nbinsY, float yMin, float yMax) const
float getIDCDeltaVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const
void dumpToTreeIDCDelta(const Side side, const char *outFileName="IDCCCDBTreeDeltaIDC.root") const
void setIDCZero(const Side side, const IDCZero &idcZero)
set the IDCZero object
std::unique_ptr< CalDet< PadFlags > > getPadStatusMap()
const PadPos & padPos(GlobalPadNumber padNumber) const
Definition Mapper.h:50
static GlobalPosition3D LocalToGlobal(const LocalPosition3D &pos, const double alpha)
Definition Mapper.h:461
GlobalPosition2D getPadCentre(const PadSecPos &padSec) const
Definition Mapper.h:163
void addValue(const float value, const float weight=1)
void reserve(const unsigned int maxValues)
Side side() const
Definition Sector.h:96
GLuint index
Definition glcorearb.h:781
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLboolean * data
Definition glcorearb.h:298
GLint GLenum GLint const GLfloat * coeffs
Definition glcorearb.h:5520
GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat maxZ
Definition glcorearb.h:2910
GLfloat GLfloat minZ
Definition glcorearb.h:2910
GEMstack
TPC GEM stack types.
Definition Defs.h:53
constexpr unsigned char SECTORSPERSIDE
Definition Defs.h:40
constexpr unsigned char SIDES
Definition Defs.h:41
Side
TPC readout sidE.
Definition Defs.h:35
PadFlags
Definition Defs.h:100
constexpr unsigned short GEMSTACKSPERSECTOR
Definition Defs.h:57
std::string filename()
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
struct containing the ITPC0 values (integrated TPC clusters)
std::vector< float > mIDCZero
I_0(r,\phi) = <I(r,\phi,t)>_t.
float getValueIDCZero(const unsigned int index) const
static std::vector< double > getRowBinningCM(uint32_t roc=72)
static void drawSectorsXY(Side side, int sectorLineColor=920, int sectorTextColor=1)
draw sector boundaris, side name and sector numbers
std::vector< int > row