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, o2::track::PID pid);
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 std::vector<o2::dataformats::PrimaryVertex> mPVUsed;
119 o2::HistoManager* mHMan = nullptr;
120};
121
123{
124 mDraw = true;
125 if (!mDrawOnly) {
126 mDraw = ic.options().get<bool>("draw-report");
127 mFillHistos = !ic.options().get<bool>("no-hist");
128 mFillTree = !ic.options().get<bool>("no-tree");
129 mNThreads = ic.options().get<int>("nthreads");
130 }
132 int lane = ic.services().get<const o2::framework::DeviceSpec>().inputTimesliceId;
133 int maxLanes = ic.services().get<const o2::framework::DeviceSpec>().maxInputTimeslices;
134 std::string nm = params.outname;
135 if (maxLanes > 1) {
136 o2::conf::ConfigurableParam::updateFromString(fmt::format("checkresid.outname={}_t{}", nm, lane));
137 }
138 if (mDraw) {
139 mFillHistos = true;
140 }
141 if (!mDrawOnly && mFillHistos) {
142 bookHistos();
143 }
144 if (!params.ext_hm_list.empty()) {
145 auto vecNames = o2::utils::Str::tokenize(params.ext_hm_list, ',');
146 auto vecLegends = o2::utils::Str::tokenize(params.ext_leg_list, ',');
147 bool useLeg = true;
148 if (vecNames.size() != vecLegends.size()) {
149 LOGP(warn, "{} legend names provided for {} external histomanagers, will use file names as legends", vecLegends.size(), vecNames.size());
150 useLeg = false;
151 }
152 int cntH = 0;
153 for (const auto& vn : vecNames) {
154 LOGP(info, "Loading external HistoManager {}", vn);
155 mHManV.emplace_back() = std::make_unique<o2::HistoManager>("", vn, true);
156 auto hm = mHManV.back().get();
157 if (!hm) {
158 LOGP(error, "Failed to load histograms from {}", vn);
159 mHManV.pop_back();
160 } else {
161 hm->SetName(useLeg ? vecLegends[cntH].c_str() : vn.c_str());
162 }
163 cntH++;
164 }
165 }
166 if (mDrawOnly) {
167 return;
168 }
170#ifndef WITH_OPENMP
171 if (mNThreads > 1) {
172 LOGP(warn, "No OpenMP");
173 }
174 mNThreads = 1;
175#endif
176 if (mFillTree) {
177 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>(fmt::format("{}.root", params.outname).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 mPVUsed.clear();
260 slots.resize(mNThreads);
261 int nvGood = 0, nvUse = 0, nvRefFail = 0;
262 long pvFitDuration{};
263 for (int iv = 0; iv < nv; iv++) {
264 const auto& vtref = vtxRefs[iv];
265 auto pve = pvvec[iv];
266 if (pve.getNContributors() < params.minPVContributors) {
267 continue;
268 }
269 nvGood++;
270 if (params.refitPV) {
271 LOGP(debug, "Refitting PV#{} of {} tracks", iv, pve.getNContributors());
272 auto tStartPVF = std::chrono::time_point_cast<std::chrono::microseconds>(std::chrono::system_clock::now()).time_since_epoch().count();
273 bool res = refitPV(pve, iv);
274 pvFitDuration += std::chrono::time_point_cast<std::chrono::microseconds>(std::chrono::system_clock::now()).time_since_epoch().count() - tStartPVF;
275 if (!res) {
276 nvRefFail++;
277 continue;
278 }
279 }
280 nvUse++;
281 mPVUsed.push_back(pve);
282 const auto& itsContRefs = vtref.getITSGloContributors();
283 int idMinITSGlo = 0, idMaxITSGlo = 0;
284 if (params.useITSGloContributors) {
285 if (itsContRefs.getFirstEntry() < vtref.getEntries() && itsContRefs.getEntries() == 0) {
286 LOGP(fatal, "Usage of stored ITS global contributors is requested but they are missing");
287 }
288 idMinITSGlo = itsContRefs.getFirstEntry();
289 idMaxITSGlo = idMinITSGlo + itsContRefs.getEntries();
290 }
291 int cntPVCont = 0;
292 for (int is = 0; is < GTrackID::NSources; is++) {
293 if (!params.useITSGloContributors && (!mTracksSrc[is] || !mRecoData->isTrackSourceLoaded(is))) {
294 continue;
295 }
296 int idMin = vtref.getFirstEntryOfSource(is), idMax = idMin + vtref.getEntriesOfSource(is);
297 DetID::mask_t dm = GTrackID::getSourceDetectorsMask(is);
298 if (!dm[DetID::ITS]) {
299 continue;
300 }
301#ifdef WITH_OPENMP
302#pragma omp parallel for schedule(dynamic) num_threads(mNThreads)
303#endif
304 for (int i = idMin; i < idMax; i++) {
305 auto vid = trackIndex[i];
306 bool pvCont = vid.isPVContributor();
307 if (!pvCont && (params.pvcontribOnly || params.useITSGloContributors)) {
308 continue;
309 }
310 GTrackID gidITS;
311 o2::track::PID pidITS;
312 if (params.useITSGloContributors) {
313 int id = idMinITSGlo + cntPVCont;
314 if (id >= idMaxITSGlo) {
315 LOGP(fatal, "Calculated GlobalContributor ITS track index {} exceeds number of stored indices {}", id, itsContRefs.getEntries());
316 }
317 pidITS = trackIndex[id].getSource();
318 gidITS = GTrackID(trackIndex[id].getIndex(), GTrackID::ITS);
319 cntPVCont++;
320 } else {
321 gidITS = mRecoData->getITSContributorGID(vid);
322 if (gidITS.getSource() != GTrackID::ITS) {
323 continue;
324 }
325 }
326 const auto& itsTrack = mRecoData->getITSTrack(gidITS);
327 if (itsTrack.getNClusters() < params.minITSCl) {
328 continue;
329 }
330 const auto& trc = params.useITSGloContributors ? ((o2::track::TrackParCov&)itsTrack) : mRecoData->getTrackParam(vid);
331 if (!params.useITSGloContributors) {
332 pidITS = trc.getPID();
333 }
334 auto pt = trc.getPt();
335 if (pt < params.minPt || pt > params.maxPt) {
336 continue;
337 }
338 if (std::abs(trc.getTgl()) > params.maxTgl) {
339 continue;
340 }
341#ifdef WITH_OPENMP
342 auto& accum = slots[omp_get_thread_num()];
343#else
344 auto& accum = slots[0];
345#endif
346 auto& resTrack = accum.emplace_back();
347 resTrack.gid = vid;
348 if (!processITSTrack(itsTrack, pve, resTrack, pidITS)) {
349 accum.pop_back();
350 continue;
351 }
352 }
353 }
354 }
355 // output
356 for (const auto& accum : slots) {
357 for (const auto& tr : accum) {
358 if (mDBGOut) {
359 (*mDBGOut) << "res" << "tr=" << tr << "\n";
360 }
361 if (mHMan) {
362 fillHistos(tr);
363 }
364 }
365 }
366 if (mDBGOut) {
367 (*mDBGOut) << "pvUsed" << "pv=" << mPVUsed << "\n";
368 }
369 if (mHMan) {
370 for (const auto& pv : mPVUsed) {
371 mHMan->getHisto(20000 + 0)->Fill(pv.getX());
372 mHMan->getHisto(20000 + 1)->Fill(pv.getY());
373 mHMan->getHisto(20000 + 2)->Fill(pv.getZ());
374 mHMan->getHisto(20000 + 3)->Fill(pv.getNContributors());
375 }
376 }
377 LOGP(info, "processed {} PVs out of {} good vertices (out of {} in total), PV refits took {} mus, {} refits failed", nvUse, nvGood, nv, pvFitDuration, nvRefFail);
378 TFCount++;
379}
380
381bool CheckResidSpec::processITSTrack(const o2::its::TrackITS& iTrack, const o2::dataformats::PrimaryVertex& pv, o2::checkresid::Track& resTrack, o2::track::PID pid)
382{
383 const auto itsClRefs = mRecoData->getITSTracksClusterRefs();
384 auto trFitInw = iTrack.getParamOut(); // seed for inward refit
385 auto trFitOut = iTrack.getParamIn(); // seed for outward refit
386 trFitInw.setPID(pid);
387 trFitOut.setPID(pid);
388 auto prop = o2::base::Propagator::Instance();
390 float pvAlpha = 0;
391 float bz = prop->getNominalBz();
392 std::array<const o2::BaseCluster<float>*, 8> clArr{};
393 const auto& params = CheckResidConfig::Instance();
394 std::array<o2::track::TrackParCov, 8> extrapOut, extrapInw; // 2-way Kalman extrapolations, vertex + 7 layers
395
396 auto rotateTrack = [bz](o2::track::TrackParCov& tr, float alpha, o2::track::TrackPar* refLin) {
397 return refLin ? tr.rotate(alpha, *refLin, bz) : tr.rotate(alpha);
398 };
399
400 auto accountCluster = [&](int i, std::array<o2::track::TrackParCov, 8>& extrapDest, o2::track::TrackParCov& tr, o2::track::TrackPar* refLin) {
401 if (clArr[i]) { // update with cluster
402 if (!rotateTrack(tr, i == 0 ? pvAlpha : geom->getSensorRefAlpha(clArr[i]->getSensorID()), refLin) ||
403 !prop->propagateTo(tr, refLin, clArr[i]->getX(), true)) {
404 return 0;
405 }
406 extrapDest[i] = tr; // before update
407 if (!tr.update(*clArr[i])) {
408 return 0;
409 }
410 } else {
411 extrapDest[i].invalidate();
412 return -1;
413 }
414 return 1;
415 };
416
417 auto inv2d = [](float s00, float s11, float s01) -> std::array<float, 3> {
418 auto det = s00 * s11 - s01 * s01;
419 if (det < 1e-16) {
420 LOGP(error, "Singular det {}, input: {} {} {}", det, s00, s11, s01);
421 return {0.f, 0.f, 0.f};
422 }
423 det = 1.f / det;
424 return {s11 * det, s00 * det, -s01 * det};
425 };
426
427 resTrack.points.clear();
428 if (!prop->propagateToDCA(pv, trFitOut, bz)) {
429 LOGP(debug, "Failed to propagateToDCA, {}", trFitOut.asString());
430 return false;
431 }
433 if (params.addPVAsCluster) {
434 float cosAlp, sinAlp;
435 pvAlpha = trFitOut.getAlpha();
436 o2::math_utils::sincos(trFitOut.getAlpha(), sinAlp, cosAlp); // vertex position rotated to track frame
437 bcPV.setXYZ(pv.getX() * cosAlp + pv.getY() * sinAlp, -pv.getX() * sinAlp + pv.getY() * cosAlp, pv.getZ());
438 bcPV.setSigmaY2(0.5 * (pv.getSigmaX2() + pv.getSigmaY2()));
439 bcPV.setSigmaZ2(pv.getSigmaZ2());
440 bcPV.setSensorID(-1);
441 clArr[0] = &bcPV;
442 }
443 // collect all track clusters to array, placing them to layer+1 slot
444 int nCl = iTrack.getNClusters();
445 for (int i = 0; i < nCl; i++) { // clusters are ordered from the outermost to the innermost
446 const auto& curClu = mITSClustersArray[itsClRefs[iTrack.getClusterEntry(i)]];
447
448 int llr = geom->getLayer(curClu.getSensorID());
449 if (clArr[1 + llr]) {
450 LOGP(error, "Cluster at lr {} was already assigned, old sens {}, new sens {}", llr, clArr[1 + llr]->getSensorID(), curClu.getSensorID());
451 }
452 clArr[1 + geom->getLayer(curClu.getSensorID())] = &curClu;
453 }
454 o2::track::TrackPar refLinInw0, refLinOut0, *refLinOut = nullptr, *refLinInw = nullptr;
455 o2::track::TrackPar refLinIBOut0, refLinOBInw0, *refLinOBInw = nullptr, *refLinIBOut = nullptr;
456 if (params.useStableRef) {
457 refLinOut = &(refLinOut0 = trFitOut);
458 refLinInw = &(refLinInw0 = trFitInw);
459 }
460 trFitOut.resetCovariance();
461 trFitOut.setCov(trFitOut.getQ2Pt() * trFitOut.getQ2Pt() * trFitOut.getCov()[14], 14);
462 trFitInw.resetCovariance();
463 trFitInw.setCov(trFitInw.getQ2Pt() * trFitInw.getQ2Pt() * trFitInw.getCov()[14], 14);
464 // fit in inward and outward direction
465 for (int i = 0; i <= 7; i++) {
466 int resOut, resInw;
467 // process resOut in ascending order (0-->7) and resInw in descending order (7-->0)
468 if (!(resOut = accountCluster(i, extrapOut, trFitOut, refLinOut)) || !(resInw = accountCluster(7 - i, extrapInw, trFitInw, refLinInw))) {
469 return false;
470 }
471 // at layer 3, find the IB track (trIBOut) and the OB track (trOBInw)
472 // propagate both trcaks to a common radius, RCompIBOB (12cm), and rotates
473 // them to the same reference frame for comparison
474 if (i == 3 && resOut == 1 && resInw == 1 && params.doIBOB && nCl == 7) {
475 resTrack.trIBOut = trFitOut; // outward track updated at outermost IB layer
476 resTrack.trOBInw = trFitInw; // inward track updated at innermost OB layer
477 o2::track::TrackPar refLinIBOut0, refLinIBIn0;
478 if (refLinOut) {
479 refLinIBOut = &(refLinIBOut0 = refLinOut0);
480 refLinOBInw = &(refLinOBInw0 = refLinInw0);
481 }
482 float xRref;
483 if (!resTrack.trOBInw.getXatLabR(params.rCompIBOB, xRref, bz) ||
484 !prop->propagateTo(resTrack.trOBInw, refLinOBInw, xRref, true) ||
485 !rotateTrack(resTrack.trOBInw, resTrack.trOBInw.getPhiPos(), refLinOBInw) || // propagate OB track to ref R and rotate
486 !rotateTrack(resTrack.trIBOut, resTrack.trOBInw.getAlpha(), refLinIBOut) ||
487 !prop->propagateTo(resTrack.trIBOut, refLinIBOut, resTrack.trOBInw.getX(), true)) { // rotate OB track to same frame and propagate to same X
488 // if any propagation or rotation steps fail, invalidate both tracks
489 return false;
490 }
491 }
492 }
493
494 bool innerDone = false;
495 if (params.doResid) {
496 for (int i = 0; i <= 7; i++) {
497 if (clArr[i]) {
498 // calculate interpolation as a weighted mean of inward/outward extrapolations to this layer
499 const auto &tInw = extrapInw[i], &tOut = extrapOut[i];
500 auto wInw = inv2d(tInw.getSigmaY2(), tInw.getSigmaZ2(), tInw.getSigmaZY());
501 auto wOut = inv2d(tOut.getSigmaY2(), tOut.getSigmaZ2(), tOut.getSigmaZY());
502 if (wInw[0] == 0.f || wOut[0] == 0.f) {
503 return false;
504 }
505 std::array<float, 3> wTot = {wInw[0] + wOut[0], wInw[1] + wOut[1], wInw[2] + wOut[2]};
506 auto cTot = inv2d(wTot[0], wTot[1], wTot[2]);
507 auto ywi = wInw[0] * tInw.getY() + wInw[2] * tInw.getZ() + wOut[0] * tOut.getY() + wOut[2] * tOut.getZ();
508 auto zwi = wInw[2] * tInw.getY() + wInw[1] * tInw.getZ() + wOut[2] * tOut.getY() + wOut[1] * tOut.getZ();
509 auto yw = ywi * cTot[0] + zwi * cTot[2];
510 auto zw = ywi * cTot[2] + zwi * cTot[1];
511 // posCl.push_back(clArr[i]->getXYZGlo(*o2::its::GeometryTGeo::Instance()));
512 auto phi = i == 0 ? tInw.getPhi() : tInw.getPhiPos();
513 o2::math_utils::bringTo02Pi(phi);
514 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);
515 if (!innerDone) {
516 resTrack.track = tInw;
517 innerDone = true;
518 }
519 } else {
520 LOGP(debug, "No cluster on lr {}", i);
521 }
522 }
523 }
524 return true;
525}
526
527bool CheckResidSpec::refitPV(o2::dataformats::PrimaryVertex& pv, int vid)
528{
530 std::vector<o2::track::TrackParCov> tracks;
531 std::vector<bool> useTrack;
532 std::vector<GTrackID> gidsITS;
533 int ntr = pv.getNContributors(), ntrIni = ntr;
534 tracks.reserve(ntr);
535 useTrack.reserve(ntr);
536 gidsITS.reserve(ntr);
537 const auto& vtref = mRecoData->getPrimaryVertexMatchedTrackRefs()[vid];
538 auto trackIndex = mRecoData->getPrimaryVertexMatchedTracks();
539 const auto& itsContRefs = vtref.getITSGloContributors();
540 if (params.useITSGloContributors && itsContRefs.getEntries()) {
541 int itr = itsContRefs.getFirstEntry(), itLim = itr + itsContRefs.getEntries();
542 for (; itr < itLim; itr++) {
543 auto tid = trackIndex[itr]; // these are ITS tracks, the Source part is substituted by the PID!!!
544 tracks.emplace_back().setPID(tid.getSource());
545 gidsITS.emplace_back(tid.getIndex(), GTrackID::ITS);
546 }
547 } else {
548 int itr = vtref.getFirstEntry(), itLim = itr + vtref.getEntries();
549 for (; itr < itLim; itr++) {
550 auto tid = trackIndex[itr];
551 if (tid.isPVContributor() && mRecoData->isTrackSourceLoaded(tid.getSource())) {
552 tracks.emplace_back().setPID(mRecoData->getTrackParam(tid).getPID());
553 gidsITS.push_back(mRecoData->getITSContributorGID(tid));
554 }
555 }
556 }
557 ntr = tracks.size();
558 useTrack.resize(ntr);
559#ifdef WITH_OPENMP
560#pragma omp parallel for schedule(dynamic) num_threads(mNThreads)
561#endif
562 for (int itr = 0; itr < ntr; itr++) {
563 if (!(useTrack[itr] = refitITStrack(tracks[itr], gidsITS[itr]))) {
564 tracks[itr] = mRecoData->getTrackParam(gidsITS[itr]); // this track will not be used but participates in prepareVertexRefit
565 }
566 }
567 ntr = 0;
568 for (auto v : useTrack) {
569 ntr++;
570 }
571 if (ntr < params.minPVContributors || !mVertexer.prepareVertexRefit(tracks, pv)) {
572 LOGP(warn, "Abandon vertex refit: NcontribNew = {} vs NcontribOld = {}", ntr, ntrIni);
573 return false;
574 }
575 LOGP(debug, "Original vtx: Nc:{} {}, chi2={}", pv.getNContributors(), pv.asString(), pv.getChi2());
576 auto pvSave = pv;
577 pv = mVertexer.refitVertexFull(useTrack, pv);
578 LOGP(debug, "Refitted vtx: Nc:{} {}, chi2={}", ntr, pv.asString(), pv.getChi2());
579 if (pv.getChi2() < 0.f) {
580 LOGP(warn, "Failed to refit PV {}", pvSave.asString());
581 return false;
582 }
583 return true;
584}
585
586bool CheckResidSpec::refitITStrack(o2::track::TrackParCov& track, GTrackID gid)
587{
588 // destination tack might have non-default PID assigned
589 const auto& trkITS = mRecoData->getITSTrack(gid);
590 const auto itsClRefs = mRecoData->getITSTracksClusterRefs();
591 const auto& params = CheckResidConfig::Instance();
592 auto pid = track.getPID();
593 track = trkITS.getParamOut();
594 track.resetCovariance();
595 track.setCov(track.getQ2Pt() * track.getQ2Pt() * track.getCov()[14], 14);
596 track.setPID(pid);
597 auto nCl = trkITS.getNumberOfClusters();
599 auto prop = o2::base::Propagator::Instance();
600 float bz = prop->getNominalBz();
601 o2::track::TrackPar refLin{track};
602
603 for (int iCl = 0; iCl < nCl; iCl++) { // clusters are stored from outer to inner layers
604 const auto& cls = mITSClustersArray[itsClRefs[trkITS.getClusterEntry(iCl)]];
605 auto alpha = geom->getSensorRefAlpha(cls.getSensorID());
606 if (!(params.useStableRef ? track.rotate(alpha, refLin, bz) : track.rotate(alpha)) ||
607 !prop->propagateTo(track, params.useStableRef ? &refLin : nullptr, cls.getX(), true)) {
608 LOGP(debug, "refitITStrack failed on propagation to cl#{}, alpha={}, x={} | {}", iCl, alpha, cls.getX(), track.asString());
609 return false;
610 }
611 if (!track.update(cls)) {
612 LOGP(debug, "refitITStrack failed on update with cl#{}, | {}", iCl, track.asString());
613 return false;
614 }
615 }
616 return true;
617}
618
619void CheckResidSpec::fillHistos(const o2::checkresid::Track& trc)
620{
621 const auto& params = CheckResidConfig::Instance();
622 int np = trc.points.size();
623 auto pt = trc.track.getPt();
624 if (pt < params.minPt || pt > params.maxPt) {
625 return;
626 }
627 for (int ip = 0; ip < np; ip++) {
628 const auto& pnt = trc.points[ip];
629 int il = pnt.lr >= 0 ? pnt.lr + 1 : 0;
630 mHMan->getHisto2F(il * 10 + 0 * 100)->Fill(pnt.phi, pnt.dy);
631 mHMan->getHisto2F(il * 10 + 0 * 100 + 1000)->Fill(pnt.z, pnt.dy);
632 mHMan->getHisto2F(il * 10 + 0 * 100 + 2000)->Fill(pt, pnt.dy);
633 mHMan->getHisto2F(il * 10 + 0 * 100 + 3000)->Fill(trc.track.getTgl(), pnt.dy);
634 if (pnt.sig2y > 0) {
635 auto pull = pnt.dy / std::sqrt(pnt.sig2y);
636 mHMan->getHisto2F(il * 10 + 0 * 100 + 5)->Fill(pnt.phi, pull);
637 mHMan->getHisto2F(il * 10 + 0 * 100 + 5 + 1000)->Fill(pnt.z, pull);
638 mHMan->getHisto2F(il * 10 + 0 * 100 + 5 + 2000)->Fill(pt, pull);
639 mHMan->getHisto2F(il * 10 + 0 * 100 + 5 + 3000)->Fill(trc.track.getTgl(), pull);
640 }
641 mHMan->getHisto2F(il * 10 + 1 * 100)->Fill(pnt.phi, pnt.dz);
642 mHMan->getHisto2F(il * 10 + 1 * 100 + 1000)->Fill(pnt.z, pnt.dz);
643 mHMan->getHisto2F(il * 10 + 1 * 100 + 2000)->Fill(pt, pnt.dz);
644 mHMan->getHisto2F(il * 10 + 1 * 100 + 3000)->Fill(trc.track.getTgl(), pnt.dz);
645 if (pnt.sig2z > 0) {
646 auto pull = pnt.dz / std::sqrt(pnt.sig2z);
647 mHMan->getHisto2F(il * 10 + 1 * 100 + 5)->Fill(pnt.phi, pull);
648 mHMan->getHisto2F(il * 10 + 1 * 100 + 5 + 1000)->Fill(pnt.z, pull);
649 mHMan->getHisto2F(il * 10 + 1 * 100 + 5 + 2000)->Fill(pt, pull);
650 mHMan->getHisto2F(il * 10 + 1 * 100 + 5 + 3000)->Fill(trc.track.getTgl(), pull);
651 }
652 }
653 //--------------
654 if (trc.trIBOut.getX() > 1 && std::abs(trc.trIBOut.getX() - trc.trOBInw.getX()) < 0.1) {
655 for (int ip = 0; ip < 5; ip++) {
656 float d = trc.trIBOut.getParam(ip) - trc.trOBInw.getParam(ip);
657 mHMan->getHisto2F(10000 + ip * 10)->Fill(trc.trIBOut.getPhiPos(), d);
658 mHMan->getHisto2F(11000 + ip * 10)->Fill(trc.trIBOut.getZ(), d);
659 mHMan->getHisto2F(12000 + ip * 10)->Fill(pt, d);
660 mHMan->getHisto2F(13000 + ip * 10)->Fill(trc.track.getTgl(), d);
661 float sg = trc.trIBOut.getCovarElem(ip, ip) + trc.trOBInw.getCovarElem(ip, ip);
662 if (sg > 0) {
663 auto pull = d / std::sqrt(sg);
664 mHMan->getHisto2F(10000 + ip * 10 + 5)->Fill(trc.trIBOut.getPhiPos(), pull);
665 mHMan->getHisto2F(11000 + ip * 10 + 5)->Fill(trc.trIBOut.getZ(), pull);
666 mHMan->getHisto2F(12000 + ip * 10 + 5)->Fill(pt, pull);
667 mHMan->getHisto2F(13000 + ip * 10 + 5)->Fill(trc.track.getTgl(), pull);
668 }
669 }
670 }
671}
672
673void CheckResidSpec::bookHistos()
674{
676 mHManV.emplace_back() = std::make_unique<o2::HistoManager>("", fmt::format("{}_hman.root", params.outname));
677 mHMan = mHManV.back().get();
678 mHMan->SetName(params.outname.c_str());
679 auto defLogAxis = [](float xMn, float xMx, int nbin) { // get array for log axis
680 if (xMn <= 0 || xMx <= xMn || nbin < 2) {
681 LOGP(fatal, "Wrong log axis request: xmin = {} xmax = {} nbins = {}", xMn, xMx, nbin);
682 }
683 auto dx = std::log(xMx / xMn) / nbin;
684 std::vector<double> xax(nbin + 1);
685 for (int i = 0; i <= nbin; i++) {
686 xax[i] = xMn * std::exp(dx * i);
687 }
688 return xax;
689 };
690 float minPt = std::max(0.1f, params.minPt), maxPt = std::min(50.f, params.maxPt);
691 auto ptax = defLogAxis(minPt, maxPt, params.nBinsPt);
692
693 for (int il = 0; il < 8; il++) {
694 std::string lrName = il == 0 ? "Vtx" : fmt::format("Lr{}", il - 1);
695 for (int iyz = 0; iyz < 2; iyz++) {
696 std::string dname = iyz == 0 ? "dy" : "dz", dtit = iyz == 0 ? "#DeltaY" : "#DeltaZ";
697 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]);
698 mHMan->addHisto(h2, il * 10 + iyz * 100);
699 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);
700 mHMan->addHisto(h2p, il * 10 + iyz * 100 + 5);
701
702 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]);
703 mHMan->addHisto(hz2, il * 10 + iyz * 100 + 1000);
704 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);
705 mHMan->addHisto(hz2p, il * 10 + iyz * 100 + 5 + 1000);
706
707 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]);
708 mHMan->addHisto(hpt2, il * 10 + iyz * 100 + 2000);
709 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);
710 mHMan->addHisto(hpt2p, il * 10 + iyz * 100 + 5 + 2000);
711
712 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]);
713 mHMan->addHisto(htgl2, il * 10 + iyz * 100 + 3000);
714 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);
715 mHMan->addHisto(htgl2p, il * 10 + iyz * 100 + 5 + 3000);
716 }
717 }
718
719 for (int ip = 0; ip < 5; ip++) {
720 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]);
721 mHMan->addHisto(h2, 10000 + ip * 10);
722 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);
723 mHMan->addHisto(h2p, 10000 + ip * 10 + 5);
724
725 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]);
726 mHMan->addHisto(hz2, 11000 + ip * 10);
727 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);
728 mHMan->addHisto(hz2p, 11000 + ip * 10 + 5);
729
730 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]);
731 mHMan->addHisto(hpt2, 12000 + ip * 10);
732 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);
733 mHMan->addHisto(hpt2p, 12000 + ip * 10 + 5);
734
735 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]);
736 mHMan->addHisto(htgl2, 13000 + ip * 10);
737 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);
738 mHMan->addHisto(htgl2p, 13000 + ip * 10 + 5);
739 }
740 // used PV
741 mHMan->addHisto(new TH1F("pvX", "PV X;x;", params.nBinsPVXYZ, -params.maxHPVXY, params.maxHPVXY), 20000 + 0);
742 mHMan->addHisto(new TH1F("pvY", "PV Y;y;", params.nBinsPVXYZ, -params.maxHPVXY, params.maxHPVXY), 20000 + 1);
743 mHMan->addHisto(new TH1F("pvZ", "PV Z;z;", params.nBinsPVXYZ, -params.maxHPVZ, params.maxHPVZ), 20000 + 2);
744 mHMan->addHisto(new TH1F("pvN", "PV Contributors;Nc;", params.maxHPVN, 1.5, params.maxHPVN + 1.5), 20000 + 3);
745}
746
747void CheckResidSpec::postProcessHistos()
748{
749 printf("Fitting histos\n");
750 if (!mHMan) {
751 if (mHManV.empty()) {
752 LOGP(warn, "nothing to process");
753 return;
754 }
755 mHMan = mHManV[0].get();
756 }
758 auto gs = new TF1("gs", "gaus", -1, 1);
759 int maxH = mPostProcOnly ? mHManV.size() : 1;
760 TObjArray arr;
761 for (int ihm = 0; ihm < maxH; ihm++) {
762 auto* histm = mHManV[ihm].get();
763 auto fitSlices = [&](int id) {
764 auto h2 = histm->getHisto2F(id);
765 if (!h2 || h2->GetEntries() < params.minHistoStat2Fit) {
766 return;
767 }
768 h2->FitSlicesY(gs, 0, -1, 0, "QNR", &arr);
769 arr.SetOwner(true);
770 TH1* hmean = (TH1*)arr.RemoveAt(1);
771 if (hmean) {
772 hmean->SetTitle(Form("<%s>", h2->GetTitle()));
773 histm->addHisto(hmean, id + 1);
774 }
775 TH1* hsig = (TH1*)arr.RemoveAt(2);
776 if (hsig) {
777 hsig->SetTitle(Form("#sigma(%s)", h2->GetTitle()));
778 histm->addHisto(hsig, id + 2);
779 }
780 };
781 for (int ioffs = 0; ioffs <= 3; ioffs++) { // vs phi, Z, pT, tgl
782 int offs = ioffs * 1000;
783 for (int iht = 0; iht < 2; iht++) { // resid, pull
784 int offsV = iht == 0 ? 0 : 5;
785 for (int il = 0; il < 8; il++) {
786 for (int iyz = 0; iyz < 2; iyz++) {
787 fitSlices(il * 10 + iyz * 100 + offsV + offs);
788 }
789 }
790 for (int ip = 0; ip < 5; ip++) {
791 fitSlices(10000 + ip * 10 + offsV + offs);
792 }
793 }
794 }
795 histm->write();
796 }
797 delete gs;
798}
799
800void CheckResidSpec::drawHistos()
801{
802 gROOT->SetBatch(true);
803 gStyle->SetTitleX(0.2);
804 gStyle->SetTitleY(0.88);
805 gStyle->SetTitleW(0.25);
806 gStyle->SetOptStat(0);
807 int nhm = mHManV.size();
808 std::array<unsigned int, 3> hcol{EColor::kRed, EColor::kBlue, EColor::kGreen + 2};
809 std::unique_ptr<TLegend> lg;
810 lg = std::make_unique<TLegend>(0.12, 0.13, 0.9, 0.13 + std::min(0.5f, nhm * 0.2f / 3.f));
811 lg->SetFillStyle(0);
812 lg->SetBorderSize(0);
813 for (int i = 0; i < nhm; i++) {
814 auto hman = mHManV[i].get();
815 if (!hman || hman->GetLast() < 1) {
816 continue;
817 }
818 hman->setMarkerStyle(20 + i + (i % 2) * 4, 0.5);
819 hman->setColor(hcol[i % hcol.size()]);
820 auto le = lg->AddEntry(hman->getHisto(1), hman->GetName(), "lp");
821 le->SetTextColor(hcol[i % hcol.size()]);
822 }
823 TCanvas cly("cly", "", 600, 800), clz("clz", "", 600, 800), clpar("clpar", "", 600, 800);
824 TCanvas czly("czly", "", 600, 800), czlz("czlz", "", 600, 800), czlpar("czlpar", "", 600, 800);
826
827 auto AddLabel = [](const char* txt, float x = 0.1, float y = 0.9, int color = kBlack, float size = 0.04) {
828 TLatex* lt = new TLatex(x, y, txt);
829 lt->SetNDC();
830 lt->SetTextColor(color);
831 lt->SetTextSize(size);
832 lt->Draw();
833 return lt;
834 };
835
836 auto drawResLr = [this](TCanvas& canv, int offs, const float resMM[8], bool logX) {
837 canv.Clear();
838 canv.Divide(2, 4);
839 int nh = this->mHManV.size();
840 for (int i = 0; i < 8; i++) {
841 canv.cd(i + 1);
842 bool same = false;
843 for (int j = 0; j < nh; j++) {
844 auto hman = this->mHManV[j].get();
845 if (!hman || hman->GetLast() < 1) {
846 continue;
847 }
848 if (auto histo = hman->getHisto(10 * i + offs)) {
849 histo->Draw(same ? "same" : "");
850 if (!same) {
851 histo->SetMinimum(-resMM[i]);
852 histo->SetMaximum(resMM[i]);
853 same = true;
854 }
855 }
856 }
857 gPad->SetGrid();
858 gPad->SetLogx(logX);
859 }
860 };
861
862 auto drawResPar = [this](TCanvas& canv, int offs, const float resMM[8], bool logX) {
863 canv.Clear();
864 canv.Divide(2, 3);
865 int nh = this->mHManV.size();
866 for (int i = 0; i < 5; i++) {
867 canv.cd(i + 1);
868 bool same = false;
869 for (int j = 0; j < nh; j++) {
870 auto hman = this->mHManV[j].get();
871 if (!hman || hman->GetLast() < 1) {
872 continue;
873 }
874 if (auto histo = hman->getHisto(10 * i + offs)) {
875 histo->Draw(same ? "same" : "");
876 if (!same) {
877 histo->SetMinimum(-resMM[i]);
878 histo->SetMaximum(resMM[i]);
879 same = true;
880 }
881 }
882 }
883 gPad->SetGrid();
884 gPad->SetLogx(logX);
885 }
886 };
887
888 cly.Print(Form("%s_hman.pdf[", params.outname.c_str()));
889 drawResLr(cly, 1, params.resMMLrY, false);
890 cly.cd(2);
891 lg->Draw();
892 AddLabel("Y residuals", 0.1, 0.95);
893 cly.Print(Form("%s_hman.pdf", params.outname.c_str()));
894
895 drawResLr(clz, 101, params.resMMLrZ, false);
896 clz.cd(2);
897 lg->Draw();
898 AddLabel("Z residuals", 0.1, 0.95);
899 clz.Print(Form("%s_hman.pdf", params.outname.c_str()));
900
901 drawResLr(czly, 1001, params.resMMLrY, false);
902 czly.cd(2);
903 lg->Draw();
904 AddLabel("Y residuals", 0.1, 0.95);
905 czly.Print(Form("%s_hman.pdf", params.outname.c_str()));
906
907 drawResLr(czlz, 1101, params.resMMLrZ, false);
908 czlz.cd(2);
909 lg->Draw();
910 AddLabel("Z residuals", 0.1, 0.95);
911 czlz.Print(Form("%s_hman.pdf", params.outname.c_str()));
912
913 drawResLr(czly, 2001, params.resMMLrY, true);
914 czly.cd(2);
915 lg->Draw();
916 AddLabel("Y residuals", 0.1, 0.95);
917 czly.Print(Form("%s_hman.pdf", params.outname.c_str()));
918
919 drawResLr(czlz, 2101, params.resMMLrZ, true);
920 czlz.cd(2);
921 lg->Draw();
922 AddLabel("Z residuals", 0.1, 0.95);
923 czlz.Print(Form("%s_hman.pdf", params.outname.c_str()));
924
925 drawResLr(czly, 3001, params.resMMLrY, false);
926 czly.cd(2);
927 lg->Draw();
928 AddLabel("Y residuals", 0.1, 0.95);
929 czly.Print(Form("%s_hman.pdf", params.outname.c_str()));
930
931 drawResLr(czlz, 3101, params.resMMLrZ, false);
932 czlz.cd(2);
933 lg->Draw();
934 AddLabel("Z residuals", 0.1, 0.95);
935 czlz.Print(Form("%s_hman.pdf", params.outname.c_str()));
936
937 drawResPar(clpar, 10001, params.resMMPar, false);
938 clpar.cd(6);
939 lg->Draw();
940 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
941 clpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
942
943 drawResPar(czlpar, 11001, params.resMMPar, false);
944 czlpar.cd(6);
945 lg->Draw();
946 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
947 czlpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
948
949 drawResPar(czlpar, 12001, params.resMMPar, true);
950 czlpar.cd(6);
951 lg->Draw();
952 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
953 czlpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
954
955 drawResPar(czlpar, 13001, params.resMMPar, false);
956 czlpar.cd(6);
957 lg->Draw();
958 AddLabel("IB-OB tracks params differences at R = 12 cm", 0.2, 0.8);
959 czlpar.Print(Form("%s_hman.pdf", params.outname.c_str()));
960
961 cly.Print(Form("%s_hman.pdf]", params.outname.c_str()));
962}
963
965{
966 mDBGOut.reset();
967 if (mHManV.size()) {
968 postProcessHistos();
969 }
970 if (mDraw) {
971 drawHistos();
972 }
973}
974
976{
978 return;
979 }
980 if (matcher == ConcreteDataMatcher("GLO", "MEANVERTEX", 0)) {
981 LOG(info) << "Imposing new MeanVertex: " << ((const o2::dataformats::MeanVertexObject*)obj)->asString();
982 mMeanVtx = *(const o2::dataformats::MeanVertexObject*)obj;
983 mMeanVertexUpdated = true;
984 return;
985 }
986 if (matcher == ConcreteDataMatcher("ITS", "CLUSDICT", 0)) {
987 LOG(info) << "cluster dictionary updated";
988 mITSDict = (const o2::itsmft::TopologyDictionary*)obj;
989 return;
990 }
991}
992
993DataProcessorSpec getCheckResidSpec(GTrackID::mask_t srcTracks, GTrackID::mask_t srcClusters, bool drawOnly, bool postProcOnly)
994{
995 std::vector<OutputSpec> outputs;
996 auto dataRequest = std::make_shared<DataRequest>();
997 if (!drawOnly && !postProcOnly) {
998 bool useMC = false;
999 dataRequest->requestTracks(srcTracks, useMC);
1000 dataRequest->requestClusters(srcClusters, useMC);
1001 dataRequest->requestPrimaryVertices(useMC);
1002 dataRequest->inputs.emplace_back("meanvtx", "GLO", "MEANVERTEX", 0, Lifetime::Condition, ccdbParamSpec("GLO/Calib/MeanVertex", {}, 1));
1003 }
1004 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
1005 true, // GRPECS=true
1006 true, // GRPLHCIF
1007 true, // GRPMagField
1008 true, // askMatLUT
1010 dataRequest->inputs, true);
1011 Options opts;
1012 if (!drawOnly) {
1013 opts = Options{
1014 {"nthreads", VariantType::Int, 1, {"number of threads"}},
1015 {"no-tree", VariantType::Bool, false, {"do not fill residuals tree"}},
1016 {"no-hist", VariantType::Bool, false, {"do not fill residuals histograms"}},
1017 {"draw-report", VariantType::Bool, false, {"fill residuals report"}},
1018 };
1019 }
1020
1021 return DataProcessorSpec{
1022 "check-resid",
1023 dataRequest->inputs,
1024 outputs,
1025 AlgorithmSpec{adaptFromTask<CheckResidSpec>(dataRequest, ggRequest, srcTracks, drawOnly, postProcOnly)},
1026 opts};
1027}
1028
1029} // 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
TH1 * getHisto(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
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
o2::dataformats::GlobalTrackID GTrackID
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::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"