CbmRoot
Loading...
Searching...
No Matches
CbmBbaAlignmentTask.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023-2023 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: S.Gorbunov[committer], N.Bluhme */
4
11
12// Cbm Headers ----------------------
13#include "CbmBbaAlignmentTask.h"
14
15
16// ROOT headers
17
18#include "FairRootManager.h"
19#include "FairRunAna.h"
20#include "TClonesArray.h"
21#include "TFile.h"
22#include "TGeoPhysicalNode.h"
23#include "TH1F.h"
24#include "TRandom.h"
25
26//
28#include "CbmGlobalTrack.h"
29#include "CbmL1.h"
30#include "CbmMuchTrack.h"
32#include "CbmMvdHit.h"
34#include "CbmStsHit.h"
35#include "CbmStsSetup.h"
36#include "CbmStsTrack.h"
38#include "CbmTofTrack.h"
40#include "CbmTrdTrack.h"
42#include "TCanvas.h"
43#include "TGeoPhysicalNode.h"
44#include "bba/BBA.h"
45
46#include <cmath>
47#include <iostream>
48#include <thread>
49
50namespace
51{
52 using namespace cbm::algo;
53}
54
55std::pair<std::string, TGeoHMatrix> AlignNode(std::string path, double shiftX, double shiftY, double shiftZ,
56 double rotX, double rotY, double rotZ)
57{
58
59 TGeoHMatrix result;
60 result.SetDx(shiftX);
61 result.SetDy(shiftY);
62 result.SetDz(shiftZ);
63 result.RotateX(rotX);
64 result.RotateY(rotY);
65 result.RotateZ(rotZ);
66
67 //std::cout << "Alignment matrix for node " << path << " is: " << std::endl;
68 //result.Print();
69 //std::cout << std::endl;
70
71 return std::pair<std::string, TGeoHMatrix>(path, result);
72}
73
75{
76 fNmvdHits = 0;
77 fNstsHits = 0;
78 fNmuchHits = 0;
79 fNtrd1dHits = 0;
80 fNtrd2dHits = 0;
81 fNtofHits = 0;
82 for (const auto& n : fUnalignedTrack.fNodes) {
83 switch (n.fHitSystemId) {
84 case ECbmModuleId::kMvd: fNmvdHits++; break;
85 case ECbmModuleId::kSts: fNstsHits++; break;
86 case ECbmModuleId::kMuch: fNmuchHits++; break;
87 case ECbmModuleId::kTrd: fNtrd1dHits++; break;
88 case ECbmModuleId::kTrd2d: fNtrd2dHits++; break;
89 case ECbmModuleId::kTof: fNtofHits++; break;
90 default: break;
91 }
92 }
93}
94
95
96CbmBbaAlignmentTask::CbmBbaAlignmentTask(const char* name, Int_t iVerbose, TString histoFileName)
97 : FairTask(name, iVerbose)
98 , fHistoFileName(histoFileName)
99{
100 TFile* curFile = gFile;
101 TDirectory* curDirectory = gDirectory;
102
103 if (!(fHistoFileName == "")) {
104 fHistoFile = new TFile(fHistoFileName.Data(), "RECREATE");
105 }
106 else {
107 fHistoFile = gFile;
108 }
109
110 fHistoFile->cd();
111
112 fHistoDir = fHistoFile->mkdir("Bba");
113 fHistoDir->cd();
114
115 gFile = curFile;
116 gDirectory = curDirectory;
117}
118
119
121
122
124{
125
126 {
127 const char* env = std::getenv("OMP_NUM_THREADS");
128 if (env) {
129 fNthreads = TString(env).Atoi();
130 LOG(info) << "BBA: found environment variable OMP_NUM_THREADS = \"" << env
131 << "\", read as integer: " << fNthreads;
132 }
133 }
134
135 // ensure that at least one thread is set
136
137 if (fNthreads < 1) {
138 fNthreads = 1;
139 }
140
141 // fNthreads = 1; // enforce n threads to one
142
143 fFitter.Init();
145
146 //Get ROOT Manager
147 FairRootManager* ioman = FairRootManager::Instance();
148
149 if (!ioman) {
150 LOG(error) << "CbmBbaAlignmentTask::Init :: RootManager not instantiated!";
151 return kERROR;
152 }
153
154 // Get global tracks & matches
155
157
158 fInputGlobalTracks = static_cast<TClonesArray*>(ioman->GetObject("GlobalTrack"));
159 if (!fInputGlobalTracks) {
160 LOG(error) << "CbmBbaAlignmentTask::Init: Global track array not found!";
161 return kERROR;
162 }
163
164 fInputGlobalTrackMatches = static_cast<TClonesArray*>(ioman->GetObject("GlobalTrackMatch"));
165
167 LOG(error) << "CbmBbaAlignmentTask::Init: Global track matches array not found!";
168 //return kERROR;
169 }
170 }
171
172 fInputStsTracks = static_cast<TClonesArray*>(ioman->GetObject("StsTrack"));
173 if (!fInputStsTracks) {
174 LOG(error) << "CbmBbaAlignmentTask::Init: Sts track array not found!";
175 return kERROR;
176 }
177
178 fInputStsTrackMatches = static_cast<TClonesArray*>(ioman->GetObject("StsTrackMatch"));
179
181 LOG(error) << "CbmBbaAlignmentTask::Init: Sts track matches array not found!";
182 //return kERROR;
183 }
184
185 // MC tracks
186 fInputMcTracks = static_cast<TClonesArray*>(ioman->GetObject("MCTrack"));
187 if (!fInputMcTracks) {
188 Warning("CbmBbaAlignmentTask::Init", "mc track array not found!");
189 //return kERROR;
190 }
191
192 fTracks.clear();
193
194
196 const ca::Parameters<ca::fvec>& l1Param = l1->GetParameters();
198 // TO get rid of weird station counting:
199 // fNtrackingStations = l1Param.GetNstationsActive() - l1Param. GetNstationsActive(ca::EDetectorID::kMvd);
200
201 TDirectory* curDirectory = gDirectory;
202
203 fHistoDir->cd();
204
205 for (int i = -1; i < fNtrackingStations; i++) {
206 const char* n1 = Form("Station%d", i);
207 const char* n2 = Form(", station %d", i);
208 if (i == -1) {
209 n1 = "All";
210 n2 = "";
211 }
212 fHistoDir->mkdir(n1);
213 fHistoDir->cd(n1);
214
215 int nBins = 301;
216
217 double sizeX = 0.1;
218 double sizeY = 0.1;
219 double sizePX = 10.;
220 double sizePY = 10.;
221 if (fTrackingMode == kSts) {
222 sizeX = 0.1;
223 sizeY = 0.1;
224 }
225 else if (fTrackingMode == kMcbm) {
226
227 if (i == -1) {
228 sizeX = 5.;
229 sizeY = 5.;
230 sizePX = 10.;
231 sizePY = 10.;
232 }
233
234 if (i == 0 || i == 1) {
235 sizeX = 1.0;
236 sizeY = 1.0;
237 }
238 if (i == 2) {
239 sizeX = 5.;
240 sizeY = 5.;
241 }
242 if (i == 3) {
243 sizeX = 5.;
244 sizeY = 20.;
245 sizePX = 10.;
246 }
247 if (i == 4) {
248 sizeX = 20.;
249 sizeY = 5.;
250 sizePY = 10.;
251 }
252 if (i >= 5) {
253 sizeX = 5.;
254 sizeY = 5.;
255 }
256 }
257
259 new TH1F(Form("ResBeforeAli%s_X", n1), Form("X-Residuals Before Alignment%s", n2), nBins, -sizeX, sizeX));
261 new TH1F(Form("ResBeforeAli%s_Y", n1), Form("Y-Residuals Before Alignment%s", n2), nBins, -sizeY, sizeY));
263 new TH1F(Form("ResAfterAli%s_X", n1), Form("X-Residuals After Alignment%s", n2), nBins, -sizeX, sizeX));
265 new TH1F(Form("ResAfterAli%s_Y", n1), Form("Y-Residuals After Alignment%s", n2), nBins, -sizeY, sizeY));
266
267 hPullsBeforeAlignmentX.push_back(
268 new TH1F(Form("PullBeforeAli%s_X", n1), Form("X-Pulls Before Alignment%s", n2), nBins, -sizePX, sizePX));
269 hPullsBeforeAlignmentY.push_back(
270 new TH1F(Form("PullBeforeAli%s_Y", n1), Form("Y-Pulls Before Alignment%s", n2), nBins, -sizePY, sizePY));
271 hPullsAfterAlignmentX.push_back(
272 new TH1F(Form("PullAfterAli%s_X", n1), Form("X-Pulls After Alignment%s", n2), nBins, -sizePX, sizePX));
273 hPullsAfterAlignmentY.push_back(
274 new TH1F(Form("PullAfterAli%s_Y", n1), Form("Y-Pulls After Alignment%s", n2), nBins, -sizePY, sizePY));
275 }
276
277 for (int i = 0; i < fNtrackingStations + 1; i++) {
278 hResidualsAfterAlignmentX[i]->SetLineColor(kRed);
279 hResidualsAfterAlignmentY[i]->SetLineColor(kRed);
280 hPullsAfterAlignmentX[i]->SetLineColor(kRed);
281 hPullsAfterAlignmentY[i]->SetLineColor(kRed);
282 }
283
284 gDirectory = curDirectory;
285
286 return kSUCCESS;
287}
288
289void CbmBbaAlignmentTask::Exec(Option_t* /*opt*/)
290{
291
292 LOG(info) << "BBA: exec event N " << fNEvents;
293
294 fNEvents++;
295
296 if (static_cast<int>(fTracks.size()) >= fMaxNtracks) {
297 return;
298 }
299
301 const ca::Parameters<ca::fvec>& l1Par = l1->GetParameters();
302
303 fFitter.SetDoSmooth(true);
304
305 // select tracks for alignment and store them
306 unsigned num_rejected_fit = 0;
308
310
311 for (int iTr = 0; iTr < fInputStsTracks->GetEntriesFast(); iTr++) {
312 if (static_cast<int>(fTracks.size()) >= fMaxNtracks) {
313 break;
314 }
315 const CbmStsTrack* stsTrack = dynamic_cast<CbmStsTrack*>(fInputStsTracks->At(iTr));
316 if (!stsTrack) continue;
318 if (!fFitter.CreateMvdStsTrack(t.fUnalignedTrack, iTr)) continue;
319 t.MakeConsistent();
320
321 for (auto& n : t.fUnalignedTrack.fNodes) {
322 n.fIsTimeSet = false;
323 // attempt to switch off MVD
324 // if (n.fHitSystemId == ECbmModuleId::kMvd) {
325 // n.fIsXySet = false;
326 // }
327 }
328
330 LOG(debug) << "failed to fit the track! ";
331 num_rejected_fit++;
332 continue;
333 }
334
335 const auto& par = t.fUnalignedTrack.fNodes.front().fParamUp;
336
337 if (t.fNstsHits < l1Par.GetNstationsActive(ca::EDetectorID::kSts)) continue;
338 if (!std::isfinite(par.GetQp())) continue;
339 if (fabs(par.GetQp()) > 1.) continue; // select tracks with min 1 GeV momentum
340
342 fTracks.push_back(t);
343 }
344 }
345 else if (fTrackingMode == kMcbm && fInputGlobalTracks) {
346
349
350 for (int iTr = 0; iTr < fInputGlobalTracks->GetEntriesFast(); iTr++) {
351 if (static_cast<int>(fTracks.size()) >= fMaxNtracks) {
352 break;
353 }
354 const CbmGlobalTrack* globalTrack = dynamic_cast<const CbmGlobalTrack*>(fInputGlobalTracks->At(iTr));
355 if (!globalTrack) {
356 LOG(fatal) << "BBA: null pointer to the global track!";
357 break;
358 }
360 if (!fFitter.CreateGlobalTrack(t.fUnalignedTrack, *globalTrack)) {
361 LOG(fatal) << "BBA: can not create the global track for the fit! ";
362 break;
363 }
364 t.MakeConsistent();
365
366 for (auto& n : t.fUnalignedTrack.fNodes) {
367 n.fIsTimeSet = false;
368 if (n.fHitSystemId == ECbmModuleId::kTrd) {
369 if (n.fMxy.Dx2() < n.fMxy.Dy2()) {
370 n.fMxy.SetNdfY(0);
371 }
372 else {
373 n.fMxy.SetNdfX(0);
374 }
375 }
376 }
377
378 if (t.fNstsHits < 2) continue;
379 if (t.fNtofHits < 2) continue;
380 if (t.fNtrd1dHits + t.fNtrd2dHits < 2) continue;
382 LOG(debug) << "failed to fit the track! ";
383 num_rejected_fit++;
384 continue;
385 }
387 fTracks.push_back(t);
388 }
389 } // end of mcbm tracking mode
390
391 if (num_rejected_fit != 0) {
392 LOG(warn) << "failed to fit " << num_rejected_fit << " tracks";
393 }
394
395 // fix the material radiation thickness to avoid the chi2 rattling at the material map bin boundaries
396 // (to be improved)
397 unsigned num_fix_rad_max = 0;
398 unsigned num_fix_rad_min = 0;
399 unsigned num_rejected_fit2 = 0;
400 for (TrackContainer& t : fTracks) {
401 if (!t.fIsActive) continue;
402 for (unsigned int in = 0; in < t.fUnalignedTrack.fNodes.size(); in++) {
403 CbmKfTrackFitter::TrajectoryNode& node = t.fUnalignedTrack.fNodes[in];
404 if (!node.fIsFitted) {
405 t.fIsActive = false;
406 num_rejected_fit2++;
407 break;
408 }
409 node.fIsRadThickFixed = true;
410 if (node.fRadThick > 0.01) {
411 LOG(debug) << "CbmBbaAlignmentTask: radiation thickness is too large: " << node.fRadThick;
412 num_fix_rad_max++;
413 node.fRadThick = 0.01;
414 }
415 if (node.fRadThick < 0.001) {
416 LOG(debug) << "CbmBbaAlignmentTask: radiation thickness is too small: " << node.fRadThick;
417 num_fix_rad_min++;
418 node.fRadThick = 0.001;
419 }
420 }
421 }
422 if (num_fix_rad_max != 0) {
423 LOG(warn) << "CbmBbaAlignmentTask: radiation thickness is too large " << num_fix_rad_max << " times";
424 }
425 if (num_fix_rad_min != 0) {
426 LOG(warn) << "CbmBbaAlignmentTask: radiation thickness is too small " << num_fix_rad_min << " times";
427 }
428 if (num_rejected_fit2 != 0) {
429 LOG(warn) << "Rejected fit 2 " << num_rejected_fit2 << " tracks";
430 }
431
432 // ensure that all the hits are not too much deviated from the trajectory
433 // (to be improved)
434 unsigned num_rejected_traj = 0;
435 for (TrackContainer& t : fTracks) {
436 if (!t.fIsActive) continue;
437 const double cutX = 8.;
438 const double cutY = 8.;
439 for (unsigned int in = 0; in < t.fUnalignedTrack.fNodes.size(); in++) {
440 CbmKfTrackFitter::Trajectory tr = t.fUnalignedTrack;
442 if (!node.fIsXySet) {
443 continue;
444 }
445 node.fIsXySet = false; // exclude the node from the fit
446 if (!fFitter.FitTrajectory(tr)) {
447 t.fIsActive = false;
448 num_rejected_traj++;
449 break;
450 }
451 double x_residual = node.fMxy.X() - node.fParamDn.X(); // this should be the difference vector
452 double y_residual = node.fMxy.Y() - node.fParamDn.Y();
453 double x_pull = x_residual / sqrt(node.fMxy.Dx2() + node.fParamDn.GetCovariance(0, 0));
454 double y_pull = y_residual / sqrt(node.fMxy.Dy2() + node.fParamDn.GetCovariance(1, 1));
455 if (!node.fIsFitted || (node.fMxy.NdfX() > 0. && fabs(x_pull) > cutX)
456 || (node.fMxy.NdfY() > 0. && fabs(y_pull) > cutY)) {
457 t.fIsActive = false;
458 num_rejected_traj++;
459 break;
460 }
461 }
462 }
463 if (num_rejected_traj != 0) {
464 LOG(warn) << "Rejected " << num_rejected_traj << " tracks for hits deviating from the trajectory";
465 }
466
467 TClonesArray* inputTracks = nullptr;
468 if (fTrackingMode == kSts) {
469 inputTracks = fInputStsTracks;
470 }
471 else if (fTrackingMode == kMcbm) {
472 inputTracks = fInputGlobalTracks;
473 }
474 static int statTracks = 0;
475 statTracks += inputTracks->GetEntriesFast();
476
477 unsigned active_tracks = 0;
478 for (auto& t : fTracks) {
479 if (t.fIsActive) active_tracks++;
480 }
481
482 LOG(info) << "Bba: " << inputTracks->GetEntriesFast() << " tracks in event, " << statTracks << " tracks in total, "
483 << fTracks.size() << " tracks selected for alignment, " << active_tracks << " tracks active";
484}
485
486
487void CbmBbaAlignmentTask::ConstrainStation(std::vector<double>& par, int iSta, int ixyz)
488{
489 // set overall shifts of the station to zero
490
491 double mean = 0.;
492 int n = 0;
493 for (unsigned int i = 0; i < fAlignmentBodies.size(); i++) {
494 if (fAlignmentBodies[i].fTrackingStation == iSta) {
495 mean += par[3 * i + ixyz];
496 n++;
497 }
498 }
499 if (n <= 0) return;
500 mean /= n;
501 for (unsigned int i = 0; i < fAlignmentBodies.size(); i++) {
502 if (fAlignmentBodies[i].fTrackingStation == iSta) {
503 par[3 * i + ixyz] -= mean;
504 }
505 }
506}
507
508void CbmBbaAlignmentTask::ApplyConstraints(std::vector<double>& par)
509{
510 // apply alignment parameters to the hits
511 ConstrainStation(par, 0, 2);
515}
516
517
518void CbmBbaAlignmentTask::ApplyAlignment(const std::vector<double>& par)
519{
520 // apply alignment parameters to the hits
521
522 std::vector<double> parConstrained = par;
523
524 ApplyConstraints(parConstrained);
525
526 for (TrackContainer& t : fTracks) {
527
528 for (unsigned int in = 0; in < t.fUnalignedTrack.fNodes.size(); in++) {
529 CbmKfTrackFitter::TrajectoryNode& node = t.fUnalignedTrack.fNodes[in];
530 if (!node.fIsXySet) {
531 continue;
532 }
533 CbmKfTrackFitter::TrajectoryNode& nodeAligned = t.fAlignedTrack.fNodes[in];
534 int iSensor = node.fReference1;
535 assert(iSensor >= 0 && iSensor < (int) fSensors.size());
536 assert(nodeAligned.fReference1 == node.fReference1);
537
538 Sensor& s = fSensors[iSensor];
539 if (s.fAlignmentBody < 0) continue;
540
541 double dx = parConstrained[3 * s.fAlignmentBody + 0];
542 double dy = parConstrained[3 * s.fAlignmentBody + 1];
543 double dz = parConstrained[3 * s.fAlignmentBody + 2];
544
545 // shift the hit
546 nodeAligned.fMxy.SetX(node.fMxy.X() + dx);
547 nodeAligned.fMxy.SetY(node.fMxy.Y() + dy);
548 nodeAligned.fZ = node.fZ + dz;
549
550 // shift the fitted track to the Z position of the hit
551 {
552 auto& p = nodeAligned.fParamDn;
553 p.SetX(p.X() + p.Tx() * dz);
554 p.SetY(p.Y() + p.Ty() * dz);
555 p.SetZ(nodeAligned.fZ);
556 }
557 {
558 auto& p = nodeAligned.fParamUp;
559 p.SetX(p.X() + p.Tx() * dz);
560 p.SetY(p.Y() + p.Ty() * dz);
561 p.SetZ(nodeAligned.fZ);
562 }
563 }
564 }
565
566 static int statNCalls = 0;
567 statNCalls++;
568 if (statNCalls % 100 == 0) {
569 std::vector<double>& r = parConstrained;
570 LOG(info) << "Current parameters: ";
571 for (int is = 0; is < fNalignmentBodies; is++) {
572 auto& b = fAlignmentBodies[is];
573 LOG(info) << "Body " << is << " sta " << b.fTrackingStation << ": x " << r[3 * is + 0] << " y " << r[3 * is + 1]
574 << " z " << r[3 * is + 2];
575 }
576 }
577}
578
579
580double CbmBbaAlignmentTask::CostFunction(const std::vector<double>& par)
581{
582 // apply new parameters to the hits
583
584 for (auto& t : fTracks) {
585 t.fAlignedTrack = t.fUnalignedTrack;
586 }
587
588 ApplyAlignment(par);
589
590 std::vector<double> chi2Thread(fNthreads, 0.);
591 std::vector<long> ndfThread(fNthreads, 0);
592
593 auto fitThread = [&](int iThread) {
595 fit.SetDoSmooth(false);
597
598 for (unsigned int itr = iThread; itr < fTracks.size(); itr += fNthreads) {
599 auto& t = fTracks[itr];
600 if (!t.fIsActive) continue;
601
602 // fit.SetDebugInfo(" track " + std::to_string(itr));
603
604 for (unsigned int in = 0; in < t.fAlignedTrack.fNodes.size(); in++) {
605 CbmKfTrackFitter::TrajectoryNode& node = t.fAlignedTrack.fNodes[in];
606 if (!node.fIsFitted) {
607 LOG(fatal) << "BBA: Cost function: aligned node is not fitted! ";
608 }
609 }
610 bool ok = fit.FitTrajectory(t.fAlignedTrack);
611 double chi2 = t.fAlignedTrack.fNodes.back().fParamDn.GetChiSq();
612 double ndf = t.fAlignedTrack.fNodes.back().fParamDn.GetNdf();
613 if (!ok || !std::isfinite(chi2) || chi2 <= 0. || ndf <= 0.) {
614 // TODO: deactivate the track and continue
615 LOG(fatal) << "BBA: fit failed! ";
616 chi2 = 0.;
617 ndf = 0.;
618 }
619 chi2Thread[iThread] += chi2;
620 ndfThread[iThread] += (int) ndf;
621 }
622 };
623
624 std::vector<std::thread> threads(fNthreads);
625
626 // run n threads
627 for (int i = 0; i < fNthreads; i++) {
628 threads[i] = std::thread(fitThread, i);
629 }
630
631 // wait for the threads to finish
632 for (auto& th : threads) {
633 th.join();
634 }
635
636 fChi2Total = 0.;
637 fNdfTotal = 0;
638 for (int i = 0; i < fNthreads; i++) {
639 fChi2Total += chi2Thread[i];
640 fNdfTotal += ndfThread[i];
641 }
642
643 double cost = fChi2Total / (fNdfTotal + 1);
644 if (fFixedNdf > 0 && fNdfTotal != fFixedNdf) {
645 cost = 1.e30;
646 }
647 //LOG(info) << "BBA: calculate cost function: ndf " << fNdfTotal << ", cost " << cost
648 // << ", diff to ideal cost: " << cost - fCostIdeal << ", diff to initial cost: " << cost - fCostInitial;
649
650 return cost;
651}
652
654{
655 // perform the alignment
656
657 LOG(info) << "BBA: start the alignment procedure with " << fTracks.size() << " tracks ...";
658
659 // collect the sensors from the hits
660 // TODO: read it from the detector setup
661
662 if (fTrackingMode == kSts) {
664 }
665
666 std::set<Sensor> sensorSet;
667 for (auto& t : fTracks) {
668 for (auto& n : t.fUnalignedTrack.fNodes) {
669 if (!n.fIsXySet) {
670 continue;
671 }
672 Sensor s;
673 s.fSystemId = n.fHitSystemId;
674 s.fSensorId = n.fHitAddress;
675 // TODO: get the station index from n.fHitSystemId, n.fHitAddress
676 s.fTrackingStation = n.fMaterialLayer;
677
680 }
681
684 // hardcode path to physical node in geometry
685 if (s.fTrackingStation == 2) {
686 s.fNodePath = "/cave_1/trd_v22h_mcbm_0/layer01_20101/module9_101001001";
687 }
688 else if (s.fTrackingStation == 3) {
689 s.fNodePath = "/cave_1/trd_v22h_mcbm_0/layer02_10202/module8_101002001";
690 }
691 else if (s.fTrackingStation == 4) {
692 s.fNodePath = "/cave_1/trd_v22h_mcbm_0/layer03_11303/module8_101303001";
693 }
694 }
695 else if (s.fSystemId == ECbmModuleId::kTof) {
696 s.fSensorId = CbmTofAddress::GetRpcFullId(n.fHitAddress);
697 }
698
699 sensorSet.insert(s);
700 }
701 }
702
703 fSensors.clear();
704 for (auto& s : sensorSet) { // convert set to a vector
705 fSensors.push_back(s);
706 }
707 std::sort(fSensors.begin(), fSensors.end(), [](const Sensor& a, const Sensor& b) { return a < b; });
708
709 for (auto& t : fTracks) {
710 for (auto& n : t.fUnalignedTrack.fNodes) {
711 if (!n.fIsXySet) {
712 continue;
713 }
714 Sensor s;
715 s.fSystemId = n.fHitSystemId;
716 s.fSensorId = n.fHitAddress;
717 s.fTrackingStation = n.fMaterialLayer;
718
721 }
722
725 }
726 else if (s.fSystemId == ECbmModuleId::kTof) {
727 s.fSensorId = CbmTofAddress::GetRpcFullId(n.fHitAddress);
728 }
729 auto iter = std::lower_bound( // find the sensor in the vector
730 fSensors.begin(), fSensors.end(), s, [](const Sensor& a, const Sensor& b) { return a < b; });
731 assert(iter != fSensors.end());
732 assert(*iter == s);
733 int iSensor = std::distance(fSensors.begin(), iter);
734 n.fReference1 = iSensor;
735 }
736 t.fAlignedTrack = t.fUnalignedTrack;
737 }
738
739 if (1) { // one alignment body per tracking station
740
743 for (int i = 0; i < fNalignmentBodies; i++) {
744 fAlignmentBodies[i].fTrackingStation = i;
745 }
746 for (auto& s : fSensors) {
747 s.fAlignmentBody = s.fTrackingStation;
748 }
749 }
750 else { // one alignment body per sensor
751
754 for (int unsigned i = 0; i < fSensors.size(); i++) {
755 fSensors[i].fAlignmentBody = i;
756 fAlignmentBodies[i].fTrackingStation = fSensors[i].fTrackingStation;
757 }
758 }
759
760 LOG(info) << "BBA: " << fSensors.size() << " sensors, " << fNalignmentBodies << " alignment bodies";
761
762 LOG(info) << "\n Sensors: ";
763 {
764 for (unsigned int is = 0; is < fSensors.size(); is++) {
765 auto& s = fSensors[is];
766 LOG(info) << "Sensor " << is << " sys " << s.fSystemId << " sta " << s.fTrackingStation << " body "
767 << s.fAlignmentBody;
768 }
769 }
770
771 LOG(info) << "\n Alignment bodies: ";
772 {
773 for (int is = 0; is < fNalignmentBodies; is++) {
774 auto& b = fAlignmentBodies[is];
775 LOG(info) << "Body: " << is << " station: " << b.fTrackingStation;
776 }
777 }
778
779 int nParameters = 3 * fNalignmentBodies; // x, y, z
780
781 std::vector<bba::Parameter> par(nParameters);
782
783 for (int ip = 0; ip < nParameters; ip++) {
784 bba::Parameter& p = par[ip];
785 p.SetActive(0);
786 p.SetValue(0.);
787 p.SetBoundaryMin(-20.); // +-20 cm alignment range
788 p.SetBoundaryMax(20.);
789 p.SetStepMin(1.e-4);
790 p.SetStepMax(0.1);
791 p.SetStep(1.e-2);
792 }
793
794 // set active parameters
795
796 for (unsigned int ib = 0; ib < fAlignmentBodies.size(); ib++) {
797 auto& b = fAlignmentBodies[ib];
798 if (b.fTrackingStation == 0) { // the first station
799 //continue;
800 }
801 if (b.fTrackingStation == fNtrackingStations - 1) { // the last station
802 //continue;
803 }
804 if (fTrackingMode == kSts && b.fTrackingStation == fNtrackingStations / 2) { // station in the middle for STS mode
805 par[3 * ib + 0].SetActive(0);
806 par[3 * ib + 1].SetActive(1);
807 par[3 * ib + 2].SetActive(0);
808 }
809 else {
810 par[3 * ib + 0].SetActive(1);
811 par[3 * ib + 1].SetActive(1);
812 par[3 * ib + 2].SetActive(1);
813 }
814 }
815
816 // set parameter ranges for different subsystems
817 for (const auto& s : fSensors) {
818 int ib = s.fAlignmentBody;
819 if (ib < 0 || ib >= fNalignmentBodies) continue;
820 //auto& b = fAlignmentBodies[ib];
821 for (int j = 0; j < 3; j++) {
822 bba::Parameter& p = par[ib * 3 + j];
823 p.SetStepMin(1.e-4);
824 if (s.fSystemId == ECbmModuleId::kTrd || s.fSystemId == ECbmModuleId::kTrd2d) {
825 p.SetStepMin(10.e-4);
826 }
827 else if (s.fSystemId == ECbmModuleId::kTof) {
828 p.SetStepMin(10.e-4);
829 }
830 }
831 }
832
833 gRandom->SetSeed(1);
834
835 for (int is = 0; is < fNalignmentBodies; is++) {
836 bba::Parameter& px = par[3 * is + 0];
837 bba::Parameter& py = par[3 * is + 1];
838 bba::Parameter& pz = par[3 * is + 2];
839
840 // +- d cm misalignment
842 if (px.IsActive()) {
843 px.SetValue(gRandom->Uniform(2. * d) - d);
844 }
845 if (py.IsActive()) {
846 py.SetValue(gRandom->Uniform(2. * d) - d);
847 }
848 if (pz.IsActive()) {
849 pz.SetValue(gRandom->Uniform(2. * d) - d);
850 }
851 }
852
853 if (0) { // test
854 LOG(info) << "STS sensor paths: ";
855 for (auto& s : fSensors) {
856 if (s.fSystemId != ECbmModuleId::kSts) {
857 continue;
858 }
859 const CbmStsElement* el = CbmStsSetup::Instance()->GetElement(s.fSensorId, kStsSensor);
860 if (!el) {
861 LOG(error) << "Element not found: " << s.fSensorId;
862 continue;
863 }
864 el->Print();
865 const auto* n = el->GetPnode();
866 if (n) {
867 LOG(info) << "Node: ";
868 n->Print();
869 }
870 }
871 }
872
873 bba::BBA alignment;
874
875 alignment.setParameters(par);
876
877 auto lambdaCostFunction = [this](const std::vector<double>& p) { return CostFunction(p); };
878
879 alignment.setChi2(lambdaCostFunction);
880 alignment.printInfo();
881
882 std::vector<double> parIdeal;
883 std::vector<double> parMisaligned;
884 {
885 const std::vector<double>& r = alignment.getResult();
886 for (unsigned int i = 0; i < r.size(); i++) {
887 parMisaligned.push_back(r[i]);
888 parIdeal.push_back(0.);
889 }
890 }
891
892 {
893 std::cout << "initial cost function..." << std::endl;
894
895 fCostInitial = CostFunction(parMisaligned);
896
897 unsigned tracks_rejected = 0;
898 for (auto& t : fTracks) {
899 double ndf = t.fAlignedTrack.fNodes.back().fParamDn.GetNdf();
900 double chi2 = t.fAlignedTrack.fNodes.back().fParamDn.GetChiSq();
901 if (ndf < 0. || chi2 < 0. || !std::isfinite(ndf) || !std::isfinite(chi2)) {
902 t.fIsActive = false;
903 tracks_rejected++;
904 }
905 }
906 std::cout << "reject " << tracks_rejected << " bad tracks and recalculate the initial cost function..."
907 << std::endl;
908 fCostInitial = CostFunction(parMisaligned);
909 fFixedNdf = -1; //fNdfTotal;
910
911 // plot the residuals before alignment
912
913 for (const auto& t : fTracks) {
914 if (!t.fIsActive) continue;
915 // calculate the residuals for all tracks at each TrajectoryNode before alignment
916 for (unsigned int in = 0; in < t.fAlignedTrack.fNodes.size(); in++) {
917 CbmKfTrackFitter::Trajectory tr = t.fAlignedTrack;
919
920 if (!node.fIsXySet) {
921 continue;
922 }
923 node.fIsXySet = false;
924
926
927 if (!node.fIsFitted) {
928 continue;
929 }
930
931 Sensor& s = fSensors[node.fReference1];
932
933 double x_residual = node.fMxy.X() - node.fParamDn.X(); // this should be the difference vector
934 double y_residual = node.fMxy.Y() - node.fParamDn.Y();
935 double x_pull = x_residual / sqrt(node.fMxy.Dx2() + node.fParamDn.GetCovariance(0, 0));
936 double y_pull = y_residual / sqrt(node.fMxy.Dy2() + node.fParamDn.GetCovariance(1, 1));
937
938 if (node.fMxy.NdfX() > 0) {
939 hResidualsBeforeAlignmentX[0]->Fill(x_residual);
940 hPullsBeforeAlignmentX[0]->Fill(x_pull);
941 }
942 if (node.fMxy.NdfY() > 0) {
943 hResidualsBeforeAlignmentY[0]->Fill(y_residual);
944 hPullsBeforeAlignmentY[0]->Fill(y_pull);
945 }
946
947 hResidualsBeforeAlignmentX[s.fTrackingStation + 1]->Fill(x_residual);
948 hResidualsBeforeAlignmentY[s.fTrackingStation + 1]->Fill(y_residual);
949
950 hPullsBeforeAlignmentX[s.fTrackingStation + 1]->Fill(x_pull);
951 hPullsBeforeAlignmentY[s.fTrackingStation + 1]->Fill(y_pull);
952 }
953 }
954 }
955
956 std::cout << "calculate ideal cost function..." << std::endl;
957 fCostIdeal = CostFunction(parIdeal);
958
959
960 LOG(info) << " Cost function for the true parameters: " << fCostIdeal;
961 LOG(info) << " Cost function for the initial parameters: " << fCostInitial;
962
963 // run the alignment
964 alignment.align();
965
966 double costResult = CostFunction(alignment.getResult());
967
968 CostFunction(alignment.getResult());
969 CostFunction(alignment.getResult());
970
971 // Create alignment matrices:
972 std::vector<double> result = alignment.getResult();
973 std::map<std::string, TGeoHMatrix> alignmentMatrices;
974
975 for (auto& s : fSensors) {
976 int i = s.fAlignmentBody;
977 assert(i < fNalignmentBodies);
978 if (i < 0) continue;
979
980 // For STS sensors
981 if (s.fSystemId == ECbmModuleId::kSts) {
982 // get sensor for sensorID
983 const CbmStsElement* el = CbmStsSetup::Instance()->GetElement(s.fSensorId, kStsSensor);
984 if (!el) {
985 LOG(error) << "Element not found: " << s.fSensorId;
986 continue;
987 }
988 const auto* n = el->GetPnode();
989 LOG(info) << "Node: ";
990 LOG(info) << n->GetName();
991 alignmentMatrices.insert(
992 AlignNode(n->GetName(), result[i * 3 + 0], result[i * 3 + 1], result[i * 3 + 2], 0., 0., 0.));
993 }
994 // for Trd stations
995 else if (s.fSystemId == ECbmModuleId::kTrd) {
996 LOG(info) << "Node: ";
997 LOG(info) << s.fNodePath;
998 alignmentMatrices.insert(
999 AlignNode(s.fNodePath, result[i * 3 + 0], result[i * 3 + 1], result[i * 3 + 2], 0., 0., 0.));
1000 }
1001 // for Trd2D station
1002 else if (s.fSystemId == ECbmModuleId::kTrd2d) {
1003 LOG(info) << "Node: ";
1004 LOG(info) << s.fNodePath;
1005 alignmentMatrices.insert(
1006 AlignNode(s.fNodePath, result[i * 3 + 0], result[i * 3 + 1], result[i * 3 + 2], 0., 0., 0.));
1007 }
1008 // TODO: for TOF stations
1009 // else {
1010 // LOG(info) << "fSystemId: " << s.fSystemId << "\tfTrackongStation: " << s.fTrackingStation
1011 // << "\t SensorId: " << s.fSensorId;
1012 // alignmentMatrices.insert(AlignNode(Form("fTrackingStation %i", s.fTrackingStation), result[i * 3 + 0],
1013 // result[i * 3 + 1], result[i * 3 + 2], 0., 0., 0.));
1014 // }
1015 } // sensors
1016
1017 // save matrices to disk
1018 TFile* misalignmentMatrixRootfile =
1019 new TFile("AlignmentMatrices_mcbm_beam_2022_05_23_nickel_finetuning.root", "RECREATE");
1020 if (misalignmentMatrixRootfile->IsOpen()) {
1021 gDirectory->WriteObject(&alignmentMatrices, "AlignmentMatrices");
1022 misalignmentMatrixRootfile->Write();
1023 misalignmentMatrixRootfile->Close();
1024 }
1025
1026
1027 LOG(info) << " Cost function for the true parameters: " << fCostIdeal;
1028 LOG(info) << " Cost function for the initial parameters: " << fCostInitial;
1029 LOG(info) << " Cost function for the aligned parameters: " << costResult;
1030
1031 LOG(info) << "\n Misaligned parameters: ";
1032 for (int is = 0; is < fNalignmentBodies; is++) {
1033 const std::vector<double>& r = parMisaligned;
1034 LOG(info) << "Body " << is << ": x " << r[3 * is + 0] << " y " << r[3 * is + 1] << " z " << r[3 * is + 2];
1035 }
1036
1037 LOG(info) << "\n Alignment results: ";
1038 {
1039 const std::vector<double>& r = alignment.getResult();
1040 for (int is = 0; is < fNalignmentBodies; is++) {
1041 auto& b = fAlignmentBodies[is];
1042 LOG(info) << "Body " << is << " sta " << b.fTrackingStation << ": x " << r[3 * is + 0] << " y " << r[3 * is + 1]
1043 << " z " << r[3 * is + 2];
1044 }
1045 }
1046
1047 LOG(info) << "\n";
1048
1049 LOG(info) << " Cost function for the true parameters: " << fCostIdeal;
1050 LOG(info) << " Cost function for the initial parameters: " << fCostInitial;
1051 LOG(info) << " Cost function for the aligned parameters: " << costResult;
1052
1053 LOG(info) << "\n Alignment results, constrained: ";
1054 {
1055 std::vector<double> r = alignment.getResult();
1057 for (int is = 0; is < fNalignmentBodies; is++) {
1058 auto& b = fAlignmentBodies[is];
1059 LOG(info) << "Body " << is << " sta " << b.fTrackingStation << ": x " << r[3 * is + 0] << " y " << r[3 * is + 1]
1060 << " z " << r[3 * is + 2];
1061 }
1062 }
1063
1064 unsigned active_tracks = 0;
1065 for (auto& t : fTracks) {
1066 if (t.fIsActive) active_tracks++;
1067 }
1068
1069 LOG(info) << "\n";
1070
1071 LOG(info) << "Bba: " << fTracks.size() << " tracks used for alignment\n";
1072 LOG(info) << " " << active_tracks << " tracks active\n";
1073
1074 LOG(info) << " Cost function for the true parameters: " << fCostIdeal;
1075 LOG(info) << " Cost function for the initial parameters: " << fCostInitial;
1076 LOG(info) << " Cost function for the aligned parameters: " << costResult;
1077
1078 for (auto& t : fTracks) {
1079 if (!t.fIsActive) continue;
1080 // calculate the residuals for all tracks at each TrajectoryNode after alignment
1081 // TODO: put into function
1082 for (unsigned int in = 0; in < t.fAlignedTrack.fNodes.size(); in++) {
1083 CbmKfTrackFitter::Trajectory tr = t.fAlignedTrack;
1085 if (!node.fIsXySet) {
1086 continue;
1087 }
1088 node.fIsXySet = false;
1090
1091 Sensor& s = fSensors[node.fReference1];
1092
1093 double x_residual = node.fMxy.X() - node.fParamDn.X();
1094 double y_residual = node.fMxy.Y() - node.fParamDn.Y();
1095 double x_pull = x_residual / sqrt(node.fMxy.Dx2() + node.fParamDn.GetCovariance(0, 0));
1096 double y_pull = y_residual / sqrt(node.fMxy.Dy2() + node.fParamDn.GetCovariance(1, 1));
1097
1098 if (node.fMxy.NdfX() > 0) {
1099 hResidualsAfterAlignmentX[0]->Fill(x_residual);
1100 hPullsAfterAlignmentX[0]->Fill(x_pull);
1101 }
1102 if (node.fMxy.NdfY() > 0) {
1103 hResidualsAfterAlignmentY[0]->Fill(y_residual);
1104 hPullsAfterAlignmentY[0]->Fill(y_pull);
1105 }
1106
1107 hResidualsAfterAlignmentX[s.fTrackingStation + 1]->Fill(x_residual);
1108 hResidualsAfterAlignmentY[s.fTrackingStation + 1]->Fill(y_residual);
1109
1110 hPullsAfterAlignmentX[s.fTrackingStation + 1]->Fill(x_pull);
1111 hPullsAfterAlignmentY[s.fTrackingStation + 1]->Fill(y_pull);
1112 }
1113 }
1114
1115 // store the histograms
1116
1117 TDirectory* curr = gDirectory;
1118 TFile* currentFile = gFile;
1119 // Open output file and write histograms
1120
1121 fHistoFile->cd();
1123 if (!(fHistoFileName == "")) {
1124 fHistoFile->Close();
1125 fHistoFile->Delete();
1126 }
1127 gFile = currentFile;
1128 gDirectory = curr;
1129}
1130
1132{
1133 if (!obj->IsFolder())
1134 obj->Write();
1135 else {
1136 TDirectory* cur = gDirectory;
1137 TFile* currentFile = gFile;
1138
1139 TDirectory* sub = cur->GetDirectory(obj->GetName());
1140 sub->cd();
1141 TList* listSub = (static_cast<TDirectory*>(obj))->GetList();
1142 TIter it(listSub);
1143 while (TObject* obj1 = it()) {
1144 WriteHistosCurFile(obj1);
1145 }
1146 cur->cd();
1147 gFile = currentFile;
1148 gDirectory = cur;
1149 }
1150}
1151
ClassImp(CbmBbaAlignmentTask)
std::pair< std::string, TGeoHMatrix > AlignNode(std::string path, double shiftX, double shiftY, double shiftZ, double rotX, double rotY, double rotZ)
Task class for alignment.
Time-slice/event reader for CA tracker in CBM (header)
@ kMvd
Micro-Vertex Detector.
@ kTrd
Transition Radiation Detector.
@ kTof
Time-of-flight Detector.
@ kSts
Silicon Tracking System.
@ kTrd2d
TRD-FASP Detector (FIXME)
@ kMuch
Muon detection system.
@ kStsSensor
Data class for a reconstructed hit in the STS.
Data class for STS tracks.
friend fvec sqrt(const fvec &a)
TClonesArray * fInputGlobalTrackMatches
void WriteHistosCurFile(TObject *obj)
std::vector< TrackContainer > fTracks
CbmBbaAlignmentTask(const char *name="CbmBbaAlignmentTask", Int_t iVerbose=0, TString histoFileName="CbmBbaAlignmentHisto.root")
std::vector< TH1F * > hPullsBeforeAlignmentX
TClonesArray * fInputStsTrackMatches
std::vector< TH1F * > hPullsAfterAlignmentX
TClonesArray * fInputMcTracks
TClonesArray * fInputStsTracks
std::vector< TH1F * > hResidualsAfterAlignmentX
std::vector< AlignmentBody > fAlignmentBodies
std::vector< TH1F * > hPullsBeforeAlignmentY
void ApplyAlignment(const std::vector< double > &par)
void ApplyConstraints(std::vector< double > &par)
TClonesArray * fInputGlobalTracks
std::vector< TH1F * > hResidualsAfterAlignmentY
std::vector< Sensor > fSensors
std::vector< TH1F * > hResidualsBeforeAlignmentX
void Exec(Option_t *opt)
std::vector< TH1F * > hResidualsBeforeAlignmentY
std::vector< TH1F * > hPullsAfterAlignmentY
double CostFunction(const std::vector< double > &par)
void ConstrainStation(std::vector< double > &par, int iSta, int ixyz)
void SetSkipUnmeasuredCoordinates(bool skip=true)
skip unmeasured coordinates
bool CreateMvdStsTrack(Trajectory &kfTrack, int stsTrackIndex)
set the input data arrays
bool CreateGlobalTrack(Trajectory &kfTrack, int globalTrackIndex)
void FixMomentumForMs(bool fix=true)
fix the inverse momentum for the Multiple Scattering calculation
void SetDefaultMomentumForMs(double p)
set the default inverse momentum for the Multiple Scattering calculation
void SetDoSmooth(bool doSmooth)
do the KF-smoothing to define track pars at all the nodes
bool FitTrajectory(CbmKfTrackFitter::Trajectory &t)
fit the track
void SetNoMultipleScattering()
set the default inverse momentum for the Multiple Scattering calculation
static CbmL1 * Instance()
Pointer to CbmL1 instance.
Definition CbmL1.h:171
ca::Framework * fpAlgo
Pointer to the L1 track finder algorithm.
Definition CbmL1.h:425
Class representing an element of the STS setup.
TGeoPhysicalNode * GetPnode() const
virtual void Print(Option_t *opt="") const
CbmStsElement * GetElement(Int_t address, Int_t level)
static CbmStsSetup * Instance()
static CbmStsTrackingInterface * Instance()
Gets pointer to the instance of the CbmStsTrackingInterface class.
int GetTrackingStationIndex(const CbmHit *hit) const override
Gets a tracking station of a CbmHit.
static int32_t GetRpcFullId(uint32_t address)
int GetNtrackingStations() const
Gets actual number of stations, provided by the current geometry setup.
const Parameters< fvec > & GetParameters() const
Gets a pointer to the Framework parameters object.
Definition CaFramework.h:87
A container for all external parameters of the CA tracking algorithm.
int GetNstationsActive() const
Gets total number of active stations.
T X() const
Gets x position [cm].
T Y() const
Gets y position [cm].
void SetX(T v)
Sets x position [cm].
T GetCovariance(int i, int j) const
Get covariance matrix element.
CbmKfTrackFitter::Trajectory fUnalignedTrack
CbmKfTrackFitter::Trajectory fAlignedTrack
cbm::algo::kf::TrackParamD fParamUp
fitted track parameters upstream the node
double fZ
Z coordinate of the node.
bool fIsFitted
true if the track parameters at the node are fitted
cbm::algo::kf::MeasurementXy< double > fMxy
== Hit information ( if present )
int fReference1
some reference can be set by the user
cbm::algo::kf::TrackParamD fParamDn
fitted track parameters downstream the node
bool fIsRadThickFixed
true if the radiation thickness is fixed to the fRadThick value
A trajectory to be fitted.
std::vector< TrajectoryNode > fNodes
nodes on the trajectory