Project
Loading...
Searching...
No Matches
TrackingStudy.cxx
Go to the documentation of this file.
1// Copyright 2019-2026 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
12#include <vector>
13#include <cmath>
14
15#include <TStopwatch.h>
16#include <TF1.h>
17#include <Eigen/Dense>
18
32#include "Framework/Task.h"
34#include "ITS3Base/SpecsV2.h"
40#include "ITS3Align/TrackFit.h"
49#include "Framework/Logger.h"
50
51namespace o2::its3::study
52{
53
54using namespace o2::framework;
60using T2VMap = std::unordered_map<GTrackID, size_t>;
61
62class TrackingStudySpec : public Task
63{
64 public:
69 TrackingStudySpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, GTrackID::mask_t src, bool useMC)
70 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(src), mUseMC(useMC) {}
71 ~TrackingStudySpec() final = default;
72 void init(InitContext& ic) final;
73 void run(ProcessingContext& pc) final;
74 void endOfStream(EndOfStreamContext& ec) final;
75 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
76
77 private:
78 void process();
79 void updateTimeDependentParams(ProcessingContext& pc);
80 void prepareITSClusters();
81 bool selectTrack(GTrackID trkID, bool checkMCTruth = true) const;
82 T2VMap buildT2V(bool includeCont = false, bool requireMCMatch = true) const;
83 bool refitITSPVTrack(o2::track::TrackParCov& trFit, GTrackID gidx);
84
85 void doDCAStudy();
86 void doDCARefitStudy();
87 void doPullStudy();
88 void doMCStudy();
89 void doResidStudy();
90 void doMisalignmentStudy();
91
92 struct TrackCounter {
93 TrackCounter() = default;
94
95 void operator+=(int src)
96 {
97 if (src >= 0 && src < static_cast<int>(mSuccess.size())) {
98 ++mSuccess[src];
99 }
100 }
101
102 void operator-=(int src)
103 {
104 if (src >= 0 && src < static_cast<int>(mirrors.size())) {
105 ++mirrors[src];
106 }
107 }
108
109 void operator&=(int src)
110 {
111 if (src >= 0 && src < static_cast<int>(mRejected.size())) {
112 ++mRejected[src];
113 }
114 }
115
116 void print() const
117 {
118 LOGP(info, "\t\t\tSuccess / Error / Rejected");
119 for (int cis = 0; cis < GTrackID::NSources; ++cis) {
120 const auto cdm = GTrackID::getSourceDetectorsMask(cis);
121 if (cdm[DetID::ITS]) {
122 LOGP(info, "\t{:{}}\t{} / {} / {}", GTrackID::getSourceName(cis), 15, mSuccess[cis], mirrors[cis], mRejected[cis]);
123 }
124 }
125 }
126
127 void reset()
128 {
129 mSuccess.fill(0);
130 mirrors.fill(0);
131 mRejected.fill(0);
132 }
133
134 std::array<size_t, GTrackID::NSources> mSuccess{};
135 std::array<size_t, GTrackID::NSources> mirrors{};
136 std::array<size_t, GTrackID::NSources> mRejected{};
137 };
138 TrackCounter mTrackCounter;
139
140 using TrackingCluster = align::TrackingCluster<float>;
141 std::vector<TrackingCluster> mITScl;
142 std::span<const int> mITSclRef;
143
144 const ITS3TrackingStudyParam* mParams{nullptr};
145 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
146 std::shared_ptr<DataRequest> mDataRequest;
147 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
148 bool mUseMC{false};
149 GTrackID::mask_t mTracksSrc;
150 o2::vertexing::PVertexer mVertexer;
151 o2::steer::MCKinematicsReader mMCReader; // reader of MC information
152 const o2::its3::TopologyDictionary* mITSDict = nullptr; // cluster patterns dictionary
154 align::MisalignmentModel mMisalignment;
155};
156
158{
160 int lane = ic.services().get<const o2::framework::DeviceSpec>().inputTimesliceId;
161 int maxLanes = ic.services().get<const o2::framework::DeviceSpec>().maxInputTimeslices;
162 std::string dbgnm = maxLanes == 1 ? "its3TrackStudy.root" : fmt::format("its3TrackStudy_{}.root", lane);
163 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>(dbgnm.c_str(), "recreate");
164
166 LOGP(fatal, "initialization of MCKinematicsReader failed");
167 }
168}
169
171{
172 mRecoData.collectData(pc, *mDataRequest);
173 updateTimeDependentParams(pc);
174 process();
175}
176
177void TrackingStudySpec::updateTimeDependentParams(ProcessingContext& pc)
178{
180 if (static bool initOnceDone{false}; !initOnceDone) { // this params need to be queried only once
181 initOnceDone = true;
183 mVertexer.init();
186 if (mParams->doMisalignment) {
187 mMisalignment = {};
188 if (!mParams->misAlgJson.empty()) {
189 mMisalignment = align::loadMisalignmentModel(mParams->misAlgJson);
190 }
191 }
192 }
193}
194
196{
197 mDBGOut.reset();
198}
199
201{
203 return;
204 }
205 if (matcher == ConcreteDataMatcher("IT3", "CLUSDICT", 0)) {
206 LOG(info) << "cluster dictionary updated";
207 mITSDict = (const o2::its3::TopologyDictionary*)obj;
208 return;
209 }
210}
211
212void TrackingStudySpec::process()
213{
214 prepareITSClusters();
215 if (mParams->doDCA) {
216 doDCAStudy();
217 }
218 if (mParams->doDCARefit) {
219 doDCARefitStudy();
220 }
221 if (mUseMC && mParams->doPull) {
222 doPullStudy();
223 }
224 if (mUseMC && mParams->doMC) {
225 doMCStudy();
226 }
227 if (mParams->doResid) {
228 doResidStudy();
229 }
230 if (mUseMC && mParams->doMisalignment) {
231 doMisalignmentStudy();
232 }
233}
234
235void TrackingStudySpec::prepareITSClusters()
236{
237 const auto& clusITS = mRecoData.getITSClusters();
238 LOGP(info, "Preparing {} measurments", clusITS.size());
239 const auto& patterns = mRecoData.getITSClustersPatterns();
240 mITScl.reserve(clusITS.size());
241 auto pattIt = patterns.begin();
242 auto geom = its::GeometryTGeo::Instance();
243 mITSclRef = mRecoData.getITSTracksClusterRefs();
244 mITScl.clear();
245 mITScl.reserve(clusITS.size());
246 for (const auto& cls : clusITS) {
247 const auto sens = cls.getSensorID();
248 float sigmaY2{0}, sigmaZ2{0};
249 math_utils::Point3D<float> locXYZ = o2::its3::ioutils::extractClusterData(cls, pattIt, mITSDict, sigmaY2, sigmaZ2);
250 // Transformation to the local --> global
251 const auto gloXYZ = geom->getMatrixL2G(sens) * locXYZ;
252 // Inverse transformation to the local --> tracking
253 o2::math_utils::Point3D<float> trkXYZ = geom->getMatrixT2L(sens) ^ locXYZ;
254 // Tracking alpha angle
255 // We want that each cluster rotates its tracking frame to the clusters phi
256 // that way the track linearization around the measurement is less biases to the arc
257 // this means automatically that the measurement on the arc is at 0 for the curved layers
258 float alpha = geom->getSensorRefAlpha(sens);
259 if (constants::detID::isDetITS3(sens)) {
260 trkXYZ.SetY(0.f);
261 // alpha&x always have to be defined wrt to the global Z axis!
262 trkXYZ.SetX(std::hypot(gloXYZ.x(), gloXYZ.y()));
263 alpha = std::atan2(gloXYZ.y(), gloXYZ.x());
264 }
265 auto& cl3d = mITScl.emplace_back(sens, trkXYZ);
266 cl3d.setErrors(sigmaY2, sigmaZ2, 0.f);
267 cl3d.alpha = alpha;
268 math_utils::detail::bringToPMPi(cl3d.alpha); // alpha is defined on -Pi,Pi
269 }
270}
271
272bool TrackingStudySpec::selectTrack(GTrackID trkID, bool checkMCTruth) const
273{
274 if (!trkID.includesDet(GTrackID::ITS)) {
275 return false;
276 }
277 if (!mRecoData.isTrackSourceLoaded(trkID.getSource())) {
278 return false;
279 }
280 auto contributorsGID = mRecoData.getSingleDetectorRefs(trkID);
281 if (!contributorsGID[GTrackID::ITS].isIndexSet()) { // we need of course ITS
282 return false;
283 }
284 // ITS specific
285 const auto& itsTrk = mRecoData.getITSTrack(contributorsGID[GTrackID::ITS]);
286 if (itsTrk.getChi2() > mParams->maxChi2 || itsTrk.getNClusters() < mParams->minITSCls) {
287 return false;
288 }
289 // TPC specific
290 if (contributorsGID[GTrackID::TPC].isIndexSet()) {
291 const auto& tpcTrk = mRecoData.getTPCTrack(contributorsGID[GTrackID::TPC]);
292 if (tpcTrk.getNClusters() < mParams->minTPCCls) {
293 return false;
294 }
295 }
296 // general
297 const auto& gTrk = mRecoData.getTrackParam(trkID);
298 if (gTrk.getPt() < mParams->minPt || gTrk.getPt() > mParams->maxPt) {
299 return false;
300 }
301 if (std::abs(gTrk.getEta()) > mParams->maxEta) {
302 return false;
303 }
304 if (mUseMC && checkMCTruth) {
305 const auto& itsLbl = mRecoData.getTrackMCLabel(contributorsGID[GTrackID::ITS]);
306 if (!itsLbl.isValid()) {
307 return false;
308 }
309 if (contributorsGID[GTrackID::TPC].isIndexSet()) {
310 const auto& tpcLbl = mRecoData.getTrackMCLabel(contributorsGID[GTrackID::TPC]);
311 if (itsLbl != tpcLbl) {
312 return false;
313 }
314 }
315 if (contributorsGID[GTrackID::TRD].isIndexSet()) {
316 // TODO
317 }
318 if (contributorsGID[GTrackID::TOF].isIndexSet()) {
319 const auto& tofLbls = mRecoData.getTOFClustersMCLabels()->getLabels(contributorsGID[GTrackID::TOF]);
320 for (const auto& lbl : tofLbls) {
321 if (lbl.isValid()) {
322 return true;
323 }
324 }
325 }
326 }
327 return true;
328}
329
330T2VMap TrackingStudySpec::buildT2V(bool includeCont, bool requireMCMatch) const
331{
332 // build track->vertex assoc., maybe including contributor tracks
333 auto pvvec = mRecoData.getPrimaryVertices();
334 auto trackIndex = mRecoData.getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
335 auto vtxRefs = mRecoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
336 auto nv = vtxRefs.size() - 1; // last entry is for unassigned tracks, ignore them
337 T2VMap t2v;
338 for (size_t iv = 0; iv < nv; ++iv) {
339 const auto& pv = pvvec[iv];
340 if (pv.getNContributors() - 1 < mParams->minPVCont) {
341 continue;
342 }
343 if (requireMCMatch) {
344 auto pvl = mRecoData.getPrimaryVertexMCLabel(iv);
345 }
346 const auto& vtxRef = vtxRefs[iv];
347 int it = vtxRef.getFirstEntry(), itLim = it + vtxRef.getEntries();
348 for (; it < itLim; it++) {
349 const auto& tvid = trackIndex[it];
350 if (tvid.isAmbiguous()) {
351 continue;
352 }
353 if (!mRecoData.isTrackSourceLoaded(tvid.getSource())) {
354 continue;
355 }
356 if (mUseMC && requireMCMatch) {
357 const auto& pvlbl = mRecoData.getPrimaryVertexMCLabel(iv);
358 if (pvlbl.getEventID() != mRecoData.getTrackMCLabel(tvid).getEventID()) {
359 continue;
360 }
361 }
362 t2v[tvid] = iv;
363 if (includeCont) {
364 auto contributorsGID = mRecoData.getSingleDetectorRefs(tvid);
365 for (int cis = 0; cis < GTrackID::NSources; cis++) {
366 const auto cdm = GTrackID::getSourceDetectorsMask(cis);
367 if (!mRecoData.isTrackSourceLoaded(cis) || !cdm[DetID::ITS] || !contributorsGID[cis].isIndexSet()) {
368 continue;
369 }
370 if (mUseMC && requireMCMatch) {
371 const auto& pvlbl = mRecoData.getPrimaryVertexMCLabel(iv);
372 if (pvlbl.getEventID() != mRecoData.getTrackMCLabel(contributorsGID[cis]).getEventID()) {
373 continue;
374 }
375 }
376 t2v[contributorsGID[cis]] = iv;
377 }
378 }
379 }
380 }
381 return std::move(t2v);
382}
383
384bool TrackingStudySpec::refitITSPVTrack(o2::track::TrackParCov& trFit, GTrackID gidx)
385{
386 if (gidx.getSource() != GTrackID::ITS) {
387 return false;
388 }
389 static auto pvvec = mRecoData.getPrimaryVertices();
390 static auto t2v = buildT2V(true, true);
391 static std::vector<unsigned int> itsTracksROF;
392 if (static bool done{false}; !done) {
393 done = true;
394 const auto& itsTracksROFRec = mRecoData.getITSTracksROFRecords();
395 itsTracksROF.resize(mRecoData.getITSTracks().size());
396 for (unsigned irf = 0, cnt = 0; irf < itsTracksROFRec.size(); irf++) {
397 int ntr = itsTracksROFRec[irf].getNEntries();
398 for (int itr = 0; itr < ntr; itr++) {
399 itsTracksROF[cnt++] = irf;
400 }
401 }
402 }
403 auto prop = o2::base::Propagator::Instance();
404 std::array<const TrackingCluster*, 8> clArr{nullptr};
405 const auto trkIn = mRecoData.getTrackParam(gidx);
406 const auto trkOut = mRecoData.getTrackParamOut(gidx);
407 const auto& itsTrOrig = mRecoData.getITSTrack(gidx);
408 int ncl = itsTrOrig.getNumberOfClusters(), rof = itsTracksROF[gidx.getIndex()];
409 const auto& itsTrackClusRefs = mRecoData.getITSTracksClusterRefs();
410 int clEntry = itsTrOrig.getFirstClusterEntry();
411 const auto propagator = o2::base::Propagator::Instance();
412 // convert PV to a fake cluster in the track DCA frame
413 const auto& pv = pvvec[t2v[gidx]];
414 auto trkPV = trkIn;
415 if (!prop->propagateToDCA(pv, trkPV, prop->getNominalBz(), 2.0, mParams->CorrType)) {
416 mTrackCounter -= gidx.getSource();
417 return false;
418 }
419 // create base cluster from the PV, with the alpha corresponding to the track at DCA
420 float cosAlp = NAN, sinAlp = NAN;
421 o2::math_utils::sincos(trkPV.getAlpha(), sinAlp, cosAlp);
422 // vertex position rotated to track frame
423 TrackingCluster pvCls;
424 pvCls.setXYZ((pv.getX() * cosAlp) + (pv.getY() * sinAlp), (-pv.getX() * sinAlp) + (pv.getY() * cosAlp), pv.getZ());
425 pvCls.setSigmaY2(0.5f * (pv.getSigmaX2() + pv.getSigmaY2()));
426 pvCls.setSigmaZ2(pv.getSigmaZ2());
427 clArr[0] = &pvCls;
428 for (int icl = 0; icl < ncl; ++icl) { // ITS clusters are referred in layer decreasing order
429 clArr[ncl - icl] = &mITScl[itsTrackClusRefs[clEntry + icl]];
430 }
431 // start refit
432 trFit = trkOut;
433 trFit.resetCovariance(1'000);
434 float chi2{0};
435 for (int icl = ncl; icl >= 0; --icl) { // go backwards
436 if (!trFit.rotate(clArr[icl]->alpha) || !prop->propagateToX(trFit, clArr[icl]->getX(), prop->getNominalBz(), 0.85, 2.0, mParams->CorrType)) {
437 mTrackCounter -= gidx.getSource();
438 return false;
439 }
440 chi2 += trFit.getPredictedChi2(*clArr[icl]);
441 if (!trFit.update(*clArr[icl])) {
442 mTrackCounter -= gidx.getSource();
443 return false;
444 }
445 }
446 // chi2 < conf.maxChi2; should I cut here?
447 return true;
448};
449
450void TrackingStudySpec::doDCAStudy()
451{
453 LOGP(info, "Doing DCA study");
454 mTrackCounter.reset();
455 auto prop = o2::base::Propagator::Instance();
456 TStopwatch sw;
457 sw.Start();
458 int nDCAFits{0}, nDCAFitsFail{0};
459 auto pvvec = mRecoData.getPrimaryVertices();
460 auto trackIndex = mRecoData.getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
461 auto vtxRefs = mRecoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
462 auto nv = vtxRefs.size() - 1; // last entry is for unassigned tracks, ignore them
463 auto& stream = (*mDBGOut) << "dca";
464 for (int iv = 0; iv < nv; iv++) {
465 const auto& pv = pvvec[iv];
466 const auto& vtref = vtxRefs[iv];
467 for (int is = 0; is < GTrackID::NSources; is++) {
468 const auto dm = GTrackID::getSourceDetectorsMask(is);
469 if (!mRecoData.isTrackSourceLoaded(is) || !dm[DetID::ITS]) {
470 mTrackCounter &= is;
471 continue;
472 }
473 int idMin = vtref.getFirstEntryOfSource(is), idMax = idMin + vtref.getEntriesOfSource(is);
474 for (int i = idMin; i < idMax; i++) {
475 const auto vid = trackIndex[i];
476 if (!vid.isPVContributor()) {
477 mTrackCounter &= vid.getSource();
478 continue;
479 }
480
481 // we fit each different sub-track type, that include ITS, e.g.
482 // ITS,ITS-TPC,ITS-TPC-TRD,ITS-TPC-TOF,ITS-TPC-TRD-TOF
483 auto contributorsGID = mRecoData.getSingleDetectorRefs(vid);
484 for (int cis = 0; cis < GTrackID::NSources && cis <= is; cis++) {
485 const auto cdm = GTrackID::getSourceDetectorsMask(cis);
486 if (!mRecoData.isTrackSourceLoaded(cis) || !cdm[DetID::ITS] || !contributorsGID[cis].isIndexSet()) {
487 mTrackCounter &= cis;
488 continue;
489 }
490 if (!selectTrack(contributorsGID[cis])) {
491 mTrackCounter &= vid.getSource();
492 continue;
493 }
494
495 o2::dataformats::DCA dcaInfo;
496 const auto& trk = mRecoData.getTrackParam(contributorsGID[cis]);
497 auto trkRefit = trk;
498 // for ITS standalone tracks instead of having the trk at the pv we refit with the pv
499 if (mParams->refitITS && cis == GTrackID::ITS && !refitITSPVTrack(trkRefit, contributorsGID[cis])) {
500 mTrackCounter -= cis;
501 continue;
502 } else {
503 trkRefit.invalidate();
504 };
505
506 auto trkDCA = trk;
507 if (!prop->propagateToDCABxByBz(pv, trkDCA, 2.f, mParams->CorrType, &dcaInfo)) {
508 mTrackCounter -= cis;
509 ++nDCAFitsFail;
510 continue;
511 }
512
513 stream << "src=" << cis
514 << "pv=" << pv
515 << "trk=" << trk
516 << "trkRefit=" << trkRefit
517 << "trkAtPV=" << trkDCA
518 << "dca=" << dcaInfo;
519
520 if (mUseMC) {
521 const auto& lbl = mRecoData.getTrackMCLabel(contributorsGID[cis]);
522 lbl.print();
523 o2::dataformats::DCA dcaInfoMC;
524 const auto& eve = mMCReader.getMCEventHeader(lbl.getSourceID(), lbl.getEventID());
526 mcEve.setPos({(float)eve.GetX(), (float)eve.GetY(), (float)eve.GetZ()});
527 auto trkC = trk;
528 if (!prop->propagateToDCABxByBz(mcEve, trkC, 2.f, mParams->CorrType, &dcaInfoMC)) {
529 mTrackCounter -= cis;
530 ++nDCAFitsFail;
531 continue;
532 }
533 const auto& mcTrk = mMCReader.getTrack(lbl);
534 if (mcTrk == nullptr) {
535 LOGP(fatal, "mcTrk is null did selection fail?");
536 }
537 stream << "mcTrk=" << *mcTrk
538 << "dca2MC=" << dcaInfoMC
539 << "lbl=" << lbl;
540 }
541 stream << "\n";
542
543 ++nDCAFits;
544 mTrackCounter += cis;
545 }
546 }
547 }
548 }
549 sw.Stop();
550 LOGP(info, "doDCAStudy: accepted {} fits, failed {} (in {:.2f} seconds)", nDCAFits, nDCAFitsFail, sw.RealTime());
551 mTrackCounter.print();
552}
553
554void TrackingStudySpec::doDCARefitStudy()
555{
557 LOGP(info, "Doing DCARefit study");
558 mTrackCounter.reset();
559 auto prop = o2::base::Propagator::Instance();
560 TStopwatch sw;
561 sw.Start();
562
563 // build track->vertex assoc.
564 auto pvvec = mRecoData.getPrimaryVertices();
565 auto vtxRefs = mRecoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
566 auto nv = vtxRefs.size() - 1; // last entry is for unassigned tracks, ignore them
567 auto t2v = buildT2V();
568 std::vector<std::vector<GTrackID>> v2t;
569 v2t.resize(nv);
570 auto creator = [&](const auto& trk, GTrackID trkID, float _t0, float terr) -> bool {
571 if constexpr (!isBarrelTrack<decltype(trk)>()) {
572 mTrackCounter &= trkID.getSource();
573 return false;
574 }
575 if (!trkID.includesDet(GTrackID::ITS)) {
576 mTrackCounter &= trkID.getSource();
577 return false;
578 }
579 // general
580 if constexpr (isBarrelTrack<decltype(trk)>()) {
581 if (trk.getPt() < mParams->minPt || trk.getPt() > mParams->maxPt) {
582 mTrackCounter &= trkID.getSource();
583 return false;
584 }
585 if (std::abs(trk.getEta()) > mParams->maxEta) {
586 mTrackCounter &= trkID.getSource();
587 return false;
588 }
589 if (!t2v.contains(trkID)) {
590 mTrackCounter &= trkID.getSource();
591 return false;
592 }
593 if (!selectTrack(trkID, mUseMC)) {
594 mTrackCounter &= trkID.getSource();
595 return false;
596 }
597 }
598 v2t[t2v[trkID]].push_back(trkID);
599 return true;
600 };
601 mRecoData.createTracksVariadic(creator);
602
603 int nDCAFits{0}, nDCAFitsFail{0};
604 auto& stream = (*mDBGOut) << "dcaRefit";
605 for (size_t iv = 0; iv < nv; ++iv) {
606 const auto& pv = pvvec[iv];
607 const auto& trkIDs = v2t[iv];
608 if (trkIDs.size() - 1 < mParams->minPVCont) {
609 continue;
610 }
611 std::vector<o2::track::TrackParCov> trks;
612 trks.reserve(trkIDs.size());
613 for (const auto& trkID : trkIDs) {
614 trks.push_back(mRecoData.getTrackParam(trkID));
615 }
616
617 if (!mVertexer.prepareVertexRefit(trks, pv)) {
618 continue;
619 }
620 std::vector<bool> trkMask(trkIDs.size(), true);
621 for (size_t it{0}; it < trkMask.size(); ++it) {
622 trkMask[it] = false; // mask current track from pv refit
623 if (it != 0) {
624 trkMask[it - 1] = true; // unmask previoustrack from pv refit
625 }
626 auto pvRefit = mVertexer.refitVertex(trkMask, pv);
627 if (pvRefit.getChi2() < 0) {
628 trkMask[it] = true;
629 continue;
630 }
631
632 // check DCA both for refitted and original PV
633 o2::dataformats::DCA dcaInfo;
634 auto trkC = trks[it];
635 if (!prop->propagateToDCABxByBz(pv, trkC, 2.f, mParams->CorrType, &dcaInfo)) {
636 mTrackCounter -= trkIDs[it].getSource();
637 ++nDCAFitsFail;
638 continue;
639 }
640 o2::dataformats::DCA dcaInfoRefit;
641 auto trkCRefit = trks[it];
642 if (!prop->propagateToDCABxByBz(pv, trkCRefit, 2.f, mParams->CorrType, &dcaInfoRefit)) {
643 mTrackCounter -= trkIDs[it].getSource();
644 ++nDCAFitsFail;
645 continue;
646 }
647
648 stream << "src=" << trkIDs[it].getSource()
649 << "pv=" << pv
650 << "trkAtPV=" << trkC
651 << "dca=" << dcaInfo
652 << "pvRefit=" << pvRefit
653 << "trkAtPVRefit=" << trkC
654 << "dcaRefit=" << dcaInfoRefit;
655 if (mUseMC) {
656 const auto& mcTrk = mMCReader.getTrack(mRecoData.getTrackMCLabel(trkIDs[it]));
657 if (mcTrk == nullptr) {
658 LOGP(fatal, "mcTrk is null did selection fail?");
659 }
660 stream << "mcTrk=" << *mcTrk;
661 }
662 stream << "\n";
663 ++nDCAFits;
664 mTrackCounter += trkIDs[it].getSource();
665 }
666 }
667 sw.Stop();
668 LOGP(info, "doDCARefitStudy: accepted {} fits, failed {} (in {:.2f} seconds)", nDCAFits, nDCAFitsFail, sw.RealTime());
669 mTrackCounter.print();
670}
671
672void TrackingStudySpec::doPullStudy()
673{
674 // check track pulls compared to mc generation
675 LOGP(info, "Doing Pull study");
676 mTrackCounter.reset();
677 TStopwatch sw;
678 sw.Start();
679 int nPulls{0}, nPullsFail{0};
680 auto prop = o2::base::Propagator::Instance();
681
682 auto checkInTrack = [&](GTrackID trkID) {
683 if (!selectTrack(trkID)) {
684 mTrackCounter &= trkID.getSource();
685 return;
686 }
687 const auto mcTrk = mMCReader.getTrack(mRecoData.getTrackMCLabel(trkID));
688 if (!mcTrk) {
689 return;
690 }
691 auto trk = mRecoData.getTrackParam(trkID);
692
693 // for ITS standalone tracks we add the PV as an additional measurement point
694 if (mParams->refitITS && trkID.getSource() == GTrackID::ITS && !refitITSPVTrack(trk, trkID)) {
695 mTrackCounter -= trkID.getSource();
696 ++nPullsFail;
697 return;
698 }
699
700 std::array<float, 3> xyz{(float)mcTrk->GetStartVertexCoordinatesX(), (float)mcTrk->GetStartVertexCoordinatesY(), (float)mcTrk->GetStartVertexCoordinatesZ()},
701 pxyz{(float)mcTrk->GetStartVertexMomentumX(), (float)mcTrk->GetStartVertexMomentumY(), (float)mcTrk->GetStartVertexMomentumZ()};
702 TParticlePDG* pPDG = TDatabasePDG::Instance()->GetParticle(mcTrk->GetPdgCode());
703 if (!pPDG) {
704 mTrackCounter -= trkID.getSource();
705 ++nPullsFail;
706 return;
707 }
708 o2::track::TrackPar mcTrkO2(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), false);
709 // propagate it to the alpha/X of the reconstructed track
710 if (!mcTrkO2.rotate(trk.getAlpha()) || !prop->PropagateToXBxByBz(mcTrkO2, trk.getX())) {
711 mTrackCounter -= trkID.getSource();
712 ++nPullsFail;
713 return;
714 }
715 const auto contTrk = mRecoData.getSingleDetectorRefs(trkID);
716 const auto& itsTrk = mRecoData.getITSTrack(contTrk[GTrackID::ITS]);
717
718 (*mDBGOut)
719 << "pull"
720 << "src=" << trkID.getSource()
721 << "itsTrk=" << itsTrk
722 << "mcTrk=" << mcTrkO2
723 << "mcPart=" << mcTrk
724 << "trk=" << trk
725 << "\n";
726 ++nPulls;
727 mTrackCounter += trkID.getSource();
728 };
729
730 for (size_t iTrk{0}; iTrk < mRecoData.getITSTracks().size(); ++iTrk) {
731 checkInTrack(GTrackID(iTrk, GTrackID::ITS));
732 }
733 for (size_t iTrk{0}; iTrk < mRecoData.getTPCITSTracks().size(); ++iTrk) {
734 checkInTrack(GTrackID(iTrk, GTrackID::ITSTPC));
735 }
736 for (size_t iTrk{0}; iTrk < mRecoData.getITSTPCTRDTracksMCLabels().size(); ++iTrk) {
737 checkInTrack(GTrackID(iTrk, GTrackID::ITSTPCTRD));
738 }
739 for (size_t iTrk{0}; iTrk < mRecoData.getITSTPCTOFMatches().size(); ++iTrk) {
740 checkInTrack(GTrackID(iTrk, GTrackID::ITSTPCTOF));
741 }
742 for (size_t iTrk{0}; iTrk < mRecoData.getITSTPCTRDTOFMatches().size(); ++iTrk) {
743 checkInTrack(GTrackID(iTrk, GTrackID::ITSTPCTRDTOF));
744 }
745 sw.Stop();
746 LOGP(info, "doPullStudy: accepted {} pulls; rejected {} (in {:.2f} seconds)", nPulls, nPullsFail, sw.RealTime());
747 mTrackCounter.print();
748}
749
750void TrackingStudySpec::doMCStudy()
751{
752 LOGP(info, "Doing MC study");
753 mTrackCounter.reset();
754 TStopwatch sw;
755 sw.Start();
756 int nTracks{0};
757
758 const int iSrc{0};
759 const int nev = mMCReader.getNEvents(iSrc);
760 std::unordered_map<o2::MCCompLabel, ParticleInfoExt> info;
761
762 LOGP(info, "** Filling particle table ... ");
763 for (int iEve{0}; iEve < nev; ++iEve) {
764 const auto& mcTrks = mMCReader.getTracks(iSrc, iEve);
765 for (int iTrk{0}; iTrk < mcTrks.size(); ++iTrk) {
766 const auto& mcTrk = mcTrks[iTrk];
767 const auto pdg = mcTrk.GetPdgCode();
768 if (o2::O2DatabasePDG::Instance()->GetParticle(pdg) == nullptr) {
769 continue;
770 }
771 const auto apdg = std::abs(pdg);
772 if (apdg != 11 && apdg != 211 && apdg != 321 && apdg != 2212) {
773 continue;
774 }
775 o2::MCCompLabel lbl(iTrk, iEve, iSrc);
776 auto& part = info[lbl];
777 part.mcTrack = mcTrk;
778 }
779 }
780 LOGP(info, "** Creating particle/clusters correspondence ... ");
781 const auto& clusters = mRecoData.getITSClusters();
782 const auto& clustersMCLCont = mRecoData.getITSClustersMCLabels();
783 for (auto iCluster{0}; iCluster < clusters.size(); ++iCluster) {
784 auto labs = clustersMCLCont->getLabels(iCluster);
785 for (auto& lab : labs) {
786 if (!lab.isValid() || lab.getSourceID() != 0 || !lab.isCorrect()) {
787 continue;
788 }
789 int trackID = 0, evID = 0, srcID = 0;
790 bool fake = false;
791 lab.get(trackID, evID, srcID, fake);
792 auto& cluster = clusters[iCluster];
793 auto layer = o2::its::GeometryTGeo::Instance()->getLayer(cluster.getSensorID());
794 auto& part = info[{trackID, evID, srcID}];
795 part.clusters |= (1 << layer);
796 if (fake) {
797 part.fakeClusters |= (1 << layer);
798 }
799 }
800 }
801 LOGP(info, "** Analysing tracks ... ");
802 auto accountLbl = [&](const globaltracking::RecoContainer::GlobalIDSet& contributorsGID, DetID::ID det) {
803 if (contributorsGID[det].isIndexSet()) {
804 const auto& lbl = mRecoData.getTrackMCLabel(contributorsGID[det]);
805 if (lbl.isValid()) {
806 o2::MCCompLabel iLbl(lbl.getTrackID(), lbl.getEventID(), lbl.getSourceID());
807 if (info.contains(iLbl)) {
808 auto& part = info[iLbl];
809 SETBIT(part.recoTracks, det);
810 if (lbl.isFake()) {
811 SETBIT(part.fakeTracks, det);
812 }
813 }
814 }
815 }
816 };
817 auto creator = [&](const auto& trk, GTrackID trkID, float _t0, float terr) -> bool {
818 if constexpr (!isBarrelTrack<decltype(trk)>()) {
819 return false;
820 }
821 if (!trkID.includesDet(GTrackID::ITS)) {
822 return false;
823 }
824 // general
825 auto contributorsGID = mRecoData.getSingleDetectorRefs(trkID);
826 if (!contributorsGID[GTrackID::ITS].isIndexSet()) { // we need of course ITS
827 return false;
828 }
829 const auto& gLbl = mRecoData.getTrackMCLabel(trkID);
830 if (!gLbl.isValid()) {
831 return false;
832 }
833 o2::MCCompLabel iLbl(gLbl.getTrackID(), gLbl.getEventID(), gLbl.getSourceID());
834 if (!info.contains(iLbl)) {
835 return false;
836 }
837 auto& part = info[iLbl];
838 part.recoTrack = mRecoData.getTrackParam(trkID);
839
840 accountLbl(contributorsGID, DetID::ITS);
841 accountLbl(contributorsGID, DetID::TPC);
842 accountLbl(contributorsGID, DetID::TRD);
843 accountLbl(contributorsGID, DetID::TOF);
844
845 ++nTracks;
846 return true;
847 };
848 mRecoData.createTracksVariadic(creator);
849
850 LOGP(info, "Streaming output to tree");
851 for (const auto& [_, part] : info) {
852 (*mDBGOut) << "mc"
853 << "part=" << part
854 << "\n";
855 }
856
857 sw.Stop();
858 LOGP(info, "doMCStudy: accounted {} MCParticles and {} tracks (in {:.2f} seconds)", info.size(), nTracks, sw.RealTime());
859}
860
861void TrackingStudySpec::doResidStudy()
862{
863 LOGP(info, "Doing residual study");
864 const auto geom = o2::its::GeometryTGeo::Instance();
865 const auto prop = o2::base::Propagator::Instance();
866 const float bz = prop->getNominalBz();
867
868 int goodRefit{0}, notPassedSel{0}, fitFail{0};
869
870 auto doRefits = [&](const o2::its::TrackITS& iTrack, const o2::MCCompLabel& lbl) {
871 std::array<TrackingCluster, 8> cl;
872 std::array<const TrackingCluster*, 8> clArr{nullptr};
873 if (mParams->addPVAsCluster) {
874 const auto& eve = mMCReader.getMCEventHeader(lbl.getSourceID(), lbl.getEventID());
876 auto trFitOut = iTrack.getParamIn();
877 pv.setXYZ(eve.GetX(), eve.GetY(), eve.GetZ());
878 if (!prop->propagateToDCA(pv, trFitOut, bz, base::Propagator::MAX_STEP, mParams->CorrType)) {
879 return;
880 }
881 pv.setSigmaX(20e-4f);
882 pv.setSigmaY(20e-4f);
883 pv.setSigmaZ(20e-4f);
884 float cosAlp = NAN, sinAlp = NAN;
885 o2::math_utils::sincos(trFitOut.getAlpha(), sinAlp, cosAlp);
886 cl[0].alpha = trFitOut.getAlpha();
887 cl[0].setXYZ((pv.getX() * cosAlp) + (pv.getY() * sinAlp), (-pv.getX() * sinAlp) + (pv.getY() * cosAlp), pv.getZ());
888 cl[0].setSigmaY2(0.5f * (pv.getSigmaX2() + pv.getSigmaY2()));
889 cl[0].setSigmaZ2(pv.getSigmaZ2());
890 cl[0].setSensorID(-1);
891 clArr[0] = &cl[0];
892 }
893
894 // collect track clusters into layer slots
895 int nCl = iTrack.getNClusters();
896 for (int i = 0; i < nCl; i++) {
897 const auto& curClu = mITScl[mITSclRef[iTrack.getClusterEntry(i)]];
898 int sens = curClu.getSensorID();
899 int llr = geom->getLayer(sens);
900 if (clArr[1 + llr]) {
901 LOGP(fatal, "Cluster at lr {} was already assigned, old sens {}, new sens {}", llr, clArr[1 + llr]->getSensorID(), sens);
902 }
903 clArr[1 + llr] = &curClu;
904 }
905
906 std::array<o2::track::TrackParCov, 8> extrapOut, extrapInw;
907 float chi2{0};
908 if (!align::doBidirRefit(iTrack, clArr, extrapOut, extrapInw, chi2, mParams->useStableRef, mParams->CorrType)) {
909 ++fitFail;
910 return;
911 }
912
913 for (int i = 0; i <= 7; i++) {
914 if (clArr[i]) {
915 const auto tInt = align::interpolateTrackParCov(extrapInw[i], extrapOut[i]);
916 if (!tInt.isValid()) {
917 continue;
918 }
919 auto phi = i == 0 ? tInt.getPhi() : tInt.getPhiPos();
920 o2::math_utils::bringTo02Pi(phi);
921 (*mDBGOut) << "res"
922 << "dYInt=" << clArr[i]->getY() - tInt.getY()
923 << "dZInt=" << clArr[i]->getZ() - tInt.getZ()
924 << "dYIn=" << clArr[i]->getY() - extrapInw[i].getY()
925 << "dZIn=" << clArr[i]->getZ() - extrapInw[i].getZ()
926 << "dYOut=" << clArr[i]->getY() - extrapOut[i].getY()
927 << "dZOut=" << clArr[i]->getZ() - extrapOut[i].getZ()
928 << "chi2=" << chi2
929 << "clY=" << clArr[i]->getY()
930 << "clZ=" << clArr[i]->getZ()
931 << "clX=" << clArr[i]->getX()
932 << "alpha=" << clArr[i]->alpha
933 << "sens=" << clArr[i]->getSensorID()
934 << "phi=" << phi
935 << "pt=" << tInt.getPt()
936 << "chip=" << constants::detID::getSensorID(clArr[i]->getSensorID())
937 << "lay=" << i - 1
938 << "\n";
939 }
940 }
941 ++goodRefit;
942 };
943
944 const auto itsTracks = mRecoData.getITSTracks();
945 const auto itsMC = mRecoData.getITSTracksMCLabels();
946 for (size_t iTrk{0}; iTrk < itsTracks.size(); ++iTrk) {
947 const auto& iTrack = itsTracks[iTrk];
948 const auto& lbl = itsMC[iTrk];
949 const auto& mc = mMCReader.getTrack(lbl);
950 if (std::abs(iTrack.getEta()) > mParams->maxEta || iTrack.getChi2() > mParams->maxChi2 || iTrack.getNClusters() < mParams->minITSCls || iTrack.getPt() < mParams->minPt || !lbl.isCorrect() || !mc->isPrimary()) {
951 ++notPassedSel;
952 continue;
953 }
954 doRefits(iTrack, lbl);
955 }
956
957 LOGP(info, "\trefitted {} out of {} tracks ({} !sel, {} !fit)", goodRefit, itsTracks.size(), notPassedSel, fitFail);
958}
959
960void TrackingStudySpec::doMisalignmentStudy()
961{
962 LOGP(info, "Doing misalignment study");
963 const auto prop = o2::base::Propagator::Instance();
964 const auto geom = o2::its::GeometryTGeo::Instance();
965
966 int goodRefit{0}, notPassedSel{0}, fitFail{0}, fitFailMis{0};
967
968 float chi2{0};
969 auto writeTree = [&](const char* treeName,
970 const std::array<const TrackingCluster*, 8>& clArr,
971 const std::array<o2::track::TrackParCov, 8>& extrapOut,
972 const std::array<o2::track::TrackParCov, 8>& extrapInw,
973 const o2::MCCompLabel& lbl) {
974 for (int i = 0; i <= 7; i++) {
975 if (!clArr[i]) {
976 continue;
977 }
978 // interpolated result
979 auto tInt = align::interpolateTrackParCov(extrapInw[i], extrapOut[i]);
980 if (!tInt.isValid()) {
981 continue;
982 }
983 float dY = clArr[i]->getY() - tInt.getY();
984 float dZ = clArr[i]->getZ() - tInt.getZ();
985 // MC truth at same (alpha, x)
986 o2::track::TrackPar mcTrkAtX;
987 const auto mcTrk = mMCReader.getTrack(lbl);
988 if (mcTrk) {
989 std::array<float, 3> xyz{(float)mcTrk->GetStartVertexCoordinatesX(), (float)mcTrk->GetStartVertexCoordinatesY(), (float)mcTrk->GetStartVertexCoordinatesZ()};
990 std::array<float, 3> pxyz{(float)mcTrk->GetStartVertexMomentumX(), (float)mcTrk->GetStartVertexMomentumY(), (float)mcTrk->GetStartVertexMomentumZ()};
991 TParticlePDG* pPDG = TDatabasePDG::Instance()->GetParticle(mcTrk->GetPdgCode());
992 if (pPDG) {
993 mcTrkAtX = o2::track::TrackPar(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), false);
994 if (mcTrkAtX.rotate(tInt.getAlpha()) && prop->PropagateToXBxByBz(mcTrkAtX, tInt.getX())) {
995 auto phi = i == 0 ? tInt.getPhi() : tInt.getPhiPos();
996 o2::math_utils::bringTo02Pi(phi);
997 (*mDBGOut) << treeName
998 << "trk=" << tInt
999 << "mcTrk=" << mcTrkAtX
1000 << "chi2=" << chi2
1001 << "dY=" << dY
1002 << "dZ=" << dZ
1003 << "phi=" << phi
1004 << "eta=" << tInt.getEta()
1005 << "lay=" << i - 1
1006 << "\n";
1007 }
1008 }
1009 }
1010 }
1011 };
1012
1013 const auto itsTracks = mRecoData.getITSTracks();
1014 const auto itsMC = mRecoData.getITSTracksMCLabels();
1015 for (size_t iTrk{0}; iTrk < itsTracks.size(); ++iTrk) {
1016 const auto& iTrack = itsTracks[iTrk];
1017 if (std::abs(iTrack.getEta()) > mParams->maxEta || iTrack.getChi2() > mParams->maxChi2 || iTrack.getNClusters() < mParams->minITSCls || iTrack.getPt() < mParams->minPt) {
1018 ++notPassedSel;
1019 continue;
1020 }
1021 const auto& lbl = itsMC[iTrk];
1022 if (!lbl.isCorrect() || !lbl.isValid()) {
1023 ++notPassedSel;
1024 continue;
1025 }
1026 const auto& mc = mMCReader.getTrack(lbl);
1027 if (!mc->isPrimary()) {
1028 ++notPassedSel;
1029 continue;
1030 }
1031
1032 // ideal clusters
1033 std::array<TrackingCluster, 8> cl;
1034 std::array<const TrackingCluster*, 8> clArr{nullptr};
1035 if (mParams->addPVAsCluster) {
1036 const auto& eve = mMCReader.getMCEventHeader(lbl.getSourceID(), lbl.getEventID());
1038 auto trFitOut = iTrack.getParamIn();
1039 pv.setXYZ(eve.GetX(), eve.GetY(), eve.GetZ());
1040 if (!prop->propagateToDCA(pv, trFitOut, prop->getNominalBz(), base::Propagator::MAX_STEP, mParams->CorrType)) {
1041 return;
1042 }
1043 pv.setSigmaX(20e-4f);
1044 pv.setSigmaY(20e-4f);
1045 pv.setSigmaZ(20e-4f);
1046 float cosAlp = NAN, sinAlp = NAN;
1047 o2::math_utils::sincos(trFitOut.getAlpha(), sinAlp, cosAlp);
1048 cl[0].alpha = trFitOut.getAlpha();
1049 cl[0].setXYZ((pv.getX() * cosAlp) + (pv.getY() * sinAlp), (-pv.getX() * sinAlp) + (pv.getY() * cosAlp), pv.getZ());
1050 cl[0].setSigmaY2(0.5f * (pv.getSigmaX2() + pv.getSigmaY2()));
1051 cl[0].setSigmaZ2(pv.getSigmaZ2());
1052 cl[0].setSensorID(-1);
1053 clArr[0] = &cl[0];
1054 }
1055
1056 // collect track clusters into layer slots
1057 int nCl = iTrack.getNClusters();
1058 for (int i = 0; i < nCl; i++) {
1059 const auto& curClu = mITScl[mITSclRef[iTrack.getClusterEntry(i)]];
1060 int sens = curClu.getSensorID();
1061 int llr = geom->getLayer(sens);
1062 if (clArr[1 + llr]) {
1063 LOGP(fatal, "Cluster at lr {} was already assigned, old sens {}, new sens {}", llr, clArr[1 + llr]->getSensorID(), sens);
1064 }
1065 clArr[1 + llr] = &curClu;
1066 }
1067 std::array<o2::track::TrackParCov, 8> extrapOut, extrapInw;
1068 chi2 = 0;
1069 if (!align::doBidirRefit(iTrack, clArr, extrapOut, extrapInw, chi2, mParams->useStableRef, mParams->CorrType)) {
1070 ++fitFail;
1071 continue;
1072 }
1073 writeTree("idealRes", clArr, extrapOut, extrapInw, lbl);
1074
1075 // Propagate MC truth to each cluster's (alpha, x) to get true track direction.
1076 // The shared misalignment evaluators then provide the tracking-frame dy/dz shift.
1077 const auto mcTrk = mMCReader.getTrack(lbl);
1078 if (!mcTrk) {
1079 continue;
1080 }
1081 std::array<float, 3> xyz{(float)mcTrk->GetStartVertexCoordinatesX(), (float)mcTrk->GetStartVertexCoordinatesY(), (float)mcTrk->GetStartVertexCoordinatesZ()};
1082 std::array<float, 3> pxyz{(float)mcTrk->GetStartVertexMomentumX(), (float)mcTrk->GetStartVertexMomentumY(), (float)mcTrk->GetStartVertexMomentumZ()};
1083 TParticlePDG* pPDG = TDatabasePDG::Instance()->GetParticle(mcTrk->GetPdgCode());
1084 if (!pPDG) {
1085 continue;
1086 }
1087 o2::track::TrackPar mcPar(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), false);
1088
1089 std::array<TrackingCluster, 3> misClArr; // shifted copies for up to 3 IT3 layers
1090 std::array<const TrackingCluster*, 8> clArrMis{};
1091 for (int i = 0; i <= 7; i++) {
1092 clArrMis[i] = clArr[i]; // PV and OB clusters stay the same
1093 }
1094 for (int iLay = 0; iLay < 3; ++iLay) {
1095 if (!clArr[1 + iLay]) {
1096 continue;
1097 }
1098 const auto& orig = *clArr[1 + iLay];
1099 const int sens = orig.getSensorID();
1100 if (!constants::detID::isDetITS3(sens)) {
1101 continue;
1102 }
1103 const int sensorID = constants::detID::getSensorID(sens);
1104 const int layerID = constants::detID::getDetID2Layer(sens);
1105 const auto& sensorMis = mMisalignment[sensorID];
1106
1107 // propagate MC track to cluster's tracking frame to get true slopes
1108 auto mcAtCl = mcPar;
1109 if (!mcAtCl.rotate(orig.alpha) || !prop->PropagateToXBxByBz(mcAtCl, orig.getX())) {
1110 clArrMis[1 + iLay] = nullptr; // can't compute slopes -> drop cluster
1111 continue;
1112 }
1113 const align::MisalignmentFrame misFrame{
1114 .sensorID = sensorID,
1115 .layerID = layerID,
1116 .x = orig.getX(),
1117 .alpha = orig.alpha,
1118 .z = orig.getZ()};
1119 const auto slopes = align::computeTrackSlopes(mcAtCl.getSnp(), mcAtCl.getTgl());
1120
1121 align::MisalignmentShift totalShift;
1122 if (sensorMis.hasLegendre) {
1123 const auto shift = align::evaluateLegendreShift(sensorMis, misFrame, slopes);
1124 if (!shift.accepted) {
1125 clArrMis[1 + iLay] = nullptr; // shifted outside acceptance
1126 continue;
1127 }
1128 totalShift += shift;
1129 }
1130 if (sensorMis.hasInextensional) {
1131 totalShift += align::evaluateInextensionalShift(sensorMis, misFrame, slopes);
1132 }
1133
1134 // create shifted copy: keep x=r (nominal), shift y and z
1135 misClArr[iLay] = orig;
1136 misClArr[iLay].setY(orig.getY() + totalShift.dy);
1137 misClArr[iLay].setZ(orig.getZ() + totalShift.dz);
1138 misClArr[iLay].setSigmaY2(orig.getSigmaY2() + (mParams->misAlgExtCY[sensorID] * mParams->misAlgExtCY[sensorID]));
1139 misClArr[iLay].setSigmaZ2(orig.getSigmaZ2() + (mParams->misAlgExtCZ[sensorID] * mParams->misAlgExtCZ[sensorID]));
1140 clArrMis[1 + iLay] = &misClArr[iLay];
1141 }
1142
1143 // refit with shifted clusters
1144 chi2 = 0;
1145 if (!align::doBidirRefit(iTrack, clArrMis, extrapOut, extrapInw, chi2, mParams->useStableRef, mParams->CorrType)) {
1146 ++fitFailMis;
1147 ++goodRefit; // ideal still succeeded
1148 continue;
1149 }
1150 writeTree("misRes", clArrMis, extrapOut, extrapInw, lbl);
1151
1152 ++goodRefit;
1153 }
1154
1155 LOGP(info, "\tdoMisalignmentStudy: refitted {} out of {} tracks ({} !sel, {} !fit, {} !fitMis)", goodRefit, itsTracks.size(), notPassedSel, fitFail, fitFailMis);
1156}
1157
1158DataProcessorSpec getTrackingStudySpec(GTrackID::mask_t srcTracks, GTrackID::mask_t srcClusters, bool useMC, bool withPV)
1159{
1160 std::vector<OutputSpec> outputs;
1161 auto dataRequest = std::make_shared<DataRequest>();
1162
1163 dataRequest->requestTracks(srcTracks, useMC);
1164 dataRequest->requestIT3Clusters(useMC);
1165 dataRequest->requestClusters(srcClusters, useMC);
1166 if (withPV) {
1167 dataRequest->requestPrimaryVertices(useMC);
1168 }
1169 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
1170 true, // GRPECS=true
1171 true, // GRPLHCIF
1172 true, // GRPMagField
1173 true, // askMatLUT
1175 dataRequest->inputs,
1176 true);
1177
1178 return DataProcessorSpec{
1179 .name = "its3-track-study",
1180 .inputs = dataRequest->inputs,
1181 .outputs = outputs,
1182 .algorithm = AlgorithmSpec{adaptFromTask<TrackingStudySpec>(dataRequest, ggRequest, srcTracks, useMC)},
1183 .options = {}};
1184}
1185
1186} // namespace o2::its3::study
Wrapper container for different reconstructed object types.
Definition of the ITSMFT digit.
Definition of the GeometryManager class.
Definition of the ITSMFT Hit class.
Definition of the BuildTopologyDictionary class for ITS3.
void print() const
int32_t i
bool done
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.
Primary vertex finder.
Wrapper container for different reconstructed object types.
o2::track::TrackParCov TrackParCov
Definition Recon.h:39
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
void print() const
int getEventID() const
Double_t GetStartVertexCoordinatesX() const
Definition MCTrack.h:82
static TDatabasePDG * Instance()
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
static std::string getCollisionContextFileName(const std::string_view prefix="")
Definition NameConf.cxx:52
static constexpr float MAX_STEP
Definition Propagator.h:73
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
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 TRD
Definition DetID.h:65
static constexpr ID TPC
Definition DetID.h:64
static constexpr ID TOF
Definition DetID.h:66
ServiceRegistryRef services()
Definition InitContext.h:34
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
TrackingStudySpec & operator=(TrackingStudySpec &&)=delete
void run(ProcessingContext &pc) final
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
TrackingStudySpec(const TrackingStudySpec &)=delete
void init(InitContext &ic) final
TrackingStudySpec(TrackingStudySpec &&)=delete
TrackingStudySpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, GTrackID::mask_t src, bool useMC)
TrackingStudySpec & operator=(const TrackingStudySpec &)=delete
static GeometryTGeo * Instance()
int getLayer(int index) const final
Get chip layer, from 0.
void fillMatrixCache(int mask) override
int getClusterEntry(int i) const
Definition TrackITS.h:67
bool initFromDigitContext(std::string_view filename)
MCTrack const * getTrack(o2::MCCompLabel const &) const
size_t getNEvents(int source) const
Get number of events.
o2::dataformats::MCEventHeader const & getMCEventHeader(int source, int event) const
retrieves the MCEventHeader for a given eventID and sourceID
std::vector< MCTrack > const & getTracks(int source, int event) const
variant returning all tracks for source and event at once
bool prepareVertexRefit(const TR &tracks, const o2d::VertexBase &vtxSeed)
Definition PVertexer.h:359
PVertex refitVertex(const std::vector< bool > useTrack, const o2d::VertexBase &vtxSeed)
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLenum src
Definition glcorearb.h:1767
GLsizeiptr size
Definition glcorearb.h:659
GLenum GLuint GLint GLint layer
Definition glcorearb.h:1310
GLuint GLuint stream
Definition glcorearb.h:1806
Defining PrimaryVertex explicitly as messageable.
Definition Cartesian.h:288
MisalignmentShift evaluateInextensionalShift(const SensorMisalignment &sensor, const MisalignmentFrame &frame, const TrackSlopes &slopes)
bool doBidirRefit(const o2::its::TrackITS &iTrack, std::array< const TrackingCluster< T > *, 8 > &clArr, std::array< o2::track::TrackParametrizationWithError< T >, 8 > &extrapOut, std::array< o2::track::TrackParametrizationWithError< T >, 8 > &extrapInw, T &chi2, bool useStableRef, typename o2::base::PropagatorImpl< T >::MatCorrType corrType)
Definition TrackFit.h:108
MisalignmentModel loadMisalignmentModel(const std::string &jsonPath)
o2::track::TrackParametrizationWithError< T > interpolateTrackParCov(const o2::track::TrackParametrizationWithError< T > &tA, const o2::track::TrackParametrizationWithError< T > &tB)
Definition TrackFit.h:59
MisalignmentShift evaluateLegendreShift(const SensorMisalignment &sensor, const MisalignmentFrame &frame, const TrackSlopes &slopes)
TrackSlopes computeTrackSlopes(double snp, double tgl)
bool isDetITS3(T detID)
Definition SpecsV2.h:210
o2::math_utils::Point3D< T > extractClusterData(const itsmft::CompClusterExt &c, iterator &iter, const o2::its3::TopologyDictionary *dict, T &sig2y, T &sig2z)
Definition IOUtils.h:37
o2::framework::DataProcessorSpec getTrackingStudySpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool useMC, bool withPV)
std::unordered_map< GTrackID, size_t > T2VMap
o2::dataformats::GlobalTrackID GTrackID
TrackParCovF TrackParCov
Definition Track.h:33
TrackParF TrackPar
Definition Track.h:29
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
GlobalIDSet getSingleDetectorRefs(GTrackID gidx) const
const o2::tpc::TrackTPC & getTPCTrack(GTrackID id) const
void createTracksVariadic(T creator, GTrackID::mask_t srcSel=GTrackID::getSourcesMask("all")) const
o2::MCCompLabel getTrackMCLabel(GTrackID id) const
const o2::track::TrackParCov & getTrackParam(GTrackID gidx) const
const o2::MCEventLabel & getPrimaryVertexMCLabel(int i) const
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
const o2::track::TrackParCov & getTrackParamOut(GTrackID gidx) const
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
std::array< GTrackID, GTrackID::NSources > GlobalIDSet
o2::base::PropagatorImpl< float >::MatCorrType CorrType
static constexpr int L2G
Definition Cartesian.h:54
static constexpr int T2L
Definition Cartesian.h:55
static constexpr int T2G
Definition Cartesian.h:56
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
TStopwatch sw
std::vector< Cluster > clusters