Project
Loading...
Searching...
No Matches
AlignmentSpec.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 <cmath>
13#include <memory>
14#include <chrono>
15
16#ifdef WITH_OPENMP
17#include <omp.h>
18#endif
19
20#include <Eigen/Dense>
21#include <GblTrajectory.h>
22#include <GblData.h>
23#include <GblPoint.h>
24#include <GblMeasurement.h>
25#include <MilleBinary.h>
26#include <nlohmann/json.hpp>
27
30#include "Framework/Task.h"
43#include "ITStracking/IOUtils.h"
45#include "ITS3Align/TrackFit.h"
52
53namespace o2::its3::align
54{
55using namespace o2::framework;
63
64namespace
65{
66// compute normalized (u,v) in [-1,1] from global position on a sensor
67std::pair<double, double> computeUV(double gloX, double gloY, double gloZ, int sensorID, double radius)
68{
69 const bool isTop = sensorID % 2 == 0;
70 const double phi = o2::math_utils::to02Pid(std::atan2(gloY, gloX));
71 const double phiBorder1 = o2::math_utils::to02Pid(((isTop ? 0. : 1.) * TMath::Pi()) + std::asin(constants::equatorialGap / 2. / radius));
72 const double phiBorder2 = o2::math_utils::to02Pid(((isTop ? 1. : 2.) * TMath::Pi()) - std::asin(constants::equatorialGap / 2. / radius));
73 const double u = (((phi - phiBorder1) * 2.) / (phiBorder2 - phiBorder1)) - 1.;
74 const double v = ((2. * gloZ + constants::segment::lengthSensitive) / constants::segment::lengthSensitive) - 1.;
75 return {u, v};
76}
77
78// evaluate Legendre polynomials P_0(x) through P_order(x) via recurrence
79std::vector<double> legendrePols(int order, double x)
80{
81 std::vector<double> p(order + 1);
82 p[0] = 1.;
83 if (order > 0) {
84 p[1] = x;
85 }
86 for (int n = 1; n < order; ++n) {
87 p[n + 1] = ((2 * n + 1) * x * p[n] - n * p[n - 1]) / (n + 1);
88 }
89 return p;
90}
91} // namespace
92
93class AlignmentSpec final : public Task
94{
95 public:
96 ~AlignmentSpec() final = default;
97 AlignmentSpec(const AlignmentSpec&) = delete;
99 AlignmentSpec& operator=(const AlignmentSpec&) = delete;
100 AlignmentSpec& operator=(AlignmentSpec&&) = delete;
101 AlignmentSpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, GTrackID::mask_t src, bool useMC, bool withPV, bool withITS, OutputEnum out)
102 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(src), mUseMC(useMC), mWithPV(withPV), mIsITS3(!withITS), mOutOpt(out)
103 {
104 }
105
106 void init(InitContext& ic) final;
107 void run(ProcessingContext& pc) final;
108 void endOfStream(EndOfStreamContext& ec) final;
109 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
110 void process();
111
112 private:
113 void updateTimeDependentParams(ProcessingContext& pc);
114 void buildHierarchy();
115
116 // calculate the transport jacobian for points FROM and TO numerically via ridder's method
117 // this assumes the track is already at point FROM and will be extrapolated to TO's x (xTo)
118 // method does not modify the original track
119 bool getTransportJacobian(const TrackD& track, double xTo, double alphaTo, gbl::Matrix5d& jac, gbl::Matrix5d& err);
120
121 // refit ITS track with inward/outward fit (opt. impose pv as additional constraint)
122 // after this we have the refitted track at the innermost update point
123 bool prepareITSTrack(int iTrk, const o2::its::TrackITS& itsTrack, Track& resTrack);
124
125 // prepare ITS measuremnt points
126 void prepareMeasurments(std::span<const itsmft::CompClusterExt> clusters, std::span<const unsigned char> pattIt);
127
128 // build track to vertex association
129 void buildT2V();
130
131 // apply some misalignment on inner ITS3 layers
132 // it can happen that a measurement is pushed outside of
133 // ITS3 acceptance so false is to discard track
134 bool applyMisalignment(Eigen::Vector2d& res, const FrameInfoExt& frame, const TrackD& wTrk, size_t iTrk);
135
136 OutputEnum mOutOpt;
137 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
138 std::vector<dataformats::VertexBase> mPVs;
139 std::vector<int> mT2PV;
140 bool mIsITS3{true};
141 const o2::itsmft::TopologyDictionary* mITSDict{nullptr};
142 const o2::its3::TopologyDictionary* mIT3Dict{nullptr};
143 o2::globaltracking::RecoContainer* mRecoData = nullptr;
144 std::unique_ptr<steer::MCKinematicsReader> mcReader;
145 std::vector<FrameInfoExt> mITSTrackingInfo;
146 std::shared_ptr<DataRequest> mDataRequest;
147 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
148 std::unique_ptr<AlignableVolume> mHierarchy; // tree-hiearchy
149 AlignableVolume::SensorMapping mChip2Hiearchy; // global label mapping to leaves in the tree
150 bool mUseMC{false};
151 bool mWithPV{false};
152 GTrackID::mask_t mTracksSrc;
153 int mNThreads{1};
154 const AlignmentParams* mParams{nullptr};
155 std::array<o2::math_utils::Legendre2DPolynominal, 6> mDeformations; // one per sensorID (0-5)
156 std::array<Eigen::Matrix<double, 6, 1>, 6> mRigidBodyParams; // (dx,dy,dz,rx,ry,rz) in LOC per sensorID
157};
158
160{
162 mNThreads = ic.options().get<int>("nthreads");
163 if (mOutOpt) {
164 LOG(info) << mOutOpt.pstring();
165 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>("its3_debug_alg.root", "recreate");
166 }
167 if (mUseMC) {
168 mcReader = std::make_unique<steer::MCKinematicsReader>("collisioncontext.root");
169 }
170}
171
173{
174 if (mOutOpt[OutputOpt::MilleRes]) {
175 updateTimeDependentParams(pc);
176 writeMillepedeResults(mHierarchy.get(), mParams->milleResFile, mParams->milleResOutJson, mParams->misAlgJson);
177 } else {
179 mRecoData = &recoData;
180 mRecoData->collectData(pc, *mDataRequest);
181 updateTimeDependentParams(pc);
182 process();
183 }
184 mRecoData = nullptr;
185}
186
188{
189 if (!mITSDict && !mIT3Dict) {
190 LOGP(fatal, "ITS data is not loaded");
191 }
193 const auto bz = prop->getNominalBz();
194 const auto itsTracks = mRecoData->getITSTracks();
195 const auto itsClRefs = mRecoData->getITSTracksClusterRefs();
196 const auto clusITS = mRecoData->getITSClusters();
197 const auto patterns = mRecoData->getITSClustersPatterns();
198 std::span<const o2::MCCompLabel> mcLbls;
199 if (mUseMC) {
200 mcLbls = mRecoData->getITSTracksMCLabels();
201 }
202 prepareMeasurments(clusITS, patterns);
203
204 if (mWithPV) {
205 buildT2V();
206 }
207
208 LOGP(info, "Starting fits with {} threads", mNThreads);
209
210 // Data
211 std::vector<std::vector<gbl::GblTrajectory>> gblTrajSlots(mNThreads);
212 std::vector<std::vector<Track>> resTrackSlots(mNThreads);
213
214 auto timeStart = std::chrono::high_resolution_clock::now();
215 int cFailedRefit{0}, cFailedProp{0}, cSelected{0}, cGBLFit{0}, cGBLFitFail{0}, cGBLChi2Rej{0}, cGBLConstruct{0};
216 double chi2Sum{0}, lostWeightSum{0};
217 int ndfSum{0};
218#ifdef WITH_OPENMP
219#pragma omp parallel num_threads(mNThreads) \
220 reduction(+ : cFailedRefit) \
221 reduction(+ : cFailedProp) \
222 reduction(+ : cSelected) \
223 reduction(+ : cGBLFit) \
224 reduction(+ : cGBLFitFail) \
225 reduction(+ : cGBLChi2Rej) \
226 reduction(+ : cGBLConstruct) \
227 reduction(+ : chi2Sum) \
228 reduction(+ : lostWeightSum) \
229 reduction(+ : ndfSum)
230#endif
231 {
232#ifdef WITH_OPENMP
233 const int tid = omp_get_thread_num();
234#else
235 const int tid = 0;
236#endif
237 auto& gblTrajSlot = gblTrajSlots[tid];
238 auto& resTrackSlot = resTrackSlots[tid];
239
240#ifdef WITH_OPENMP
241#pragma omp for schedule(dynamic)
242#endif
243 for (size_t iTrk = 0; iTrk < (int)itsTracks.size(); ++iTrk) {
244 const auto& trk = itsTracks[iTrk];
245 if (trk.getNClusters() < mParams->minITSCls ||
246 (trk.getChi2() / ((float)trk.getNClusters() * 2 - 5)) >= mParams->maxITSChi2Ndf ||
247 trk.getPt() < mParams->minPt ||
248 (mUseMC && (!mcLbls[iTrk].isValid() || !mcLbls[iTrk].isCorrect()))) {
249 continue;
250 }
251 ++cSelected;
252 Track& resTrack = resTrackSlot.emplace_back();
253 if (!prepareITSTrack((int)iTrk, trk, resTrack)) {
254 ++cFailedRefit;
255 resTrackSlot.pop_back();
256 continue;
257 }
258
259 o2::track::TrackParD* refLin = nullptr;
260 if (mParams->useStableRef) {
261 refLin = &resTrack.track;
262 }
263
264 // outward stepping from track IU
265 auto wTrk = resTrack.track;
266 const bool hasPV = resTrack.info[0].lr == -1;
267 std::vector<gbl::GblPoint> points;
268 bool failed = false;
269 const int np = (int)resTrack.points.size();
271 lt.setTimeNotNeeded();
272 constexpr int perm[5] = {4, 2, 3, 0, 1}; // ALICE->GBL: Q/Pt,Snp,Tgl,Y,Z
273 for (int ip{0}; ip < np; ++ip) {
274 const auto& frame = resTrack.info[ip];
275 gbl::Matrix5d err = gbl::Matrix5d::Identity(), jacALICE = gbl::Matrix5d::Identity(), jacGBL;
276 float msErr = 0.f;
277 if (ip) {
278 // numerically calculates the transport jacobian from prev. point to this point
279 // then we actually do the step to the point and accumulate the material
280 if (!getTransportJacobian(wTrk, frame.x, frame.alpha, jacALICE, err) ||
281 !prop->propagateToAlphaX(wTrk, refLin, frame.alpha, frame.x, false, mParams->maxSnp, mParams->maxStep, 1, mParams->corrType, &lt)) {
282 ++cFailedProp;
283 failed = true;
284 break;
285 }
286 msErr = its::math_utils::MSangle(trk.getPID().getMass(), trk.getP(), lt.getX2X0());
287 // after computing jac, reorder to GBL convention
288 for (int i = 0; i < 5; i++) {
289 for (int j = 0; j < 5; j++) {
290 jacGBL(i, j) = jacALICE(perm[i], perm[j]);
291 }
292 }
293 }
294
295 // wTrk is now in the measurment frame
296 gbl::GblPoint point(jacGBL);
297 // measurement
298 Eigen::Vector2d res, prec;
299 res << frame.positionTrackingFrame[0] - wTrk.getY(), frame.positionTrackingFrame[1] - wTrk.getZ();
300
301 // here we can apply some misalignment on the measurment
302 if (!applyMisalignment(res, frame, wTrk, iTrk)) {
303 failed = true;
304 break;
305 }
306
307 prec << 1. / resTrack.points[ip].sig2y, 1. / resTrack.points[ip].sig2z;
308 // the projection matrix is in the tracking frame the idendity so no need to diagonalize it
309 point.addMeasurement(res, prec);
310 if (msErr > mParams->minMS && ip < np - 1) {
311 Eigen::Vector2d scat(0., 0.), scatPrec = Eigen::Vector2d::Constant(1. / (msErr * msErr));
312 point.addScatterer(scat, scatPrec);
313 lt.clearFast(); // clear if accounted
314 }
315
316 if (frame.lr >= 0) {
317 GlobalLabel lbl(0, frame.sens, true);
318 if (mChip2Hiearchy.find(lbl) == mChip2Hiearchy.end()) {
319 LOGP(fatal, "Cannot find global label: {}", lbl.asString());
320 }
321
322 // derivatives for all sensitive volumes and their parents
323 // this is the derivative in TRK but we want to align in LOC
324 // so dr/da_(LOC) = dr/da_(TRK) * da_(TRK)/da_(LOC)
325 const auto* tileVol = mChip2Hiearchy.at(lbl);
327
328 // count rigid body columns: only volumes with real DOFs (not DOFPseudo)
329 int nColRB{0};
330 for (const auto* v = tileVol; v && !v->isRoot(); v = v->getParent()) {
331 if (v->getRigidBody()) {
332 nColRB += v->getRigidBody()->nDOFs();
333 }
334 }
335
336 // count calibration columns
337 const auto* sensorVol = tileVol->getParent();
338 const auto* calibSet = sensorVol ? sensorVol->getCalib() : nullptr;
339 const int nCalib = calibSet ? calibSet->nDOFs() : 0;
340 const int nCol = nColRB + nCalib;
341
342 std::vector<int> gLabels;
343 gLabels.reserve(nCol);
344 Eigen::MatrixXd gDer(3, nCol);
345 gDer.setZero();
346 Eigen::Index curCol{0};
347
348 // 1) tile: TRK -> LOC via precomputed T2L and J_L2T
349 const double posTrk[3] = {frame.x, 0., 0.};
350 double posLoc[3];
351 tileVol->getT2L().LocalToMaster(posTrk, posLoc);
352 Matrix66 jacL2T;
353 tileVol->computeJacobianL2T(posLoc, jacL2T);
354 der *= jacL2T;
355 if (tileVol->getRigidBody()) {
356 const int nd = tileVol->getRigidBody()->nDOFs();
357 for (int iDOF = 0; iDOF < nd; ++iDOF) {
358 gLabels.push_back(tileVol->getLabel().rawGBL(iDOF));
359 }
360 gDer.middleCols(curCol, nd) = der;
361 curCol += nd;
362 }
363
364 // 2) chain through parents: child's J_L2P
365 for (const auto* child = tileVol; child->getParent() && !child->getParent()->isRoot(); child = child->getParent()) {
366 der *= child->getJL2P();
367 const auto* parent = child->getParent();
368 if (parent->getRigidBody()) {
369 const int nd = parent->getRigidBody()->nDOFs();
370 for (int iDOF = 0; iDOF < nd; ++iDOF) {
371 gLabels.push_back(parent->getLabel().rawGBL(iDOF));
372 }
373 gDer.middleCols(curCol, nd) = der;
374 curCol += nd;
375 }
376 }
377
378 // 3) calibration derivatives (e.g. Legendre for ITS3 sensors, apply directly on the whole sensor, not on inidividual tiles)
379 if (calibSet && calibSet->type() == DOFSet::Type::Legendre) {
380 const auto* legSet = static_cast<const LegendreDOFSet*>(calibSet);
381 const int N = legSet->order();
382 const int sensorID = constants::detID::getSensorID(frame.sens);
383 const int layerID = constants::detID::getDetID2Layer(frame.sens);
384
385 const double r = frame.x;
386 const double gX = r * std::cos(frame.alpha);
387 const double gY = r * std::sin(frame.alpha);
388 const double gZ = frame.positionTrackingFrame[1];
389 auto [u, v] = computeUV(gX, gY, gZ, sensorID, constants::radii[layerID]);
390
391 const double snp = wTrk.getSnp();
392 const double tgl = wTrk.getTgl();
393 const double csci = 1. / std::sqrt(1. - (snp * snp));
394 const double dydx = snp * csci;
395 const double dzdx = tgl * csci;
396
397 auto pu = legendrePols(N, u);
398 auto pv = legendrePols(N, v);
399
400 int legIdx = 0;
401 const int legColStart = nColRB;
402 for (int i = 0; i <= N; ++i) {
403 for (int j = 0; j <= i; ++j) {
404 const double basis = pu[j] * pv[i - j];
405 gLabels.push_back(sensorVol->getLabel().asCalib().rawGBL(legIdx));
406 gDer(0, legColStart + legIdx) = dydx * basis;
407 gDer(1, legColStart + legIdx) = dzdx * basis;
408 ++legIdx;
409 }
410 }
411 }
412 point.addGlobals(gLabels, gDer);
413 }
414
415 if (mOutOpt[OutputOpt::VerboseGBL]) {
416 static Eigen::IOFormat fmt(4, 0, ", ", "\n", "[", "]");
417 LOGP(info, "WORKING-POINT {}", ip);
418 LOGP(info, "Track: {}", wTrk.asString());
419 LOGP(info, "FrameInfo: {}", frame.asString());
420 std::cout << "jacALICE:\n"
421 << jacALICE.format(fmt) << '\n';
422 std::cout << "jacGBL:\n"
423 << jacGBL.format(fmt) << '\n';
424 LOGP(info, "Point {}: GBL res=({}, {}), KF stored res=({}, {})",
425 ip, res[0], res[1], resTrack.points[ip].dy, resTrack.points[ip].dz);
426 LOGP(info, "residual: dy={} dz={}", res[0], res[1]);
427 LOGP(info, "precision: precY={} precZ={}", prec[0], prec[1]);
428 point.printPoint(5);
429 }
430 points.push_back(point);
431 }
432 if (!failed) {
433 gbl::GblTrajectory traj(points, std::abs(bz) > 0.01);
434 if (traj.isValid()) {
435 double chi2 = NAN, lostWeight = NAN;
436 int ndf = 0;
437 if (auto ierr = traj.fit(chi2, ndf, lostWeight); !ierr) {
438 if (mOutOpt[OutputOpt::VerboseGBL]) {
439 LOGP(info, "GBL FIT chi2 {} ndf {}", chi2, ndf);
440 traj.printTrajectory(5);
441 }
442 if (chi2 / ndf > mParams->maxChi2Ndf && cGBLChi2Rej++ < 10) {
443 LOGP(error, "GBL fit exceeded red chi2 {}", chi2 / ndf);
444 if (std::abs(resTrack.kfFit.chi2Ndf - 1) < 0.02) {
445 LOGP(error, "\tGBL is far away from good KF fit!!!!");
446 continue;
447 }
448 } else {
449 ++cGBLFit;
450 chi2Sum += chi2;
451 lostWeightSum += lostWeight;
452 ndfSum += ndf;
453 if (mOutOpt[OutputOpt::MilleData]) {
454 gblTrajSlot.push_back(traj);
455 }
456 FitInfo fit;
457 fit.ndf = ndf;
458 fit.chi2 = (float)chi2;
459 fit.chi2Ndf = (float)chi2 / (float)ndf;
460 resTrack.gblFit = fit;
461 }
462 } else {
463 ++cGBLFitFail;
464 }
465 } else {
466 ++cGBLConstruct;
467 }
468 }
469 }
470 }
471 auto timeEnd = std::chrono::high_resolution_clock::now();
472 auto duration = std::chrono::duration_cast<std::chrono::milliseconds>(timeEnd - timeStart);
473 LOGP(info, "Fitted {} tracks out of {} (selected {}) in {} sec", cGBLFit, itsTracks.size(), cSelected, duration.count() / 1e3);
474 LOGP(info, "\tRefit failed for {} tracks; Failed prop for {} tracks", cFailedRefit, cFailedProp);
475 LOGP(info, "\tGBL SUMMARY:");
476 LOGP(info, "\t\tGBL construction failed {}", cGBLConstruct);
477 LOGP(info, "\t\tGBL fit failed {}", cGBLFitFail);
478 LOGP(info, "\t\tGBL chi2Ndf rejected {}", cGBLChi2Rej);
479 if (!ndfSum) {
480 LOGP(info, "\t\tGBL Chi2/Ndf = NDF IS 0");
481 } else {
482 LOGP(info, "\t\tGBL Chi2/Ndf = {}", chi2Sum / ndfSum);
483 }
484 LOGP(info, "\t\tGBL LostWeight = {}", lostWeightSum);
485 LOGP(info, "Streaming results to output");
486 if (mOutOpt[OutputOpt::MilleData]) {
487 gbl::MilleBinary mille(mParams->milleBinFile, true);
488 for (auto& slot : gblTrajSlots) {
489 for (auto& traj : slot) {
490 traj.milleOut(mille);
491 }
492 }
493 }
494 if (mOutOpt[OutputOpt::Debug]) {
495 for (auto& slot : resTrackSlots) {
496 for (auto& res : slot) {
497 (*mDBGOut) << "res"
498 << "trk=" << res
499 << "\n";
500 }
501 }
502 }
503}
504
505void AlignmentSpec::updateTimeDependentParams(ProcessingContext& pc)
506{
508 if (static bool initOnce{false}; !initOnce) {
509 initOnce = true;
512 mParams = &AlignmentParams::Instance();
513 mParams->printKeyValues(true, true);
514
515 buildHierarchy();
516
517 if (mParams->doMisalignmentLeg || mParams->doMisalignmentRB) {
518 TMatrixD null(1, 1);
519 null(0, 0) = 0;
520 for (int i = 0; i < 6; ++i) {
521 mDeformations[i] = o2::math_utils::Legendre2DPolynominal(null);
522 mRigidBodyParams[i].setZero();
523 }
524 if (!mParams->misAlgJson.empty()) {
525 using json = nlohmann::json;
526 std::ifstream f(mParams->misAlgJson);
527 auto data = json::parse(f);
528 for (const auto& item : data) {
529 int id = item["id"].get<int>();
530 if (mParams->doMisalignmentLeg && item.contains("matrix")) {
531 auto v = item["matrix"].get<std::vector<std::vector<double>>>();
532 TMatrixD m(v.size(), v[v.size() - 1].size());
533 for (size_t r{0}; r < v.size(); ++r) {
534 for (size_t c{0}; c < v[r].size(); ++c) {
535 m(r, c) = v[r][c];
536 }
537 }
539 }
540 if (mParams->doMisalignmentRB && item.contains("rigidBody")) {
541 auto rb = item["rigidBody"].get<std::vector<double>>();
542 for (int k = 0; k < 6 && k < (int)rb.size(); ++k) {
543 mRigidBodyParams[id](k) = rb[k];
544 }
545 }
546 }
547 }
548 }
549 }
550}
551
552void AlignmentSpec::buildHierarchy()
553{
554 if (mIsITS3) {
555 mHierarchy = buildHierarchyIT3(mChip2Hiearchy);
556 } else {
557 mHierarchy = buildHierarchyITS(mChip2Hiearchy);
558 }
559
560 if (!mParams->dofConfigJson.empty()) {
561 applyDOFConfig(mHierarchy.get(), mParams->dofConfigJson);
562 }
563
564 mHierarchy->finalise();
565 if (mOutOpt[OutputOpt::MilleSteer]) {
566 std::ofstream tree(mParams->milleTreeFile);
567 mHierarchy->writeTree(tree);
568 std::ofstream cons(mParams->milleConFile);
569 mHierarchy->writeRigidBodyConstraints(cons);
570 std::ofstream par(mParams->milleParamFile);
571 mHierarchy->writeParameters(par);
572 }
573}
574
575bool AlignmentSpec::getTransportJacobian(const TrackD& track, double xTo, double alphaTo, gbl::Matrix5d& jac, gbl::Matrix5d& err)
576{
578 const auto bz = prop->getNominalBz();
579 const auto minStep = std::sqrt(std::numeric_limits<double>::epsilon());
580 const gbl::Vector5d x0(track.getParams());
581 auto trackC = track;
582 o2::track::TrackParD* refLin{nullptr};
583 if (mParams->useStableRef) {
584 refLin = &trackC;
585 }
586
587 auto propagate = [&](gbl::Vector5d& p) -> bool {
588 TrackD tmp(track);
589 for (int i{0}; i < track::kNParams; ++i) {
590 tmp.setParam(p[i], i);
591 }
592 if (!prop->propagateToAlphaX(tmp, refLin, alphaTo, xTo, false, mParams->maxSnp, mParams->maxStep, 1, mParams->corrType)) {
593 return false;
594 }
595 p = gbl::Vector5d(tmp.getParams());
596 return true;
597 };
598
599 for (int iPar{0}; iPar < track::kNParams; ++iPar) {
600 // step size
601 double h = std::min(mParams->ridderMaxIniStep[iPar], std::max(minStep, std::abs(track.getParam(iPar)) * mParams->ridderRelIniStep[iPar]) * std::pow(mParams->ridderShrinkFac, mParams->ridderMaxExtrap / 2));
602 ;
603 // romberg tableu
604 Eigen::MatrixXd cur(track::kNParams, mParams->ridderMaxExtrap);
605 Eigen::MatrixXd pre(track::kNParams, mParams->ridderMaxExtrap);
606 double normErr = std::numeric_limits<double>::max();
607 gbl::Vector5d bestDeriv = gbl::Vector5d::Constant(std::numeric_limits<double>::max());
608 for (int iExt{0}; iExt < mParams->ridderMaxExtrap; ++iExt) {
609 gbl::Vector5d xPlus = x0, xMinus = x0;
610 xPlus(iPar) += h;
611 xMinus(iPar) -= h;
612 if (!propagate(xPlus) || !propagate(xMinus)) {
613 return false;
614 }
615 cur.col(0) = (xPlus - xMinus) / (2.0 * h);
616 if (!iExt) {
617 bestDeriv = cur.col(0);
618 }
619 // shrink step in next iteration
620 h /= mParams->ridderShrinkFac;
621 // richardson extrapolation
622 double fac = mParams->ridderShrinkFac * mParams->ridderShrinkFac;
623 for (int k{1}; k <= iExt; ++k) {
624 cur.col(k) = (fac * cur.col(k - 1) - pre.col(k - 1)) / (fac - 1.0);
625 fac *= mParams->ridderShrinkFac * mParams->ridderShrinkFac;
626 double e = std::max((cur.col(k) - cur.col(k - 1)).norm(), (cur.col(k) - pre.col(k - 1)).norm());
627 if (e <= normErr) {
628 normErr = e;
629 bestDeriv = cur.col(k);
630 if (normErr < mParams->ridderEps) {
631 break;
632 }
633 }
634 }
635 if (normErr < mParams->ridderEps) {
636 break;
637 }
638 // check stability
639 if (iExt > 0) {
640 double tableauErr = (cur.col(iExt) - pre.col(iExt - 1)).norm();
641 if (tableauErr >= 2.0 * normErr) {
642 break;
643 }
644 }
645 std::swap(cur, pre);
646 }
647 if (bestDeriv.isApproxToConstant(std::numeric_limits<double>::max())) {
648 return false;
649 }
650 jac.col(iPar) = bestDeriv;
651 err.col(iPar) = gbl::Vector5d::Constant(normErr);
652 }
653
654 if (jac.isIdentity(1e-8)) {
655 LOGP(error, "Near jacobian idendity for taking track from {} to {}", track.getX(), xTo);
656 return false;
657 }
658
659 return true;
660}
661
662bool AlignmentSpec::prepareITSTrack(int iTrk, const o2::its::TrackITS& itsTrack, align::Track& resTrack)
663{
664 const auto itsClRefs = mRecoData->getITSTracksClusterRefs();
665 auto trFit = convertTrack<double>(itsTrack.getParamOut()); // take outer track fit as start of refit
668 const auto bz = prop->getNominalBz();
669 std::array<const FrameInfoExt*, 8> frameArr{};
670 o2::track::TrackParD trkOut, *refLin = nullptr;
671 if (mParams->useStableRef) {
672 refLin = &(trkOut = trFit);
673 }
674
675 auto accountCluster = [&](int i, TrackD& tr, float& chi2, Measurement& meas, o2::track::TrackParD* refLin) {
676 if (frameArr[i]) { // update with cluster
677 if (!prop->propagateToAlphaX(tr, refLin, frameArr[i]->alpha, frameArr[i]->x, false, mParams->maxSnp, mParams->maxStep, 1, mParams->corrType)) {
678 return 2;
679 }
680 meas.dy = frameArr[i]->positionTrackingFrame[0] - tr.getY();
681 meas.dz = frameArr[i]->positionTrackingFrame[1] - tr.getZ();
682 meas.sig2y = frameArr[i]->covarianceTrackingFrame[0];
683 meas.sig2z = frameArr[i]->covarianceTrackingFrame[2];
684 meas.z = tr.getZ();
685 meas.phi = tr.getPhi();
686 o2::math_utils::bringTo02Pid(meas.phi);
687 chi2 += (float)tr.getPredictedChi2Quiet(frameArr[i]->positionTrackingFrame, frameArr[i]->covarianceTrackingFrame);
688 if (!tr.update(frameArr[i]->positionTrackingFrame, frameArr[i]->covarianceTrackingFrame)) {
689 return 2;
690 }
691 if (refLin) { // displace the reference to the last updated cluster
692 refLin->setY(frameArr[i]->positionTrackingFrame[0]);
693 refLin->setZ(frameArr[i]->positionTrackingFrame[1]);
694 }
695 return 0;
696 }
697 return 1;
698 };
699
700 FrameInfoExt pvInfo;
701 if (mWithPV) { // add PV as constraint
702 const int iPV = mT2PV[iTrk];
703 if (iPV < 0) {
704 return false;
705 }
706 const auto& pv = mPVs[iPV];
707 auto tmp = convertTrack<double>(itsTrack.getParamIn());
708 if (!prop->propagateToDCA(pv, tmp, bz)) {
709 return false;
710 }
711 pvInfo.alpha = (float)tmp.getAlpha();
712 double ca{0}, sa{0};
713 o2::math_utils::bringToPMPid(pvInfo.alpha);
714 o2::math_utils::sincosd(pvInfo.alpha, sa, ca);
715 pvInfo.x = tmp.getX();
716 pvInfo.positionTrackingFrame[0] = -pv.getX() * sa + pv.getY() * ca;
717 pvInfo.positionTrackingFrame[1] = pv.getZ();
718 pvInfo.covarianceTrackingFrame[0] = 0.5 * (pv.getSigmaX2() + pv.getSigmaY2());
719 pvInfo.covarianceTrackingFrame[2] = pv.getSigmaY2();
720 pvInfo.sens = -1;
721 pvInfo.lr = -1;
722 frameArr[0] = &pvInfo;
723 }
724
725 // collect all track clusters to array, placing them to layer+1 slot
726 int nCl = itsTrack.getNClusters();
727 for (int i = 0; i < nCl; i++) { // clusters are ordered from the outermost to the innermost
728 const auto& curInfo = mITSTrackingInfo[itsClRefs[itsTrack.getClusterEntry(i)]];
729 frameArr[1 + curInfo.lr] = &curInfo;
730 }
731
732 // start refit
733 resTrack.points.clear();
734 resTrack.info.clear();
735 trFit.resetCovariance();
736 trFit.setCov(trFit.getQ2Pt() * trFit.getQ2Pt() * trFit.getCov()[14], 14);
737 float chi2{0};
738 for (int i{7}; i >= 0; --i) {
739 Measurement point;
740 int res = accountCluster(i, trFit, chi2, point, refLin);
741 if (res == 2) {
742 return false;
743 } else if (res == 0) {
744 resTrack.points.push_back(point);
745 resTrack.info.push_back(*frameArr[i]);
746 resTrack.track = trFit; // put track to whatever the IU is
747 }
748 }
749 // reverse inserted points so they are in the same order as the track
750 std::reverse(resTrack.info.begin(), resTrack.info.end());
751 std::reverse(resTrack.points.begin(), resTrack.points.end());
752 resTrack.kfFit.chi2 = chi2;
753 resTrack.kfFit.ndf = (int)resTrack.info.size() * 2 - 5;
754 resTrack.kfFit.chi2Ndf = chi2 / (float)resTrack.kfFit.ndf;
755
756 return true;
757}
758
759void AlignmentSpec::prepareMeasurments(std::span<const itsmft::CompClusterExt> clusters, std::span<const unsigned char> patterns)
760{
761 LOGP(info, "Preparing {} measurments", clusters.size());
762 auto geom = its::GeometryTGeo::Instance();
763 geom->fillMatrixCache(o2::math_utils::bit2Mask(o2::math_utils::TransformType::T2L, o2::math_utils::TransformType::L2G));
764 mITSTrackingInfo.clear();
765 mITSTrackingInfo.reserve(clusters.size());
766 auto pattIt = patterns.begin();
767 for (const auto& cls : clusters) {
768 const auto sens = cls.getSensorID();
769 const auto lay = geom->getLayer(sens);
770 double sigmaY2{0}, sigmaZ2{0};
772 if (mIsITS3) {
773 locXYZ = o2::its3::ioutils::extractClusterData(cls, pattIt, mIT3Dict, sigmaY2, sigmaZ2);
774 } else {
775 locXYZ = o2::its::ioutils::extractClusterData(cls, pattIt, mITSDict, sigmaY2, sigmaZ2);
776 }
777 sigmaY2 += mParams->extraClsErrY[lay] * mParams->extraClsErrY[lay];
778 sigmaZ2 += mParams->extraClsErrZ[lay] * mParams->extraClsErrZ[lay];
779 // Transformation to the local --> global
780 const auto gloXYZ = geom->getMatrixL2G(sens) * locXYZ;
781 // Inverse transformation to the local --> tracking
782 auto trkXYZf = geom->getMatrixT2L(sens) ^ locXYZ;
784 trkXYZ.SetCoordinates(trkXYZf.X(), trkXYZf.Y(), trkXYZf.Z());
785 // Tracking alpha angle
786 // We want that each cluster rotates its tracking frame to the clusters phi
787 // that way the track linearization around the measurement is less biases to the arc
788 // this means automatically that the measurement on the arc is at 0 for the curved layers
789 double alpha = geom->getSensorRefAlpha(sens);
790 double x = trkXYZ.x();
791 if (mIsITS3 && constants::detID::isDetITS3(sens)) {
792 trkXYZ.SetY(0.f);
793 // alpha&x always have to be defined wrt to the global Z axis!
794 x = std::hypot(gloXYZ.x(), gloXYZ.y());
795 trkXYZ.SetX(x);
796 alpha = std::atan2(gloXYZ.y(), gloXYZ.x());
797 auto chip = constants::detID::getSensorID(sens);
798 sigmaY2 += mParams->extraClsErrY[chip] * mParams->extraClsErrY[chip];
799 sigmaZ2 += mParams->extraClsErrZ[chip] * mParams->extraClsErrZ[chip];
800 }
802 mITSTrackingInfo.emplace_back(sens, lay, x, alpha,
803 std::array<double, 2>{trkXYZ.y(), trkXYZ.z()},
804 std::array<double, 3>{sigmaY2, 0., sigmaZ2});
805 }
806}
807
808void AlignmentSpec::buildT2V()
809{
810 const auto& itsTracks = mRecoData->getITSTracks();
811 mT2PV.clear();
812 mT2PV.resize(itsTracks.size(), -1);
813 if (mUseMC) {
814 mPVs.reserve(mcReader->getNEvents(0));
815 for (int iEve{0}; iEve < mcReader->getNEvents(0); ++iEve) {
816 const auto& eve = mcReader->getMCEventHeader(0, iEve);
818 constexpr float err{22e-4f};
819 vtx.setX((float)eve.GetX());
820 vtx.setY((float)eve.GetY());
821 vtx.setZ((float)eve.GetZ());
822 vtx.setSigmaX(err);
823 vtx.setSigmaY(err);
824 vtx.setSigmaZ(err);
825 mPVs.push_back(vtx);
826 }
827 const auto& mcLbls = mRecoData->getITSTracksMCLabels();
828 for (size_t iTrk{0}; iTrk < mcLbls.size(); ++iTrk) {
829 const auto& lbl = mcLbls[iTrk];
830 if (!lbl.isValid() || !lbl.isCorrect()) {
831 continue;
832 }
833 const auto& mcTrk = mcReader->getTrack(lbl);
834 if (mcTrk->isPrimary()) {
835 mT2PV[iTrk] = lbl.getEventID();
836 }
837 }
838 } else {
839 LOGP(fatal, "Data PV to track TODO");
840 }
841}
842
843bool AlignmentSpec::applyMisalignment(Eigen::Vector2d& res, const FrameInfoExt& frame, const TrackD& wTrk, size_t iTrk)
844{
845 if (!constants::detID::isDetITS3(frame.sens)) {
846 return true;
847 }
848
849 const int sensorID = constants::detID::getSensorID(frame.sens);
850 const int layerID = constants::detID::getDetID2Layer(frame.sens);
851
852 // --- Legendre deformation (non-rigid-body) ---
853 if (mParams->doMisalignmentLeg && mIsITS3 && mUseMC) {
854 const auto prop = o2::base::PropagatorD::Instance();
855
856 const auto lbl = mRecoData->getITSTracksMCLabels()[iTrk];
857 const auto mcTrk = mcReader->getTrack(lbl);
858 if (!mcTrk) {
859 return false;
860 }
861 std::array<double, 3> xyz{mcTrk->GetStartVertexCoordinatesX(), mcTrk->GetStartVertexCoordinatesY(), mcTrk->GetStartVertexCoordinatesZ()};
862 std::array<double, 3> pxyz{mcTrk->GetStartVertexMomentumX(), mcTrk->GetStartVertexMomentumY(), mcTrk->GetStartVertexMomentumZ()};
863 TParticlePDG* pPDG = TDatabasePDG::Instance()->GetParticle(mcTrk->GetPdgCode());
864 if (!pPDG) {
865 return false;
866 }
867 o2::track::TrackParD mcPar(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), false);
868
869 const double r = frame.x;
870 const double gloX = r * std::cos(frame.alpha);
871 const double gloY = r * std::sin(frame.alpha);
872 const double gloZ = frame.positionTrackingFrame[1];
873 auto [u, v] = computeUV(gloX, gloY, gloZ, sensorID, constants::radii[layerID]);
874 const double h = mDeformations[sensorID](u, v);
875
876 auto mcAtCl = mcPar;
877 if (!mcAtCl.rotate(frame.alpha) || !prop->PropagateToXBxByBz(mcAtCl, frame.x)) {
878 return false;
879 }
880
881 const double snp = mcAtCl.getSnp();
882 const double tgl = mcAtCl.getTgl();
883 const double csci = 1. / std::sqrt(1. - (snp * snp));
884 const double dydx = snp * csci;
885 const double dzdx = tgl * csci;
886 const double dy = dydx * h;
887 const double dz = dzdx * h;
888
889 const double newGloY = (r * std::sin(frame.alpha)) + (dy * std::cos(frame.alpha));
890 const double newGloX = (r * std::cos(frame.alpha)) - (dy * std::sin(frame.alpha));
891 const double newGloZ = gloZ + dz;
892 auto [uNew, vNew] = computeUV(newGloX, newGloY, newGloZ, sensorID, constants::radii[layerID]);
893 if (std::abs(uNew) > 1. || std::abs(vNew) > 1.) {
894 return false;
895 }
896
897 res[0] += dy;
898 res[1] += dz;
899 }
900
901 // --- Rigid body misalignment ---
902 // Must use the same derivative chain as GBL:
903 // dres/da_parent = dres/da_TRK * J_L2T_tile * J_L2P_tile
904 // The tile is a pseudo-volume; Millepede fits at the halfBarrel (parent) level.
905 if (mParams->doMisalignmentRB) {
906 GlobalLabel lbl(0, frame.sens, true);
907 if (mChip2Hiearchy.find(lbl) == mChip2Hiearchy.end()) {
908 return true; // sensor not in hierarchy, skip
909 }
910 const auto* tileVol = mChip2Hiearchy.at(lbl);
911
912 // derivative in TRK frame (3x6: rows = dy, dz, dsnp)
914
915 // TRK -> tile LOC
916 const double posTrk[3] = {frame.x, 0., 0.};
917 double posLoc[3];
918 tileVol->getT2L().LocalToMaster(posTrk, posLoc);
919 Matrix66 jacL2T;
920 tileVol->computeJacobianL2T(posLoc, jacL2T);
921 der *= jacL2T;
922
923 // tile LOC -> halfBarrel LOC (same chain as GBL hierarchy walk)
924 der *= tileVol->getJL2P();
925
926 // apply: delta_res = der * delta_a_halfBarrel
927 Eigen::Vector3d shift = der * mRigidBodyParams[sensorID];
928 res[0] += shift[0]; // dy
929 res[1] += shift[1]; // dz
930 }
931
932 return true;
933}
934
936{
937 mDBGOut->Close();
938 mDBGOut.reset();
939}
940
942{
944 return;
945 }
946 if (matcher == ConcreteDataMatcher("ITS", "CLUSDICT", 0)) {
947 LOG(info) << "its cluster dictionary updated";
948 mITSDict = (const o2::itsmft::TopologyDictionary*)obj;
949 return;
950 }
951 if (matcher == ConcreteDataMatcher("IT3", "CLUSDICT", 0)) {
952 LOG(info) << "it3 cluster dictionary updated";
953 mIT3Dict = (const o2::its3::TopologyDictionary*)obj;
954 return;
955 }
956}
957
958DataProcessorSpec getAlignmentSpec(GTrackID::mask_t srcTracks, GTrackID::mask_t srcClusters, bool useMC, bool withPV, bool withITS, OutputEnum out)
959{
960 auto dataRequest = std::make_shared<DataRequest>();
961 std::shared_ptr<o2::base::GRPGeomRequest> ggRequest{nullptr};
962 if (!out[OutputOpt::MilleRes]) {
963 dataRequest->requestTracks(srcTracks, useMC);
964 if (!withITS) {
965 dataRequest->requestIT3Clusters(useMC);
966 } else {
967 dataRequest->requestClusters(srcClusters, useMC);
968 }
969 if (withPV && !useMC) {
970 dataRequest->requestPrimaryVertices(useMC);
971 }
972 ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
973 false, // GRPECS=true
974 true, // GRPLHCIF
975 true, // GRPMagField
976 true, // askMatLUT
978 dataRequest->inputs, // inputs
979 true, // askOnce
980 true); // propagatorD
981 } else {
982 dataRequest->inputs.emplace_back("dummy", "GLO", "DUMMY_OUT", 0);
983 ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
984 false, // GRPECS=true
985 false, // GRPLHCIF
986 false, // GRPMagField
987 false, // askMatLUT
989 dataRequest->inputs);
990 }
991
992 Options opts{
993 {"nthreads", VariantType::Int, 1, {"number of threads"}},
994 };
995
996 return DataProcessorSpec{
997 .name = "its3-alignment",
998 .inputs = dataRequest->inputs,
999 .outputs = {},
1000 .algorithm = AlgorithmSpec{adaptFromTask<AlignmentSpec>(dataRequest, ggRequest, srcTracks, useMC, withPV, withITS, out)},
1001 .options = opts};
1002}
1003} // namespace o2::its3::align
Wrapper container for different reconstructed object types.
Definition of the ClusterTopology class.
Definition of the BuildTopologyDictionary class for ITS3.
int32_t i
#define failed(...)
Definition Utils.h:42
Helper for geometry and GRP related CCDB requests.
Definition of the GeometryTGeo class.
Definition of the NDim Legendre Polynominals.
uint32_t j
Definition RawData.h:0
uint32_t res
Definition RawData.h:0
uint32_t c
Definition RawData.h:2
nlohmann::json json
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
Measurement
Class for time synchronization of RawReader instances.
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 printKeyValues(bool showProv=true, bool useLogger=false, bool withPadding=true, bool showHash=true) const final
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
ConfigParamRegistry const & options()
Definition InitContext.h:33
std::map< GlobalLabel, AlignableVolume * > SensorMapping
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
void run(ProcessingContext &pc) final
void init(InitContext &ic) final
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
static GeometryTGeo * Instance()
void fillMatrixCache(int mask) override
int getClusterEntry(int i) const
Definition TrackITS.h:67
std::string pstring(bool withNewline=false) const
Definition EnumFlags.h:510
GLdouble n
Definition glcorearb.h:1982
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
const GLfloat * m
Definition glcorearb.h:4066
GLenum src
Definition glcorearb.h:1767
const GLdouble * v
Definition glcorearb.h:832
GLdouble f
Definition glcorearb.h:310
GLboolean * data
Definition glcorearb.h:298
GLuint GLfloat x0
Definition glcorearb.h:5034
GLboolean r
Definition glcorearb.h:1233
GLuint id
Definition glcorearb.h:650
Node par(int index)
Parameters.
Defining PrimaryVertex explicitly as messageable.
std::vector< ConfigParamSpec > Options
o2::track::TrackParCovD TrackD
Definition TrackFit.h:27
Eigen::Matrix< double, 3, 6 > Matrix36
AlignableVolume::Ptr buildHierarchyITS(AlignableVolume::SensorMapping &sensorMap)
Eigen::Matrix< double, 6, 6 > Matrix66
o2::framework::DataProcessorSpec getAlignmentSpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool useMC, bool withPV, bool withITS3, OutputEnum out)
void writeMillepedeResults(AlignableVolume *root, const std::string &milleResPath, const std::string &outJsonPath, const std::string &injectedJsonPath="")
void applyDOFConfig(AlignableVolume *root, const std::string &jsonPath)
AlignableVolume::Ptr buildHierarchyIT3(AlignableVolume::SensorMapping &sensorMap)
auto getRigidBodyDerivatives(const auto &trk)
bool isDetITS3(T detID)
Definition SpecsV2.h:210
constexpr std::array< double, nLayers > radii
Definition SpecsV2.h:134
constexpr double equatorialGap
Definition SpecsV2.h:131
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::math_utils::Point3D< T > extractClusterData(const itsmft::CompClusterExt &c, iterator &iter, const itsmft::TopologyDictionary *dict, T &sig2y, T &sig2z)
Definition IOUtils.h:41
void bringToPMPid(double &phi)
Definition Utils.h:105
constexpr int kNParams
TrackParametrizationWithError< double > TrackParCovD
Definition Track.h:32
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Defining DataPointCompositeObject explicitly as copiable.
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
o2::base::PropagatorD::MatCorrType corrType
std::vector< Measurement > points
o2::track::TrackParCovD track
std::vector< FrameInfoExt > info
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"
std::vector< Cluster > clusters
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))