Project
Loading...
Searching...
No Matches
CheckResidSpec.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
15#include <vector>
17#include <TStopwatch.h>
35#include "ITStracking/IOUtils.h"
43#include <TROOT.h>
44#include <TStyle.h>
45#include <TLatex.h>
46#include <TCanvas.h>
47#include <TLegend.h>
48#include <TLegendEntry.h>
49#include <TH1F.h>
50#include <TH2F.h>
51#include <TProfile.h>
52#include <TGraph.h>
53#include <TF1.h>
54#ifdef WITH_OPENMP
55#include <omp.h>
56#endif
57
58// Attention: in case the residuals are checked with geometry different from the one used for initial reconstruction,
59// pass a --configKeyValues option for vertex refit as:
60// ;pvertexer.useMeanVertexConstraint=false;pvertexer.meanVertexExtraErrSelection=0.2;pvertexer.iniScale2=100;pvertexer.acceptableScale2=10.;
61// In any case, it is better to pass ;pvertexer.useMeanVertexConstraint=false;
62
63namespace o2::checkresid
64{
65using namespace o2::framework;
68
74
75class CheckResidSpec : public Task
76{
77 public:
78 CheckResidSpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, GTrackID::mask_t src, bool drawOnly, bool postProcOnly)
79 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(src), mDrawOnly(drawOnly), mPostProcOnly(postProcOnly)
80 {
81 }
82 ~CheckResidSpec() final = default;
83 void init(InitContext& ic) final;
84 void run(ProcessingContext& pc) final;
85 void endOfStream(EndOfStreamContext& ec) final;
86 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
87 void process();
88
89 private:
90 void updateTimeDependentParams(ProcessingContext& pc);
91 bool refitPV(o2::dataformats::PrimaryVertex& pv, int vid);
92 bool refitITStrack(o2::track::TrackParCov& track, GTrackID gid);
93 bool processITSTrack(const o2::its::TrackITS& iTrack, const o2::dataformats::PrimaryVertex& pv, o2::checkresid::Track& resTrack);
94 void bookHistos();
95 void fillHistos(const o2::checkresid::Track& trc);
96 void postProcessHistos();
97 void drawHistos();
98
99 o2::globaltracking::RecoContainer* mRecoData = nullptr;
100 int mNThreads = 1;
101 bool mMeanVertexUpdated = false;
102 float mITSROFrameLengthMUS = 0.f;
103 o2::dataformats::MeanVertexObject mMeanVtx{};
104 std::vector<o2::BaseCluster<float>> mITSClustersArray;
105 const o2::itsmft::TopologyDictionary* mITSDict = nullptr;
106 o2::vertexing::PVertexer mVertexer;
107 std::shared_ptr<DataRequest> mDataRequest;
108 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
109 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
110 GTrackID::mask_t mTracksSrc{};
111
112 bool mDrawOnly = false;
113 bool mPostProcOnly = false;
114 bool mDraw = false;
115 bool mFillHistos = true;
116 bool mFillTree = true;
117 std::vector<std::unique_ptr<o2::HistoManager>> mHManV{};
118 o2::HistoManager* mHMan = nullptr;
119};
120
122{
123 mDraw = true;
124 if (!mDrawOnly) {
125 mDraw = ic.options().get<bool>("draw-report");
126 mFillHistos = !ic.options().get<bool>("no-hist");
127 mFillTree = !ic.options().get<bool>("no-tree");
128 mNThreads = ic.options().get<int>("nthreads");
129 }
131 int lane = ic.services().get<const o2::framework::DeviceSpec>().inputTimesliceId;
132 int maxLanes = ic.services().get<const o2::framework::DeviceSpec>().maxInputTimeslices;
133 std::string nm = params.outname;
134 if (maxLanes > 1) {
135 o2::conf::ConfigurableParam::updateFromString(fmt::format("checkresid.outname={}_{}", nm, lane));
136 }
137 if (mDraw) {
138 mFillHistos = true;
139 }
140 if (!mDrawOnly && mFillHistos) {
141 bookHistos();
142 }
143 if (!params.ext_hm_list.empty()) {
144 auto vecNames = o2::utils::Str::tokenize(params.ext_hm_list, ',');
145 auto vecLegends = o2::utils::Str::tokenize(params.ext_leg_list, ',');
146 bool useLeg = true;
147 if (vecNames.size() != vecLegends.size()) {
148 LOGP(warn, "{} legend names provided for {} external histomanagers, will use file names as legends", vecLegends.size(), vecNames.size());
149 useLeg = false;
150 }
151 int cntH = 0;
152 for (const auto& vn : vecNames) {
153 LOGP(info, "Loading external HistoManager {}", vn);
154 mHManV.emplace_back() = std::make_unique<o2::HistoManager>("", vn, true);
155 auto hm = mHManV.back().get();
156 if (!hm) {
157 LOGP(error, "Failed to load histograms from {}", vn);
158 mHManV.pop_back();
159 } else {
160 hm->SetName(useLeg ? vecLegends[cntH].c_str() : vn.c_str());
161 }
162 cntH++;
163 }
164 }
165 if (mDrawOnly) {
166 return;
167 }
169#ifndef WITH_OPENMP
170 if (mNThreads > 1) {
171 LOGP(warn, "No OpenMP");
172 }
173 mNThreads = 1;
174#endif
175 if (mFillTree) {
176 nm += ".root";
177 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>(nm.c_str(), "recreate");
178 }
179}
180
182{
183 bool quit = false;
184 if (mPostProcOnly) {
185
186 postProcessHistos();
187 quit = true;
188 }
189 if (mDrawOnly) {
190 drawHistos();
191 quit = true;
192 }
193 if (quit) {
195 pc.services().get<ControlService>().readyToQuit(QuitRequest::Me);
196 return;
197 }
199 mRecoData = &recoData;
200 mRecoData->collectData(pc, *mDataRequest.get()); // select tracks of needed type, with minimal cuts, the real selected will be done in the vertexer
201 mRecoData = &recoData;
202 updateTimeDependentParams(pc); // Make sure this is called after recoData.collectData, which may load some conditions
203 process();
204 mRecoData = nullptr;
205}
206
207void CheckResidSpec::updateTimeDependentParams(ProcessingContext& pc)
208{
211 static bool initOnceDone = false;
212 if (!initOnceDone) { // this params need to be queried only once
214 initOnceDone = true;
215 // Note: reading of the ITS AlpideParam needed for ITS timing is done by the RecoContainer
218 if (!grp->isDetContinuousReadOut(DetID::ITS)) {
219 mITSROFrameLengthMUS = alpParams.roFrameLengthTrig / 1.e3; // ITS ROFrame duration in \mus
220 } else {
221 mITSROFrameLengthMUS = alpParams.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingNS * 1e-3; // ITS ROFrame duration in \mus
222 }
225 o2::conf::ConfigurableParam::updateFromString("pvertexer.useTimeInChi2=false;");
226 mVertexer.init();
227 }
228 if (mMeanVertexUpdated) {
229 mMeanVertexUpdated = false;
230 mVertexer.setMeanVertex(&mMeanVtx);
231 mVertexer.initMeanVertexConstraint();
232 }
233}
234
236{
237 if (!mITSDict) {
238 LOGP(fatal, "ITS data is not loaded");
239 }
240 const auto itsTracks = mRecoData->getITSTracks();
241 // const auto itsLbls = mRecoData->getITSTracksMCLabels();
242 const auto itsClRefs = mRecoData->getITSTracksClusterRefs();
243 const auto clusITS = mRecoData->getITSClusters();
244 const auto patterns = mRecoData->getITSClustersPatterns();
246 auto pattIt = patterns.begin();
247 mITSClustersArray.clear();
248 mITSClustersArray.reserve(clusITS.size());
249
250 o2::its::ioutils::convertCompactClusters(clusITS, pattIt, mITSClustersArray, mITSDict);
251
252 auto pvvec = mRecoData->getPrimaryVertices();
253 auto trackIndex = mRecoData->getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
254 auto vtxRefs = mRecoData->getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
255 auto prop = o2::base::Propagator::Instance();
256 static int TFCount = 0;
257 int nv = vtxRefs.size() - 1;
258 std::vector<std::vector<checkresid::Track>> slots;
259 slots.resize(mNThreads);
260 int nvGood = 0, nvUse = 0, nvRefFail = 0;
261 long pvFitDuration{};
262 for (int iv = 0; iv < nv; iv++) {
263 const auto& vtref = vtxRefs[iv];
264 auto pve = pvvec[iv];
265 if (pve.getNContributors() < params.minPVContributors) {
266 continue;
267 }
268 nvGood++;
269 if (params.refitPV) {
270 LOGP(debug, "Refitting PV#{} of {} tracks", iv, pve.getNContributors());
271 auto tStartPVF = std::chrono::time_point_cast<std::chrono::microseconds>(std::chrono::system_clock::now()).time_since_epoch().count();
272 bool res = refitPV(pve, iv);
273 pvFitDuration += std::chrono::time_point_cast<std::chrono::microseconds>(std::chrono::system_clock::now()).time_since_epoch().count() - tStartPVF;
274 if (!res) {
275 nvRefFail++;
276 continue;
277 }
278 }
279 nvUse++;
280 for (int is = 0; is < GTrackID::NSources; is++) {
281 if (!mTracksSrc[is] || !mRecoData->isTrackSourceLoaded(is)) {
282 continue;
283 }
284 int idMin = vtref.getFirstEntryOfSource(is), idMax = idMin + vtref.getEntriesOfSource(is);
285 DetID::mask_t dm = GTrackID::getSourceDetectorsMask(is);
286 if (!dm[DetID::ITS]) {
287 continue;
288 }
289 if (dm[DetID::TPC] && params.minTPCCl > 0 && !mRecoData->isTrackSourceLoaded(GTrackID::TPC)) {
290 LOGP(fatal, "Cut on TPC tracks is requested by they are not loaded");
291 }
292#ifdef WITH_OPENMP
293#pragma omp parallel for schedule(dynamic) num_threads(mNThreads)
294#endif
295 for (int i = idMin; i < idMax; i++) {
296 auto vid = trackIndex[i];
297 bool pvCont = vid.isPVContributor();
298 if (!pvCont && params.pvcontribOnly) {
299 continue;
300 }
301 if (dm[DetID::TPC] && params.minTPCCl > 0 && mRecoData->getTPCTrack(mRecoData->getTPCContributorGID(vid)).getNClusters() < params.minTPCCl) {
302 continue;
303 }
304 auto gidITS = mRecoData->getITSContributorGID(vid);
305 if (gidITS.getSource() != GTrackID::ITS) {
306 continue;
307 }
308 const auto& trc = mRecoData->getTrackParam(vid);
309 const auto& itsTrack = mRecoData->getITSTrack(gidITS);
310 if (itsTrack.getNClusters() < params.minITSCl) {
311 continue;
312 }
313 auto pt = trc.getPt();
314 if (pt < params.minPt || pt > params.maxPt) {
315 continue;
316 }
317 if (std::abs(trc.getTgl()) > params.maxTgl) {
318 continue;
319 }
320
321#ifdef WITH_OPENMP
322 auto& accum = slots[omp_get_thread_num()];
323#else
324 auto& accum = slots[0];
325#endif
326 auto& resTrack = accum.emplace_back();
327 resTrack.gid = vid;
328 if (!processITSTrack(itsTrack, pve, resTrack)) {
329 accum.pop_back();
330 continue;
331 }
332 }
333 }
334 }
335 // output
336 for (const auto& accum : slots) {
337 for (const auto& tr : accum) {
338 if (mDBGOut) {
339 (*mDBGOut) << "res" << "tr=" << tr << "\n";
340 }
341 if (mHMan) {
342 fillHistos(tr);
343 }
344 }
345 }
346 LOGP(info, "processed {} PVs out of {} good vertices (out of {} in total), PV refits took {} mus, {} refits failed", nvUse, nvGood, nv, pvFitDuration, nvRefFail);
347 TFCount++;
348}
349
350bool CheckResidSpec::processITSTrack(const o2::its::TrackITS& iTrack, const o2::dataformats::PrimaryVertex& pv, o2::checkresid::Track& resTrack)
351{
352 const auto itsClRefs = mRecoData->getITSTracksClusterRefs();
353 auto trFitInw = iTrack.getParamOut(); // seed for inward refit
354 auto trFitOut = iTrack.getParamIn(); // seed for outward refit
355 auto prop = o2::base::Propagator::Instance();
357 float pvAlpha = 0;
358 float bz = prop->getNominalBz();
359 std::array<const o2::BaseCluster<float>*, 8> clArr{};
360 const auto& params = CheckResidConfig::Instance();
361 std::array<o2::track::TrackParCov, 8> extrapOut, extrapInw; // 2-way Kalman extrapolations, vertex + 7 layers
362
363 auto rotateTrack = [bz](o2::track::TrackParCov& tr, float alpha, o2::track::TrackPar* refLin) {
364 return refLin ? tr.rotate(alpha, *refLin, bz) : tr.rotate(alpha);
365 };
366
367 auto accountCluster = [&](int i, std::array<o2::track::TrackParCov, 8>& extrapDest, o2::track::TrackParCov& tr, o2::track::TrackPar* refLin) {
368 if (clArr[i]) { // update with cluster
369 if (!rotateTrack(tr, i == 0 ? pvAlpha : geom->getSensorRefAlpha(clArr[i]->getSensorID()), refLin) ||
370 !prop->propagateTo(tr, refLin, clArr[i]->getX(), true)) {
371 return 0;
372 }
373 extrapDest[i] = tr; // before update
374 if (!tr.update(*clArr[i])) {
375 return 0;
376 }
377 } else {
378 extrapDest[i].invalidate();
379 return -1;
380 }
381 return 1;
382 };
383
384 auto inv2d = [](float s00, float s11, float s01) -> std::array<float, 3> {
385 auto det = s00 * s11 - s01 * s01;
386 if (det < 1e-16) {
387 LOGP(error, "Singular det {}, input: {} {} {}", det, s00, s11, s01);
388 return {0.f, 0.f, 0.f};
389 }
390 det = 1.f / det;
391 return {s11 * det, s00 * det, -s01 * det};
392 };
393
394 resTrack.points.clear();
395 if (!prop->propagateToDCA(pv, trFitOut, bz)) {
396 LOGP(debug, "Failed to propagateToDCA, {}", trFitOut.asString());
397 return false;
398 }
400 if (params.addPVAsCluster) {
401 float cosAlp, sinAlp;
402 pvAlpha = trFitOut.getAlpha();
403 o2::math_utils::sincos(trFitOut.getAlpha(), sinAlp, cosAlp); // vertex position rotated to track frame
404 bcPV.setXYZ(pv.getX() * cosAlp + pv.getY() * sinAlp, -pv.getX() * sinAlp + pv.getY() * cosAlp, pv.getZ());
405 bcPV.setSigmaY2(0.5 * (pv.getSigmaX2() + pv.getSigmaY2()));
406 bcPV.setSigmaZ2(pv.getSigmaZ2());
407 bcPV.setSensorID(-1);
408 clArr[0] = &bcPV;
409 }
410 // collect all track clusters to array, placing them to layer+1 slot
411 int nCl = iTrack.getNClusters();
412 for (int i = 0; i < nCl; i++) { // clusters are ordered from the outermost to the innermost
413 const auto& curClu = mITSClustersArray[itsClRefs[iTrack.getClusterEntry(i)]];
414
415 int llr = geom->getLayer(curClu.getSensorID());
416 if (clArr[1 + llr]) {
417 LOGP(error, "Cluster at lr {} was already assigned, old sens {}, new sens {}", llr, clArr[1 + llr]->getSensorID(), curClu.getSensorID());
418 }
419 clArr[1 + geom->getLayer(curClu.getSensorID())] = &curClu;
420 }
421 o2::track::TrackPar refLinInw0, refLinOut0, *refLinOut = nullptr, *refLinInw = nullptr;
422 o2::track::TrackPar refLinIBOut0, refLinOBInw0, *refLinOBInw = nullptr, *refLinIBOut = nullptr;
423 if (params.useStableRef) {
424 refLinOut = &(refLinOut0 = trFitOut);
425 refLinInw = &(refLinInw0 = trFitInw);
426 }
427 trFitOut.resetCovariance();
428 trFitOut.setCov(trFitOut.getQ2Pt() * trFitOut.getQ2Pt() * trFitOut.getCov()[14], 14);
429 trFitInw.resetCovariance();
430 trFitInw.setCov(trFitInw.getQ2Pt() * trFitInw.getQ2Pt() * trFitInw.getCov()[14], 14);
431 // fit in inward and outward direction
432 for (int i = 0; i <= 7; i++) {
433 int resOut, resInw;
434 // process resOut in ascending order (0-->7) and resInw in descending order (7-->0)
435 if (!(resOut = accountCluster(i, extrapOut, trFitOut, refLinOut)) || !(resInw = accountCluster(7 - i, extrapInw, trFitInw, refLinInw))) {
436 return false;
437 }
438 // at layer 3, find the IB track (trIBOut) and the OB track (trOBInw)
439 // propagate both trcaks to a common radius, RCompIBOB (12cm), and rotates
440 // them to the same reference frame for comparison
441 if (i == 3 && resOut == 1 && resInw == 1 && params.doIBOB && nCl == 7) {
442 resTrack.trIBOut = trFitOut; // outward track updated at outermost IB layer
443 resTrack.trOBInw = trFitInw; // inward track updated at innermost OB layer
444 o2::track::TrackPar refLinIBOut0, refLinIBIn0;
445 if (refLinOut) {
446 refLinIBOut = &(refLinIBOut0 = refLinOut0);
447 refLinOBInw = &(refLinOBInw0 = refLinInw0);
448 }
449 float xRref;
450 if (!resTrack.trOBInw.getXatLabR(params.rCompIBOB, xRref, bz) ||
451 !prop->propagateTo(resTrack.trOBInw, refLinOBInw, xRref, true) ||
452 !rotateTrack(resTrack.trOBInw, resTrack.trOBInw.getPhiPos(), refLinOBInw) || // propagate OB track to ref R and rotate
453 !rotateTrack(resTrack.trIBOut, resTrack.trOBInw.getAlpha(), refLinIBOut) ||
454 !prop->propagateTo(resTrack.trIBOut, refLinIBOut, resTrack.trOBInw.getX(), true)) { // rotate OB track to same frame and propagate to same X
455 // if any propagation or rotation steps fail, invalidate both tracks
456 return false;
457 }
458 }
459 }
460
461 bool innerDone = false;
462 if (params.doResid) {
463 for (int i = 0; i <= 7; i++) {
464 if (clArr[i]) {
465 // calculate interpolation as a weighted mean of inward/outward extrapolations to this layer
466 const auto &tInw = extrapInw[i], &tOut = extrapOut[i];
467 auto wInw = inv2d(tInw.getSigmaY2(), tInw.getSigmaZ2(), tInw.getSigmaZY());
468 auto wOut = inv2d(tOut.getSigmaY2(), tOut.getSigmaZ2(), tOut.getSigmaZY());
469 if (wInw[0] == 0.f || wOut[0] == 0.f) {
470 return false;
471 }
472 std::array<float, 3> wTot = {wInw[0] + wOut[0], wInw[1] + wOut[1], wInw[2] + wOut[2]};
473 auto cTot = inv2d(wTot[0], wTot[1], wTot[2]);
474 auto ywi = wInw[0] * tInw.getY() + wInw[2] * tInw.getZ() + wOut[0] * tOut.getY() + wOut[2] * tOut.getZ();
475 auto zwi = wInw[2] * tInw.getY() + wInw[1] * tInw.getZ() + wOut[2] * tOut.getY() + wOut[1] * tOut.getZ();
476 auto yw = ywi * cTot[0] + zwi * cTot[2];
477 auto zw = ywi * cTot[2] + zwi * cTot[1];
478 // posCl.push_back(clArr[i]->getXYZGlo(*o2::its::GeometryTGeo::Instance()));
479 auto phi = i == 0 ? tInw.getPhi() : tInw.getPhiPos();
480 o2::math_utils::bringTo02Pi(phi);
481 resTrack.points.emplace_back(clArr[i]->getY() - yw, clArr[i]->getZ() - zw, cTot[0] + clArr[i]->getSigmaY2(), cTot[1] + clArr[i]->getSigmaZ2(), phi, clArr[i]->getZ(), clArr[i]->getSensorID(), i - 1);
482 if (!innerDone) {
483 resTrack.track = tInw;
484 innerDone = true;
485 }
486 } else {
487 LOGP(debug, "No cluster on lr {}", i);
488 }
489 }
490 }
491 return true;
492}
493
494bool CheckResidSpec::refitPV(o2::dataformats::PrimaryVertex& pv, int vid)
495{
497 std::vector<o2::track::TrackParCov> tracks;
498 std::vector<bool> useTrack;
499 std::vector<GTrackID> gidsITS;
500 int ntr = pv.getNContributors(), ntrIni = ntr;
501 tracks.reserve(ntr);
502 useTrack.reserve(ntr);
503 gidsITS.reserve(ntr);
504 const auto& vtref = mRecoData->getPrimaryVertexMatchedTrackRefs()[vid];
505 auto trackIndex = mRecoData->getPrimaryVertexMatchedTracks();
506 int itr = vtref.getFirstEntry(), itLim = itr + vtref.getEntries();
507 for (; itr < itLim; itr++) {
508 auto tid = trackIndex[itr];
509 if (tid.isPVContributor() && mRecoData->isTrackSourceLoaded(tid.getSource())) {
510 tracks.emplace_back().setPID(mRecoData->getTrackParam(tid).getPID());
511 gidsITS.push_back(mRecoData->getITSContributorGID(tid));
512 }
513 }
514 ntr = tracks.size();
515 useTrack.resize(ntr);
516#ifdef WITH_OPENMP
517#pragma omp parallel for schedule(dynamic) num_threads(mNThreads)
518#endif
519 for (int itr = 0; itr < ntr; itr++) {
520 if (!(useTrack[itr] = refitITStrack(tracks[itr], gidsITS[itr]))) {
521 tracks[itr] = mRecoData->getTrackParam(gidsITS[itr]); // this track will not be used but participates in prepareVertexRefit
522 }
523 }
524 ntr = 0;
525 for (auto v : useTrack) {
526 ntr++;
527 }
528 if (ntr < params.minPVContributors || !mVertexer.prepareVertexRefit(tracks, pv)) {
529 LOGP(warn, "Abandon vertex refit: NcontribNew = {} vs NcontribOld = {}", ntr, ntrIni);
530 return false;
531 }
532 LOGP(debug, "Original vtx: Nc:{} {}, chi2={}", pv.getNContributors(), pv.asString(), pv.getChi2());
533 auto pvSave = pv;
534 pv = mVertexer.refitVertexFull(useTrack, pv);
535 LOGP(debug, "Refitted vtx: Nc:{} {}, chi2={}", ntr, pv.asString(), pv.getChi2());
536 if (pv.getChi2() < 0.f) {
537 LOGP(warn, "Failed to refit PV {}", pvSave.asString());
538 return false;
539 }
540 return true;
541}
542
543bool CheckResidSpec::refitITStrack(o2::track::TrackParCov& track, GTrackID gid)
544{
545 // destination tack might have non-default PID assigned
546 const auto& trkITS = mRecoData->getITSTrack(gid);
547 const auto itsClRefs = mRecoData->getITSTracksClusterRefs();
548 const auto& params = CheckResidConfig::Instance();
549 auto pid = track.getPID();
550 track = trkITS.getParamOut();
551 track.resetCovariance();
552 track.setCov(track.getQ2Pt() * track.getQ2Pt() * track.getCov()[14], 14);
553 track.setPID(pid);
554 auto nCl = trkITS.getNumberOfClusters();
556 auto prop = o2::base::Propagator::Instance();
557 float bz = prop->getNominalBz();
558 o2::track::TrackPar refLin{track};
559
560 for (int iCl = 0; iCl < nCl; iCl++) { // clusters are stored from outer to inner layers
561 const auto& cls = mITSClustersArray[itsClRefs[trkITS.getClusterEntry(iCl)]];
562 auto alpha = geom->getSensorRefAlpha(cls.getSensorID());
563 if (!(params.useStableRef ? track.rotate(alpha, refLin, bz) : track.rotate(alpha)) ||
564 !prop->propagateTo(track, params.useStableRef ? &refLin : nullptr, cls.getX(), true)) {
565 LOGP(debug, "refitITStrack failed on propagation to cl#{}, alpha={}, x={} | {}", iCl, alpha, cls.getX(), track.asString());
566 return false;
567 }
568 if (!track.update(cls)) {
569 LOGP(debug, "refitITStrack failed on update with cl#{}, | {}", iCl, track.asString());
570 return false;
571 }
572 }
573 return true;
574}
575
576void CheckResidSpec::fillHistos(const o2::checkresid::Track& trc)
577{
578 const auto& params = CheckResidConfig::Instance();
579 int np = trc.points.size();
580 auto pt = trc.track.getPt();
581 if (pt < params.minPt || pt > params.maxPt) {
582 return;
583 }
584 for (int ip = 0; ip < np; ip++) {
585 const auto& pnt = trc.points[ip];
586 int il = pnt.lr >= 0 ? pnt.lr + 1 : 0;
587 mHMan->getHisto2F(il * 10 + 0 * 100)->Fill(pnt.phi, pnt.dy);
588 mHMan->getHisto2F(il * 10 + 0 * 100 + 1000)->Fill(pnt.z, pnt.dy);
589 mHMan->getHisto2F(il * 10 + 0 * 100 + 2000)->Fill(pt, pnt.dy);
590 mHMan->getHisto2F(il * 10 + 0 * 100 + 3000)->Fill(trc.track.getTgl(), pnt.dy);
591 if (pnt.sig2y > 0) {
592 auto pull = pnt.dy / std::sqrt(pnt.sig2y);
593 mHMan->getHisto2F(il * 10 + 0 * 100 + 5)->Fill(pnt.phi, pull);
594 mHMan->getHisto2F(il * 10 + 0 * 100 + 5 + 1000)->Fill(pnt.z, pull);
595 mHMan->getHisto2F(il * 10 + 0 * 100 + 5 + 2000)->Fill(pt, pull);
596 mHMan->getHisto2F(il * 10 + 0 * 100 + 5 + 3000)->Fill(trc.track.getTgl(), pull);
597 }
598 mHMan->getHisto2F(il * 10 + 1 * 100)->Fill(pnt.phi, pnt.dz);
599 mHMan->getHisto2F(il * 10 + 1 * 100 + 1000)->Fill(pnt.z, pnt.dz);
600 mHMan->getHisto2F(il * 10 + 1 * 100 + 2000)->Fill(pt, pnt.dz);
601 mHMan->getHisto2F(il * 10 + 1 * 100 + 3000)->Fill(trc.track.getTgl(), pnt.dz);
602 if (pnt.sig2z > 0) {
603 auto pull = pnt.dz / std::sqrt(pnt.sig2z);
604 mHMan->getHisto2F(il * 10 + 1 * 100 + 5)->Fill(pnt.phi, pull);
605 mHMan->getHisto2F(il * 10 + 1 * 100 + 5 + 1000)->Fill(pnt.z, pull);
606 mHMan->getHisto2F(il * 10 + 1 * 100 + 5 + 2000)->Fill(pt, pull);
607 mHMan->getHisto2F(il * 10 + 1 * 100 + 5 + 3000)->Fill(trc.track.getTgl(), pull);
608 }
609 }
610 //--------------
611 if (trc.trIBOut.getX() > 1 && std::abs(trc.trIBOut.getX() - trc.trOBInw.getX()) < 0.1) {
612 for (int ip = 0; ip < 5; ip++) {
613 float d = trc.trIBOut.getParam(ip) - trc.trOBInw.getParam(ip);
614 mHMan->getHisto2F(10000 + ip * 10)->Fill(trc.trIBOut.getPhiPos(), d);
615 mHMan->getHisto2F(11000 + ip * 10)->Fill(trc.trIBOut.getZ(), d);
616 mHMan->getHisto2F(12000 + ip * 10)->Fill(pt, d);
617 mHMan->getHisto2F(13000 + ip * 10)->Fill(trc.track.getTgl(), d);
618 float sg = trc.trIBOut.getCovarElem(ip, ip) + trc.trOBInw.getCovarElem(ip, ip);
619 if (sg > 0) {
620 auto pull = d / std::sqrt(sg);
621 mHMan->getHisto2F(10000 + ip * 10 + 5)->Fill(trc.trIBOut.getPhiPos(), pull);
622 mHMan->getHisto2F(11000 + ip * 10 + 5)->Fill(trc.trIBOut.getZ(), pull);
623 mHMan->getHisto2F(12000 + ip * 10 + 5)->Fill(pt, pull);
624 mHMan->getHisto2F(13000 + ip * 10 + 5)->Fill(trc.track.getTgl(), pull);
625 }
626 }
627 }
628}
629
630void CheckResidSpec::bookHistos()
631{
633 mHManV.emplace_back() = std::make_unique<o2::HistoManager>("", fmt::format("{}_hman.root", params.outname));
634 mHMan = mHManV.back().get();
635 mHMan->SetName(params.outname.c_str());
636 auto defLogAxis = [](float xMn, float xMx, int nbin) { // get array for log axis
637 if (xMn <= 0 || xMx <= xMn || nbin < 2) {
638 LOGP(fatal, "Wrong log axis request: xmin = {} xmax = {} nbins = {}", xMn, xMx, nbin);
639 }
640 auto dx = std::log(xMx / xMn) / nbin;
641 std::vector<double> xax(nbin + 1);
642 for (int i = 0; i <= nbin; i++) {
643 xax[i] = xMn * std::exp(dx * i);
644 }
645 return xax;
646 };
647 float minPt = std::max(0.1f, params.minPt), maxPt = std::min(50.f, params.maxPt);
648 auto ptax = defLogAxis(minPt, maxPt, params.nBinsPt);
649
650 for (int il = 0; il < 8; il++) {
651 std::string lrName = il == 0 ? "Vtx" : fmt::format("Lr{}", il - 1);
652 for (int iyz = 0; iyz < 2; iyz++) {
653 std::string dname = iyz == 0 ? "dy" : "dz", dtit = iyz == 0 ? "#DeltaY" : "#DeltaZ";
654 auto h2 = new TH2F(fmt::format("{}_{}_{}", dname, lrName, "phi").c_str(), fmt::format("{}_{{{}}} vs {};#phi;{}", dtit, lrName, "#phi", dtit).c_str(), params.nBinsPhi, 0, TMath::Pi() * 2, params.nBinsRes, -params.maxDYZ[il], params.maxDYZ[il]);
655 mHMan->addHisto(h2, il * 10 + iyz * 100);
656 auto h2p = new TH2F(fmt::format("{}_{}_{}_pull", dname, lrName, "phi").c_str(), fmt::format("pull {}_{{{}}} vs {};#phi; pull{}", dtit, lrName, "phi", dtit).c_str(), params.nBinsPhi, 0, TMath::Pi() * 2, params.nBinsRes, -params.maxPull, params.maxPull);
657 mHMan->addHisto(h2p, il * 10 + iyz * 100 + 5);
658
659 auto hz2 = new TH2F(fmt::format("{}_{}_{}", dname, lrName, "Z").c_str(), fmt::format("{}_{{{}}} vs {};Z;{}", dtit, lrName, "Z", dtit).c_str(), params.nBinsZ, -params.zranges[il], params.zranges[il], params.nBinsRes, -params.maxDYZ[il], params.maxDYZ[il]);
660 mHMan->addHisto(hz2, il * 10 + iyz * 100 + 1000);
661 auto hz2p = new TH2F(fmt::format("{}_{}_{}_pull", dname, lrName, "Z").c_str(), fmt::format("pull {}_{{{}}} vs {};Z; pull{}", dtit, lrName, "Z", dtit).c_str(), params.nBinsZ, -params.zranges[il], params.zranges[il], params.nBinsRes, -params.maxPull, params.maxPull);
662 mHMan->addHisto(hz2p, il * 10 + iyz * 100 + 5 + 1000);
663
664 auto hpt2 = new TH2F(fmt::format("{}_{}_{}", dname, lrName, "Pt").c_str(), fmt::format("{}_{{{}}} vs {};p_{{T}};{}", dtit, lrName, "p_{T}", dtit).c_str(), params.nBinsPt, ptax.data(), params.nBinsRes, -params.maxDYZ[il], params.maxDYZ[il]);
665 mHMan->addHisto(hpt2, il * 10 + iyz * 100 + 2000);
666 auto hpt2p = new TH2F(fmt::format("{}_{}_{}_pull", dname, lrName, "Pt").c_str(), fmt::format("pull {}_{{{}}} vs {};p_{{T}}; pull{}", dtit, lrName, "p_{T}", dtit).c_str(), params.nBinsPt, ptax.data(), params.nBinsRes, -params.maxPull, params.maxPull);
667 mHMan->addHisto(hpt2p, il * 10 + iyz * 100 + 5 + 2000);
668
669 auto htgl2 = new TH2F(fmt::format("{}_{}_{}", dname, lrName, "tgl").c_str(), fmt::format("{}_{{{}}} vs {};tg#lambda;{}", dtit, lrName, "tg#lambda", dtit).c_str(), params.nBinsTgl, -params.maxTgl, params.maxTgl, params.nBinsRes, -params.maxDYZ[il], params.maxDYZ[il]);
670 mHMan->addHisto(htgl2, il * 10 + iyz * 100 + 3000);
671 auto htgl2p = new TH2F(fmt::format("{}_{}_{}_pull", dname, lrName, "tgl").c_str(), fmt::format("pull {}_{{{}}} vs {};tg#lambda; pull{}", dtit, lrName, "tg#lambda", dtit).c_str(), params.nBinsTgl, -params.maxTgl, params.maxTgl, params.nBinsRes, -params.maxPull, params.maxPull);
672 mHMan->addHisto(htgl2p, il * 10 + iyz * 100 + 5 + 3000);
673 }
674 }
675
676 for (int ip = 0; ip < 5; ip++) {
677 auto h2 = new TH2F(fmt::format("dPar{}_IBOBphi", ip).c_str(), fmt::format("#Delta par{} IB-OB vs phi;#phi;#Delta par{}", ip, ip).c_str(), params.nBinsPhi, 0, TMath::Pi() * 2, params.nBinsRes, -params.maxDPar[ip], params.maxDPar[ip]);
678 mHMan->addHisto(h2, 10000 + ip * 10);
679 auto h2p = new TH2F(fmt::format("dPar{}_IBOBphi_pull", ip).c_str(), fmt::format("pull #Delta par{} IB-OB vs phi;#phi;pull #Delta par{}", ip, ip).c_str(), params.nBinsPhi, 0, TMath::Pi() * 2, params.nBinsRes, -params.maxPull, params.maxPull);
680 mHMan->addHisto(h2p, 10000 + ip * 10 + 5);
681
682 auto hz2 = new TH2F(fmt::format("dPar{}_IBOBz", ip).c_str(), fmt::format("#Delta par{} IB-OB vs Z;Z;#Delta par{}", ip, ip).c_str(), params.nBinsZ, -20., 20., params.nBinsRes, -params.maxDPar[ip], params.maxDPar[ip]);
683 mHMan->addHisto(hz2, 11000 + ip * 10);
684 auto hz2p = new TH2F(fmt::format("dPar{}_IBOBz_pull", ip).c_str(), fmt::format("pull #Delta par{} IB-OB vs Z;Z;pull #Delta par{}", ip, ip).c_str(), params.nBinsZ, -20., 20., params.nBinsRes, -params.maxPull, params.maxPull);
685 mHMan->addHisto(hz2p, 11000 + ip * 10 + 5);
686
687 auto hpt2 = new TH2F(fmt::format("dPar{}_IBOBpt", ip).c_str(), fmt::format("#Delta par{} IB-OB vs pT;p_{{T}};#Delta par{}", ip, ip).c_str(), params.nBinsPt, ptax.data(), params.nBinsRes, -params.maxDPar[ip], params.maxDPar[ip]);
688 mHMan->addHisto(hpt2, 12000 + ip * 10);
689 auto hpt2p = new TH2F(fmt::format("dPar{}_IBOBpt_pull", ip).c_str(), fmt::format("pull #Delta par{} IB-OB vs pT;p_{{T}};pull #Delta par{}", ip, ip).c_str(), params.nBinsPt, ptax.data(), params.nBinsRes, -params.maxPull, params.maxPull);
690 mHMan->addHisto(hpt2p, 12000 + ip * 10 + 5);
691
692 auto htgl2 = new TH2F(fmt::format("dPar{}_IBOBtgl", ip).c_str(), fmt::format("#Delta par{} IB-OB vs tg#lambda;tg#lambda;#Delta par{}", ip, ip).c_str(), params.nBinsTgl, -params.maxTgl, params.maxTgl, params.nBinsRes, -params.maxDPar[ip], params.maxDPar[ip]);
693 mHMan->addHisto(htgl2, 13000 + ip * 10);
694 auto htgl2p = new TH2F(fmt::format("dPar{}_IBOBtgl_pull", ip).c_str(), fmt::format("pull #Delta par{} IB-OB vs tg#lambda;tg#lambda;pull #Delta par{}", ip, ip).c_str(), params.nBinsTgl, -params.maxTgl, params.maxTgl, params.nBinsRes, -params.maxPull, params.maxPull);
695 mHMan->addHisto(htgl2p, 13000 + ip * 10 + 5);
696 }
697}
698
699void CheckResidSpec::postProcessHistos()
700{
701 printf("Fitting histos\n");
702 if (!mHMan) {
703 if (mHManV.empty()) {
704 LOGP(warn, "nothing to process");
705 return;
706 }
707 mHMan = mHManV[0].get();
708 }
710 auto gs = new TF1("gs", "gaus", -1, 1);
711 int maxH = mPostProcOnly ? mHManV.size() : 1;
712 TObjArray arr;
713 for (int ihm = 0; ihm < maxH; ihm++) {
714 auto* histm = mHManV[ihm].get();
715 auto fitSlices = [&](int id) {
716 auto h2 = histm->getHisto2F(id);
717 if (!h2 || h2->GetEntries() < params.minHistoStat2Fit) {
718 return;
719 }
720 h2->FitSlicesY(gs, 0, -1, 0, "QNR", &arr);
721 arr.SetOwner(true);
722 TH1* hmean = (TH1*)arr.RemoveAt(1);
723 if (hmean) {
724 hmean->SetTitle(Form("<%s>", h2->GetTitle()));
725 histm->addHisto(hmean, id + 1);
726 }
727 TH1* hsig = (TH1*)arr.RemoveAt(2);
728 if (hsig) {
729 hsig->SetTitle(Form("#sigma(%s)", h2->GetTitle()));
730 histm->addHisto(hsig, id + 2);
731 }
732 };
733 for (int ioffs = 0; ioffs <= 3; ioffs++) { // vs phi, Z, pT, tgl
734 int offs = ioffs * 1000;
735 for (int iht = 0; iht < 2; iht++) { // resid, pull
736 int offsV = iht == 0 ? 0 : 5;
737 for (int il = 0; il < 8; il++) {
738 for (int iyz = 0; iyz < 2; iyz++) {
739 fitSlices(il * 10 + iyz * 100 + offsV + offs);
740 }
741 }
742 for (int ip = 0; ip < 5; ip++) {
743 fitSlices(10000 + ip * 10 + offsV + offs);
744 }
745 }
746 }
747 histm->write();
748 }
749 delete gs;
750}
751
752void CheckResidSpec::drawHistos()
753{
754 gROOT->SetBatch(true);
755 gStyle->SetTitleX(0.2);
756 gStyle->SetTitleY(0.88);
757 gStyle->SetTitleW(0.25);
758 gStyle->SetOptStat(0);
759 int nhm = mHManV.size();
760 std::array<unsigned int, 3> hcol{EColor::kRed, EColor::kBlue, EColor::kGreen + 2};
761 std::unique_ptr<TLegend> lg;
762 lg = std::make_unique<TLegend>(0.12, 0.13, 0.9, 0.13 + std::min(0.5f, nhm * 0.2f / 3.f));
763 lg->SetFillStyle(0);
764 lg->SetBorderSize(0);
765 for (int i = 0; i < nhm; i++) {
766 auto hman = mHManV[i].get();
767 if (!hman || hman->GetLast() < 1) {
768 continue;
769 }
770 hman->setMarkerStyle(20 + i + (i % 2) * 4, 0.5);
771 hman->setColor(hcol[i % hcol.size()]);
772 auto le = lg->AddEntry(hman->getHisto(1), hman->GetName(), "lp");
773 le->SetTextColor(hcol[i % hcol.size()]);
774 }
775 TCanvas cly("cly", "", 600, 800), clz("clz", "", 600, 800), clpar("clpar", "", 600, 800);
776 TCanvas czly("czly", "", 600, 800), czlz("czlz", "", 600, 800), czlpar("czlpar", "", 600, 800);
778
779 auto AddLabel = [](const char* txt, float x = 0.1, float y = 0.9, int color = kBlack, float size = 0.04) {
780 TLatex* lt = new TLatex(x, y, txt);
781 lt->SetNDC();
782 lt->SetTextColor(color);
783 lt->SetTextSize(size);
784 lt->Draw();
785 return lt;
786 };
787
788 auto drawResLr = [this](TCanvas& canv, int offs, const float resMM[8], bool logX) {
789 canv.Clear();
790 canv.Divide(2, 4);
791 int nh = this->mHManV.size();
792 for (int i = 0; i < 8; i++) {
793 canv.cd(i + 1);
794 bool same = false;
795 for (int j = 0; j < nh; j++) {
796 auto hman = this->mHManV[j].get();
797 if (!hman || hman->GetLast() < 1) {
798 continue;
799 }
800 if (auto histo = hman->getHisto(10 * i + offs)) {
801 histo->Draw(same ? "same" : "");
802 if (!same) {
803 histo->SetMinimum(-resMM[i]);
804 histo->SetMaximum(resMM[i]);
805 same = true;
806 }
807 }
808 }
809 gPad->SetGrid();
810 gPad->SetLogx(logX);
811 }
812 };
813
814 auto drawResPar = [this](TCanvas& canv, int offs, const float resMM[8], bool logX) {
815 canv.Clear();
816 canv.Divide(2, 3);
817 int nh = this->mHManV.size();
818 for (int i = 0; i < 5; i++) {
819 canv.cd(i + 1);
820 bool same = false;
821 for (int j = 0; j < nh; j++) {
822 auto hman = this->mHManV[j].get();
823 if (!hman || hman->GetLast() < 1) {
824 continue;
825 }
826 if (auto histo = hman->getHisto(10 * i + offs)) {
827 histo->Draw(same ? "same" : "");
828 if (!same) {
829 histo->SetMinimum(-resMM[i]);
830 histo->SetMaximum(resMM[i]);
831 same = true;
832 }
833 }
834 }
835 gPad->SetGrid();
836 gPad->SetLogx(logX);
837 }
838 };
839
840 cly.Print(Form("%s_hman.pdf[", params.outname.c_str()));
841 drawResLr(cly, 1, params.resMMLrY, false);
842 cly.cd(2);
843 lg->Draw();
844 AddLabel("Y residuals", 0.1, 0.95);
845 cly.Print(Form("%s_hman.pdf", params.outname.c_str()));
846
847 drawResLr(clz, 101, params.resMMLrZ, false);
848 clz.cd(2);
849 lg->Draw();
850 AddLabel("Z residuals", 0.1, 0.95);
851 clz.Print(Form("%s_hman.pdf", params.outname.c_str()));
852
853 drawResLr(czly, 1001, params.resMMLrY, false);
854 czly.cd(2);
855 lg->Draw();
856 AddLabel("Y residuals", 0.1, 0.95);
857 czly.Print(Form("%s_hman.pdf", params.outname.c_str()));
858
859 drawResLr(czlz, 1101, params.resMMLrZ, false);
860 czlz.cd(2);
861 lg->Draw();
862 AddLabel("Z residuals", 0.1, 0.95);
863 czlz.Print(Form("%s_hman.pdf", params.outname.c_str()));
864
865 drawResLr(czly, 2001, params.resMMLrY, true);
866 czly.cd(2);
867 lg->Draw();
868 AddLabel("Y residuals", 0.1, 0.95);
869 czly.Print(Form("%s_hman.pdf", params.outname.c_str()));
870
871 drawResLr(czlz, 2101, params.resMMLrZ, true);
872 czlz.cd(2);
873 lg->Draw();
874 AddLabel("Z residuals", 0.1, 0.95);
875 czlz.Print(Form("%s_hman.pdf", params.outname.c_str()));
876
877 drawResLr(czly, 3001, params.resMMLrY, false);
878 czly.cd(2);
879 lg->Draw();
880 AddLabel("Y residuals", 0.1, 0.95);
881 czly.Print(Form("%s_hman.pdf", params.outname.c_str()));
882
883 drawResLr(czlz, 3101, params.resMMLrZ, false);
884 czlz.cd(2);
885 lg->Draw();
886 AddLabel("Z residuals", 0.1, 0.95);
887 czlz.Print(Form("%s_hman.pdf", params.outname.c_str()));
888
889 drawResPar(clpar, 10001, params.resMMPar, false);
890 clpar.cd(6);
891 lg->Draw();
892 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
893 clpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
894
895 drawResPar(czlpar, 11001, params.resMMPar, false);
896 czlpar.cd(6);
897 lg->Draw();
898 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
899 czlpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
900
901 drawResPar(czlpar, 12001, params.resMMPar, true);
902 czlpar.cd(6);
903 lg->Draw();
904 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
905 czlpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
906
907 drawResPar(czlpar, 13001, params.resMMPar, false);
908 czlpar.cd(6);
909 lg->Draw();
910 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
911 czlpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
912
913 cly.Print(Form("%s_hman.pdf]", params.outname.c_str()));
914}
915
917{
918 mDBGOut.reset();
919 if (mHManV.size()) {
920 postProcessHistos();
921 }
922 if (mDraw) {
923 drawHistos();
924 }
925}
926
928{
930 return;
931 }
932 if (matcher == ConcreteDataMatcher("GLO", "MEANVERTEX", 0)) {
933 LOG(info) << "Imposing new MeanVertex: " << ((const o2::dataformats::MeanVertexObject*)obj)->asString();
934 mMeanVtx = *(const o2::dataformats::MeanVertexObject*)obj;
935 mMeanVertexUpdated = true;
936 return;
937 }
938 if (matcher == ConcreteDataMatcher("ITS", "CLUSDICT", 0)) {
939 LOG(info) << "cluster dictionary updated";
940 mITSDict = (const o2::itsmft::TopologyDictionary*)obj;
941 return;
942 }
943}
944
945DataProcessorSpec getCheckResidSpec(GTrackID::mask_t srcTracks, GTrackID::mask_t srcClusters, bool drawOnly, bool postProcOnly)
946{
947 std::vector<OutputSpec> outputs;
948 auto dataRequest = std::make_shared<DataRequest>();
949 if (!drawOnly && !postProcOnly) {
950 bool useMC = false;
951 dataRequest->requestTracks(srcTracks, useMC);
952 dataRequest->requestClusters(srcClusters, useMC);
953 dataRequest->requestPrimaryVertices(useMC);
954 dataRequest->inputs.emplace_back("meanvtx", "GLO", "MEANVERTEX", 0, Lifetime::Condition, ccdbParamSpec("GLO/Calib/MeanVertex", {}, 1));
955 }
956 auto ggRequest = drawOnly ? std::make_shared<o2::base::GRPGeomRequest>(false, false, false, false, false, o2::base::GRPGeomRequest::None, dataRequest->inputs) : std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
957 true, // GRPECS=true
958 true, // GRPLHCIF
959 true, // GRPMagField
960 true, // askMatLUT
962 dataRequest->inputs, true);
963 Options opts;
964 if (!drawOnly) {
965 opts = Options{
966 {"nthreads", VariantType::Int, 1, {"number of threads"}},
967 {"no-tree", VariantType::Bool, false, {"do not fill residuals tree"}},
968 {"no-hist", VariantType::Bool, false, {"do not fill residuals histograms"}},
969 {"draw-report", VariantType::Bool, false, {"fill residuals report"}},
970 };
971 }
972
973 return DataProcessorSpec{
974 "check-resid",
975 dataRequest->inputs,
976 outputs,
977 AlgorithmSpec{adaptFromTask<CheckResidSpec>(dataRequest, ggRequest, srcTracks, drawOnly, postProcOnly)},
978 opts};
979}
980
981} // namespace o2::checkresid
Wrapper container for different reconstructed object types.
Base track model for the Barrel, params only, w/o covariance.
Definition of the GeometryManager class.
std::ostringstream debug
int32_t i
Helper for geometry and GRP related CCDB requests.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
Definition of the GeometryTGeo class.
Utility functions for MC particles.
Definition of the Names Generator class.
Primary vertex finder.
uint32_t j
Definition RawData.h:0
uint16_t pid
Definition RawData.h:2
uint32_t res
Definition RawData.h:0
Wrapper container for different reconstructed object types.
o2::track::TrackParCov TrackParCov
Definition Recon.h:39
Result of refitting TPC-ITS matched track.
Reference on ITS/MFT clusters set.
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
void setSigmaZ2(T v)
void setSigmaY2(T v)
void setSensorID(std::int16_t sid)
Definition BaseCluster.h:92
void setXYZ(T x, T y, T z)
int addHisto(TH1 *histo, int at=-1)
TH2F * getHisto2F(int id) const
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
void run(ProcessingContext &pc) final
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
CheckResidSpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, GTrackID::mask_t src, bool drawOnly, bool postProcOnly)
void init(InitContext &ic) final
static void updateFromString(std::string const &)
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
static constexpr ID ITS
Definition DetID.h:63
static constexpr ID TPC
Definition DetID.h:64
ServiceRegistryRef services()
Definition InitContext.h:34
ConfigParamRegistry const & options()
Definition InitContext.h:33
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
ServiceRegistryRef services()
The services registry associated with this processing context.
static GeometryTGeo * Instance()
int getClusterEntry(int i) const
Definition TrackITS.h:67
PVertex refitVertexFull(const std::vector< bool > useTrack, const o2d::VertexBase &vtxSeed)
bool prepareVertexRefit(const TR &tracks, const o2d::VertexBase &vtxSeed)
Definition PVertexer.h:362
void setMeanVertex(const o2d::MeanVertexObject *v)
Definition PVertexer.h:98
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum src
Definition glcorearb.h:1767
GLsizeiptr size
Definition glcorearb.h:659
GLuint color
Definition glcorearb.h:1272
const GLdouble * v
Definition glcorearb.h:832
GLdouble f
Definition glcorearb.h:310
GLenum const GLfloat * params
Definition glcorearb.h:272
GLuint id
Definition glcorearb.h:650
o2::framework::DataProcessorSpec getCheckResidSpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool drawOnly, bool postProcOnly)
create a processor spec
constexpr double LHCBunchSpacingNS
Defining ITS Vertex explicitly as messageable.
Definition Cartesian.h:288
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
std::vector< ConfigParamSpec > Options
void convertCompactClusters(gsl::span< const itsmft::CompClusterExt > clusters, gsl::span< const unsigned char >::iterator &pattIt, std::vector< o2::BaseCluster< float > > &output, const itsmft::TopologyDictionary *dict)
convert compact clusters to 3D spacepoints
Definition IOUtils.cxx:35
int int int float float float int nCl
const TrackingFrameInfo *const const Cluster *const const float const float bz
double * getX(double *xyDxy, int N)
double * getY(double *xyDxy, int N)
TrackParCovF TrackParCov
Definition Track.h:33
struct o2::upgrades_utils::@469 tracks
structure to keep trigger-related info
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
o2::track::TrackPar track
o2::track::TrackParCov trIBOut
o2::track::TrackParCov trOBInw
std::vector< Point > points
GTrackID getITSContributorGID(GTrackID source) const
const o2::tpc::TrackTPC & getTPCTrack(GTrackID id) const
GTrackID getTPCContributorGID(GTrackID source) const
const o2::track::TrackParCov & getTrackParam(GTrackID gidx) const
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
static constexpr int L2G
Definition Cartesian.h:54
static constexpr int T2L
Definition Cartesian.h:55
static constexpr int T2G
Definition Cartesian.h:56
static std::vector< std::string > tokenize(const std::string &src, char delim, bool trimToken=true, bool skipEmpty=true)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"