CbmRoot
Loading...
Searching...
No Matches
CbmMatchRecoToMC.cxx
Go to the documentation of this file.
1/* Copyright (C) 2013-2023 GSI/JINR-LIT, Darmstadt/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Andrey Lebedev [committer], Florian Uhlig, Volker Friese, Pierre-Alain Loizeau */
4
5/*
6 * CbmMatchRecoToMC.cxx
7 *
8 * Created on: Oct 17, 2013
9 * Author: andrey
10 */
11#include "CbmMatchRecoToMC.h"
12
13#include "CbmCluster.h" // for CbmCluster
14#include "CbmDigiManager.h" // for CbmDigiManager
15#include "CbmFsdDigi.h"
16#include "CbmFsdHit.h" // for CbmFsdHit
17#include "CbmHit.h" // for CbmHit, kMUCHPIXELHIT
18#include "CbmLink.h" // for CbmLink
19#include "CbmMCDataArray.h" // for CbmMCDataArray
20#include "CbmMCDataManager.h" // for CbmMCDataManager
21#include "CbmMCTrack.h" // for CbmMCTrack
22#include "CbmMatch.h" // for CbmMatch
23#include "CbmPixelHit.h" // for CbmPixelHit
24#include "CbmRichDigi.h" // for CbmRichDigi
25#include "CbmRichHit.h" // for CbmRichHit
26#include "CbmRichPoint.h" // for CbmRichPoint
27#include "CbmRichRing.h" // for CbmRichRing
28#include "CbmStsAddress.h" // for GetSystemId
29#include "CbmStsHit.h" // for CbmStsHit
30#include "CbmStsTrack.h" // for CbmStsTrack
31#include "CbmTofDigi.h" // for CbmTofDigi
32#include "CbmTofHit.h" // for CbmTofHit
33#include "CbmTrack.h" // for CbmTrack
34#include "CbmTrackMatchNew.h" // for CbmTrackMatchNew
35
36#include <FairMCPoint.h> // for FairMCPoint
37#include <FairRootManager.h> // for FairRootManager
38#include <FairTask.h> // for FairTask, InitStatus
39#include <Logger.h> // for LOG, Logger
40
41#include <TClonesArray.h> // for TClonesArray
42
43#include <boost/exception/exception.hpp> // for clone_impl
44#include <boost/type_index/type_index_facade.hpp> // for operator==
45
46#include <algorithm> // for find
47#include <iomanip> // for std::setw
48
49using std::pair;
50using std::vector;
51
53
54CbmMatchRecoToMC::CbmMatchRecoToMC() : FairTask("MatchRecoToMC") {}
55
57{
58
59 if (fStsClusterMatches != nullptr) {
60 fStsClusterMatches->Delete();
61 delete fStsClusterMatches;
62 }
63 if (fStsHitMatches != nullptr) {
64 fStsHitMatches->Delete();
65 delete fStsHitMatches;
66 }
67 if (fStsTrackMatches) {
68 fStsTrackMatches->Delete();
69 delete fStsTrackMatches;
70 }
71
73 fRichTrackMatches->Delete();
74 delete fRichTrackMatches;
75 }
76
77 if (fTrdClusterMatches != nullptr) {
78 fTrdClusterMatches->Delete();
79 delete fTrdClusterMatches;
80 }
81 if (fTrdHitMatches != nullptr) {
82 fTrdHitMatches->Delete();
83 delete fTrdHitMatches;
84 }
85 if (fTrdTrackMatches) {
86 fTrdTrackMatches->Delete();
87 delete fTrdTrackMatches;
88 }
89
90 if (fMuchClusterMatches != nullptr) {
91 fMuchClusterMatches->Delete();
93 }
94 if (fMuchPixelHitMatches != nullptr) {
95 fMuchPixelHitMatches->Delete();
97 }
99 fMuchTrackMatches->Delete();
100 delete fMuchTrackMatches;
101 }
102
103 if (fMvdClusterMatches != nullptr) {
104 fMvdClusterMatches->Delete();
105 delete fMvdClusterMatches;
106 }
107 if (fMvdHitMatches != nullptr) {
108 fMvdHitMatches->Delete();
109 delete fMvdHitMatches;
110 }
111
112 if (fTofHitMatches != nullptr) {
113 fTofHitMatches->Delete();
114 delete fTofHitMatches;
115 }
116
117 if (fFsdHitMatches != nullptr) {
118 fFsdHitMatches->Delete();
119 delete fFsdHitMatches;
120 }
121}
122
124{
126
127 // Notification to a user about matching suppression
129 std::stringstream msg;
130 msg << "\033[1;31mCbmMatchRecoToMC (\"" << fName << "\"): the cluster and hit matching routines for MVD, STS, ";
131 msg
132 << "MuCh, TRD, TOF, FSD will be suppressed by a request from CbmMatchRecoToMC::SuppressHitReMatching():\033[0m\n";
133 LOG(warn) << msg.str();
134 }
135
136 return kSUCCESS;
137}
138
139void CbmMatchRecoToMC::Exec(Option_t* /*opt*/)
140{
142 if (fMvdHitMatches != nullptr) fMvdHitMatches->Delete();
143 if (fMvdClusterMatches != nullptr) fMvdClusterMatches->Delete();
144 if (fStsClusterMatches != nullptr) fStsClusterMatches->Delete();
145 if (fStsHitMatches != nullptr) fStsHitMatches->Delete();
146 if (fTrdClusterMatches != nullptr) fTrdClusterMatches->Delete();
147 if (fTrdHitMatches != nullptr) fTrdHitMatches->Delete();
148 if (fMuchClusterMatches != nullptr) fMuchClusterMatches->Delete();
149 if (fMuchPixelHitMatches != nullptr) fMuchPixelHitMatches->Delete();
150 if (fTofHitMatches != nullptr) fTofHitMatches->Delete();
151 if (fFsdHitMatches != nullptr) fFsdHitMatches->Delete();
152 }
153 if (fTrdClusterMatches != nullptr) fTrdClusterMatches->Delete();
154 if (fStsTrackMatches != nullptr) fStsTrackMatches->Delete();
155 if (fRichTrackMatches != nullptr) fRichTrackMatches->Delete();
156 if (fTrdTrackMatches != nullptr) fTrdTrackMatches->Delete();
157 if (fMuchTrackMatches != nullptr) fMuchTrackMatches->Delete();
158
159
160 // ----- MVD -----
162 // MVD: digi->hit
165 }
166 else {
167 // MVD: digi->cluster
169 // MVD: cluster->hit
171 }
172 }
173
174 // ----- STS -----
176 // STS: digi->cluster
178 // STS: cluster->hit
180 }
181 // STS: hit->track
183
184 // ----- MUCH -----
186 // MUCH: digi->cluster
188 // MUCH: cluster->hit
190 }
191 // MUCH: hit->track
193
194
195 //RICH
198 }
199
200 // TRD
201 if (fTrdClusters && fTrdHits) { // MC->digi->cluster->hit->track
205 }
207 }
208 else if (fTrdHits) { // MC->hit->track
211 }
213 }
214
215 // TOF: (Digi->MC)+(Hit->Digi)=>(Hit->MC)
218 }
219
220 // ----- FSD -----
222 // FSD: digi->hit
223 if (fFsdHits && fFsdHitMatches) {
225 }
226 }
227
228
229 //static Int_t eventNo = 0;
230 LOG(info) << "CbmMatchRecoToMC::Exec eventNo=" << fEventNumber++;
231}
232
234
236{
237 FairRootManager* ioman = FairRootManager::Instance();
238 if (nullptr == ioman) {
239 LOG(fatal) << "CbmMatchRecoToMC::ReadAndCreateDataBranches() NULL FairRootManager.";
240 }
241
242 CbmMCDataManager* mcManager = (CbmMCDataManager*) ioman->GetObject("MCDataManager");
243
244 if (nullptr == mcManager) {
245 LOG(fatal) << "CbmMatchRecoToMC::ReadAndCreateDataBranches() NULL MCDataManager.";
246 }
247
248 fMCTracks = mcManager->InitBranch("MCTrack");
249
250 //fMCTracksArray= (TClonesArray*) ioman->GetObject("MCTrack");
251
252 // --- Digi Manager
255
256
257 // Register output branches only if they don't exist
258 // This should solve the problem with BranchName_1, BranchName_2
259 // when having two instances of CbmMatchRecoToMC in the same macro
260
261 // STS
262 fStsPoints = mcManager->InitBranch("StsPoint");
263
264 fStsClusters = static_cast<TClonesArray*>(ioman->GetObject("StsCluster"));
265 if (nullptr != fStsClusters) {
266 fStsClusterMatches = static_cast<TClonesArray*>(ioman->GetObject("StsClusterMatch"));
267 if (nullptr == fStsClusterMatches) {
268 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
269 fStsClusterMatches = new TClonesArray("CbmMatch", 100);
270 ioman->Register("StsClusterMatch", "STS", fStsClusterMatches, IsOutputBranchPersistent("StsClusterMatch"));
271 }
272 }
273
274 fStsHits = static_cast<TClonesArray*>(ioman->GetObject("StsHit"));
275 if (nullptr != fStsHits) {
276 fStsHitMatches = static_cast<TClonesArray*>(ioman->GetObject("StsHitMatch"));
277 if (nullptr == fStsHitMatches) {
278 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
279 fStsHitMatches = new TClonesArray("CbmMatch", 100);
280 ioman->Register("StsHitMatch", "STS", fStsHitMatches, IsOutputBranchPersistent("StsHitMatch"));
281 }
282 }
283
284 fStsTracks = static_cast<TClonesArray*>(ioman->GetObject("StsTrack"));
285 if (nullptr != fStsTracks) {
286 fStsTrackMatches = static_cast<TClonesArray*>(ioman->GetObject("StsTrackMatch"));
287 if (nullptr == fStsTrackMatches) {
288 fStsTrackMatches = new TClonesArray("CbmTrackMatchNew", 100);
289 ioman->Register("StsTrackMatch", "STS", fStsTrackMatches, IsOutputBranchPersistent("StsTrackMatch"));
290 }
291 }
292
293 //RICH
294 fRichMcPoints = mcManager->InitBranch("RichPoint");
295 fRichHits = static_cast<TClonesArray*>(ioman->GetObject("RichHit"));
296
297 fRichRings = static_cast<TClonesArray*>(ioman->GetObject("RichRing"));
298 if (nullptr != fRichRings) {
299 fRichTrackMatches = static_cast<TClonesArray*>(ioman->GetObject("RichRingMatch"));
300 if (nullptr == fRichTrackMatches) {
301 fRichTrackMatches = new TClonesArray("CbmTrackMatchNew", 100);
302 ioman->Register("RichRingMatch", "RICH", fRichTrackMatches, IsOutputBranchPersistent("RichRingMatch"));
303 }
304 }
305
306 // TRD
307 fTrdPoints = mcManager->InitBranch("TrdPoint");
308 fTrdClusters = static_cast<TClonesArray*>(ioman->GetObject("TrdCluster"));
309 if (nullptr != fTrdClusters) {
310 fTrdClusterMatches = static_cast<TClonesArray*>(ioman->GetObject("TrdClusterMatch"));
311 if (nullptr == fTrdClusterMatches) {
312 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
313 fTrdClusterMatches = new TClonesArray("CbmMatch", 100);
314 ioman->Register("TrdClusterMatch", "TRD", fTrdClusterMatches, IsOutputBranchPersistent("TrdClusterMatch"));
315 }
316 }
317
318 fTrdHits = static_cast<TClonesArray*>(ioman->GetObject("TrdHit"));
319 if (nullptr != fTrdHits) {
320 fTrdHitMatches = static_cast<TClonesArray*>(ioman->GetObject("TrdHitMatch"));
321 if (nullptr == fTrdHitMatches) {
322 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
323 fTrdHitMatches = new TClonesArray("CbmMatch", 100);
324 ioman->Register("TrdHitMatch", "TRD", fTrdHitMatches, IsOutputBranchPersistent("TrdHitMatch"));
325 }
326 }
327
328
329 fTrdTracks = static_cast<TClonesArray*>(ioman->GetObject("TrdTrack"));
330 if (nullptr != fTrdTracks) {
331 fTrdTrackMatches = static_cast<TClonesArray*>(ioman->GetObject("TrdTrackMatch"));
332 if (nullptr == fTrdTrackMatches) {
333 fTrdTrackMatches = new TClonesArray("CbmTrackMatchNew", 100);
334 ioman->Register("TrdTrackMatch", "TRD", fTrdTrackMatches, IsOutputBranchPersistent("TrdTrackMatch"));
335 }
336 }
337
338 // MUCH
339 fMuchPoints = mcManager->InitBranch("MuchPoint");
340
341 fMuchTracks = (TClonesArray*) ioman->GetObject("MuchTrack");
342 fMuchClusters = static_cast<TClonesArray*>(ioman->GetObject("MuchCluster"));
343 if (nullptr != fMuchClusters) {
344 fMuchClusterMatches = static_cast<TClonesArray*>(ioman->GetObject("MuchClusterMatch"));
345 if (nullptr == fMuchClusterMatches) {
346 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
347 fMuchClusterMatches = new TClonesArray("CbmMatch", 100);
348 ioman->Register("MuchClusterMatch", "MUCH", fMuchClusterMatches, IsOutputBranchPersistent("MuchClusterMatch"));
349 }
350 }
351
352 fMuchPixelHits = static_cast<TClonesArray*>(ioman->GetObject("MuchPixelHit"));
353 if (nullptr != fMuchPixelHits) {
354 fMuchPixelHitMatches = static_cast<TClonesArray*>(ioman->GetObject("MuchPixelHitMatch"));
355 if (nullptr == fMuchPixelHitMatches) {
356 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
357 fMuchPixelHitMatches = new TClonesArray("CbmMatch", 100);
358 ioman->Register("MuchPixelHitMatch", "MUCH", fMuchPixelHitMatches, IsOutputBranchPersistent("MuchPixelHitMatch"));
359 }
360 }
361
362 fMuchTracks = static_cast<TClonesArray*>(ioman->GetObject("MuchTrack"));
363 if (nullptr != fMuchTracks) {
364 fMuchTrackMatches = static_cast<TClonesArray*>(ioman->GetObject("MuchTrackMatch"));
365 if (nullptr == fMuchTrackMatches) {
366 fMuchTrackMatches = new TClonesArray("CbmTrackMatchNew", 100);
367 ioman->Register("MuchTrackMatch", "MUCH", fMuchTrackMatches, IsOutputBranchPersistent("MuchTrackMatch"));
368 }
369 }
370
371 // MVD
372 fMvdPoints = mcManager->InitBranch("MvdPoint");
373
374 fMvdCluster = static_cast<TClonesArray*>(ioman->GetObject("MvdCluster"));
375 if (nullptr != fMvdCluster) {
376 fMvdClusterMatches = static_cast<TClonesArray*>(ioman->GetObject("MvdClusterMatch"));
377 if (nullptr == fMvdClusterMatches) {
378 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
379 fMvdClusterMatches = new TClonesArray("CbmMatch", 100);
380 ioman->Register("MvdClusterMatch", "MVD", fMvdClusterMatches, IsOutputBranchPersistent("MvdClusterMatch"));
381 }
382 }
383
384 fMvdHits = static_cast<TClonesArray*>(ioman->GetObject("MvdHit"));
385 if (nullptr != fMvdHits) {
386 fMvdHitMatches = static_cast<TClonesArray*>(ioman->GetObject("MvdHitMatch"));
387 if (nullptr == fMvdHitMatches) {
388 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
389 fMvdHitMatches = new TClonesArray("CbmMatch", 100);
390 ioman->Register("MvdHitMatch", "MVD", fMvdHitMatches, IsOutputBranchPersistent("MvdHitMatch"));
391 }
392 }
393
395
396 // Currently in the time-based mode MVD is present but not reconstructed
397 // TODO: remove the check once the reconstruction works
398 if (fIsMvdActive && !fMvdHits) {
399 LOG(warning) << "CbmMatchRecoToMC: MVD hits are missing, MVD will not be "
400 "included to the STS track match";
401 fIsMvdActive = false;
402 }
403
404
405 // TOF
406 fTofPoints = mcManager->InitBranch("TofPoint");
407 fTofHitDigiMatches = static_cast<TClonesArray*>(ioman->GetObject("TofHitDigiMatch"));
408
411 if (nullptr == fTofHitDigiMatches) {
412 fTofHitDigiMatches = static_cast<TClonesArray*>(ioman->GetObject("TofHitCalDigiMatch"));
413 if (nullptr == fTofHitDigiMatches) {
414 LOG(warning) << "CbmMatchRecoToMC::ReadAndCreateDataBranches()"
415 << " no TOF Hit to Digi array found!";
416 } // if (nullptr == fTofHitDigiMatches)
417 else {
418 LOG(info) << "CbmMatchRecoToMC: Using alternative TOF Hit to Digi match array!";
419 } // else of if (nullptr == fTofHitDigiMatches)
420 } // if (nullptr == fTofHitDigiMatches)
421
422 fTofHits = static_cast<TClonesArray*>(ioman->GetObject("TofHit"));
423 if (nullptr != fTofHits) {
424 fTofHitMatches = static_cast<TClonesArray*>(ioman->GetObject("TofHitMatch"));
425 if (nullptr == fTofHitMatches) {
426 fTofHitMatches = new TClonesArray("CbmMatch", 100);
427 ioman->Register("TofHitMatch", "TOF", fTofHitMatches, IsOutputBranchPersistent("TofHitMatch"));
428 }
429 }
430
432 fTofDigis = ioman->InitObjectAs<std::vector<CbmTofDigi> const*>("TofCalDigi");
433 if (nullptr == fTofDigis) {
434 LOG(info) << "No calibrated tof digi vector in the input files => trying original vector";
435 fTofDigis = ioman->InitObjectAs<std::vector<CbmTofDigi> const*>("TofDigi");
436 if (nullptr == fTofDigis) {
437 LOG(info) << "No original tof digi vector in the input files! Ignore TOF!";
438 }
439 }
440 else
441 LOG(info) << "Found calibrated tof digi vector in one of the input files";
442
443 if (nullptr != fTofDigis) {
444 fTofDigiMatch = ioman->InitObjectAs<std::vector<CbmMatch> const*>("TofCalDigiMatch");
445 if (nullptr == fTofDigiMatch) {
446 LOG(info) << "No calibrated tof digi to point match vector in the input files => trying original vector";
447 fTofDigiMatch = ioman->InitObjectAs<std::vector<CbmMatch> const*>("TofDigiMatch");
448 if (nullptr == fTofDigiMatch) {
449 LOG(fatal) << "No original tof digi to point match vector in the input files!";
450 }
451 }
452 else
453 LOG(info) << "Found calibrated tof digi to point match vector in one of the input files";
454 }
455
456 // FSD
457 fFsdPoints = mcManager->InitBranch("FsdPoint");
458
459 fFsdHits = static_cast<TClonesArray*>(ioman->GetObject("FsdHit"));
460 if (nullptr != fFsdHits) {
461 fFsdHitMatches = static_cast<TClonesArray*>(ioman->GetObject("FsdHitMatch"));
462 if (nullptr == fFsdHitMatches) {
463 fbSuppressHitReMatching = false; // Ensure, that matching is done, if matches are absent
464 fFsdHitMatches = new TClonesArray("CbmMatch", 100);
465 ioman->Register("FsdHitMatch", "FSD", fFsdHitMatches, IsOutputBranchPersistent("FsdHitMatch"));
466 }
467 }
468}
469
470
471void CbmMatchRecoToMC::MatchClusters(const TClonesArray* digiMatches, const TClonesArray* clusters,
472 TClonesArray* clusterMatches)
473{
474 if (!(digiMatches && clusters && clusterMatches)) return;
475 Int_t nofClusters = clusters->GetEntriesFast();
476 for (Int_t iCluster = 0; iCluster < nofClusters; iCluster++) {
477 CbmCluster* cluster = static_cast<CbmCluster*>(clusters->At(iCluster));
478 CbmMatch* clusterMatch = new ((*clusterMatches)[iCluster]) CbmMatch();
479 Int_t nofDigis = cluster->GetNofDigis();
480 for (Int_t iDigi = 0; iDigi < nofDigis; iDigi++) {
481 const CbmMatch* digiMatch = static_cast<const CbmMatch*>(digiMatches->At(cluster->GetDigi(iDigi)));
482 clusterMatch->AddLinks(*digiMatch);
483 } //# digis in cluster
484 } //# clusters
485}
486
487void CbmMatchRecoToMC::MatchClusters(ECbmModuleId systemId, const TClonesArray* clusters, TClonesArray* clusterMatches)
488{
489 if (!fDigiManager->IsMatchPresent(systemId)) return;
490 if (!(clusters && clusterMatches)) return;
491 Int_t nofClusters = clusters->GetEntriesFast();
492 for (Int_t iCluster = 0; iCluster < nofClusters; iCluster++) {
493 CbmCluster* cluster = static_cast<CbmCluster*>(clusters->At(iCluster));
494 CbmMatch* clusterMatch = new ((*clusterMatches)[iCluster]) CbmMatch();
495 Int_t nofDigis = cluster->GetNofDigis();
496 for (Int_t iDigi = 0; iDigi < nofDigis; iDigi++) {
497 Int_t digiIndex = cluster->GetDigi(iDigi);
498 const CbmMatch* digiMatch = fDigiManager->GetMatch(systemId, digiIndex);
499 if (nullptr == digiMatch) {
500 LOG(fatal) << "CbmMatchRecoToMC::MatchClusters => no Match found for system " << systemId << " digi index "
501 << digiIndex << " (digi " << iDigi << " from cluster " << iCluster << ") !";
502 }
503 clusterMatch->AddLinks(*digiMatch);
504 } //# digis in cluster
505 } //# clusters
506}
507
508void CbmMatchRecoToMC::MatchHits(const TClonesArray* matches, const TClonesArray* hits, TClonesArray* hitMatches)
509{
510 if (!(matches && hits && hitMatches)) return;
511 Int_t nofHits = hits->GetEntriesFast();
512 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
513 CbmHit* hit = static_cast<CbmHit*>(hits->At(iHit));
514 CbmMatch* hitMatch = new ((*hitMatches)[iHit]) CbmMatch();
515 const CbmMatch* clusterMatch = static_cast<const CbmMatch*>(matches->At(hit->GetRefId()));
516 hitMatch->AddLinks(*clusterMatch);
517 }
518}
519
520void CbmMatchRecoToMC::MatchHitsSts(const TClonesArray* cluMatches, const TClonesArray* hits, TClonesArray* hitMatches)
521{
522 if (!(cluMatches && hits && hitMatches)) return;
523 Int_t nofHits = hits->GetEntriesFast();
524 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
525 CbmStsHit* hit = static_cast<CbmStsHit*>(hits->At(iHit));
526 CbmMatch* hitMatch = new ((*hitMatches)[iHit]) CbmMatch();
527 const CbmMatch* frontClusterMatch = static_cast<const CbmMatch*>(cluMatches->At(hit->GetFrontClusterId()));
528 const CbmMatch* backClusterMatch = static_cast<const CbmMatch*>(cluMatches->At(hit->GetBackClusterId()));
529
530 for (int iLinkF = 0; iLinkF < frontClusterMatch->GetNofLinks(); ++iLinkF) {
531 const auto& linkF = frontClusterMatch->GetLink(iLinkF);
532 for (int iLinkB = 0; iLinkB < backClusterMatch->GetNofLinks(); ++iLinkB) {
533 const auto& linkB = backClusterMatch->GetLink(iLinkB);
534 if (linkB == linkF) {
535 hitMatch->AddLink(linkF);
536 hitMatch->AddLink(linkB);
537 }
538 }
539 }
540 }
541}
542
543void CbmMatchRecoToMC::MatchHitsMvd(const TClonesArray* hits, TClonesArray* hitMatches)
544{
545 if (!(hits && hitMatches)) return;
546 Int_t nofHits = hits->GetEntriesFast();
547 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
548 CbmPixelHit* hit = static_cast<CbmPixelHit*>(hits->At(iHit));
549 CbmMatch* hitMatch = new ((*hitMatches)[iHit]) CbmMatch();
550 const CbmMatch* digiMatch = fDigiManager->GetMatch(ECbmModuleId::kMvd, hit->GetRefId());
551 hitMatch->AddLinks(*digiMatch);
552 }
553}
554
555void CbmMatchRecoToMC::MatchHitsFsd(const TClonesArray* hits, TClonesArray* hitMatches)
556{
557 if (!(hits && hitMatches)) return;
558 Int_t nofHits = hits->GetEntriesFast();
559 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
560 CbmPixelHit* hit = static_cast<CbmPixelHit*>(hits->At(iHit));
561 CbmMatch* hitMatch = new ((*hitMatches)[iHit]) CbmMatch();
563 for (Int_t iDigi = 0; iDigi < nofDigis; iDigi++) {
564 const CbmFsdDigi* digi = fDigiManager->Get<const CbmFsdDigi>(iDigi);
565 Int_t adresa_digi = digi->GetAddress();
566 if (adresa_digi == hit->GetAddress()) {
567 const CbmMatch* digiMatch = fDigiManager->GetMatch(ECbmModuleId::kFsd, iDigi);
568 hitMatch->AddLinks(*digiMatch);
569 }
570 }
571 }
572}
573
574void CbmMatchRecoToMC::MatchHitsTof(const TClonesArray* HitDigiMatches, const TClonesArray* hits,
575 TClonesArray* hitMatches)
576{
577 if (!(HitDigiMatches && hits && hitMatches)) return;
578
579 Int_t iNbTofDigis = fTofDigis->size();
580 Int_t iNbTofDigiMatches = fTofDigiMatch->size();
581 if (iNbTofDigis != iNbTofDigiMatches)
582 LOG(fatal) << "CbmMatchRecoToMC::MatchHitsTof => Nb digis in vector not matching nb matches in vector: "
583 << iNbTofDigis << " VS " << iNbTofDigiMatches;
584
585 Int_t nofHits = hits->GetEntriesFast();
586 CbmTofDigi digiTof;
587 CbmMatch matchDigiPnt;
588
589 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
590 CbmMatch* hitDigiMatch = static_cast<CbmMatch*>(HitDigiMatches->At(iHit));
591 CbmMatch* hitMatch = new ((*hitMatches)[iHit]) CbmMatch();
592
593 Int_t iNbDigisHit = hitDigiMatch->GetNofLinks();
594 for (Int_t iDigi = 0; iDigi < iNbDigisHit; iDigi++) {
595 CbmLink lDigi = hitDigiMatch->GetLink(iDigi);
596 if (lDigi.IsNoise()) continue;
597 Int_t iDigiIdx = lDigi.GetIndex();
598
599 if (iNbTofDigis <= iDigiIdx) {
600 LOG(error) << "CbmMatchRecoToMC::MatchHitsTof => Digi index from Hit #" << iHit << "/" << nofHits << " Digi "
601 << iDigi << "/" << iNbDigisHit << " is bigger than nb entries in Digis arrays: " << iDigiIdx
602 << " VS " << iNbTofDigis << " => ignore it!!!";
603 const CbmTofHit* pTofHit = static_cast<CbmTofHit*>(hits->At(iHit));
604 LOG(error) << " Hit position: ( " << pTofHit->GetX() << " , "
605 << pTofHit->GetY() << " , " << pTofHit->GetZ() << " ) ";
606 continue;
607 } // if( iNbTofDigis <= iDigiIdx )
608
609 digiTof = fTofDigis->at(iDigiIdx);
610 matchDigiPnt = fTofDigiMatch->at(iDigiIdx);
611 Int_t iNbPointsDigi = matchDigiPnt.GetNofLinks();
612 if (iNbPointsDigi <= 0) {
613 LOG(error) << "CbmMatchRecoToMC::MatchHitsTof => No entries in Digi to point match for Hit #" << iHit << "/"
614 << nofHits << " Digi " << iDigi << "/" << iNbDigisHit << ": " << iNbPointsDigi << " (digi index is "
615 << iDigiIdx << "/" << iNbTofDigis << ") => ignore it!!!";
616 LOG(error) << " Digi address: 0x" << std::setw(8) << std::hex
617 << digiTof.GetAddress() << std::dec;
618 continue;
619 } // if( iNbTofDigis <= iDigiIdx )
620 CbmLink lTruePoint = matchDigiPnt.GetMatchedLink(); // Point generating the Digi
621 if (lTruePoint.IsNoise()) continue;
622 Int_t iTruePointIdx = lTruePoint.GetIndex();
623 for (Int_t iPoint = 0; iPoint < iNbPointsDigi; iPoint++) {
624 CbmLink lPoint = matchDigiPnt.GetLink(iPoint);
625 if (lPoint.IsNoise()) continue;
626 Int_t iPointIdx = lPoint.GetIndex();
627
628 if (iPointIdx == iTruePointIdx)
629 hitMatch->AddLink(CbmLink(digiTof.GetTot(), iPointIdx, lPoint.GetEntry(),
630 lPoint.GetFile())); // Point generating the Digi
631 else
632 hitMatch->AddLink(CbmLink(0, iPointIdx, lPoint.GetEntry(),
633 lPoint.GetFile())); // Point whose Digi was hidden by True one
634 } // for( Int_t iPoint = 0; iPoint < iNbPointsDigi; iPoint ++)
635 } // for (Int_t iDigi = 0; iDigi < iNbDigisHit; iDigi++)
636 } // for (Int_t iHit = 0; iHit < nofHits; iHit++)
637}
638
639void CbmMatchRecoToMC::MatchHitsToPoints(CbmMCDataArray* points, const TClonesArray* hits, TClonesArray* hitMatches)
640{
641 if (!(hits && hitMatches)) return;
642 Int_t nofHits = hits->GetEntriesFast();
643 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
644 CbmHit* hit = static_cast<CbmHit*>(hits->At(iHit));
645 CbmMatch* hitMatch = new ((*hitMatches)[iHit]) CbmMatch();
646 const FairMCPoint* point = static_cast<const FairMCPoint*>(points->Get(0, fEventNumber, hit->GetRefId()));
647 hitMatch->AddLink(CbmLink(point->GetEnergyLoss(), hit->GetRefId(), fEventNumber));
648 }
649}
650
651void CbmMatchRecoToMC::MatchTracks(const TClonesArray* hitMatches, CbmMCDataArray* points, const TClonesArray* tracks,
652 TClonesArray* trackMatches)
653{
654 if (!(hitMatches && points && tracks && trackMatches)) return;
655
656 Bool_t editMode = (trackMatches->GetEntriesFast() != 0);
657
658 Int_t nofTracks = tracks->GetEntriesFast();
659 for (Int_t iTrack = 0; iTrack < nofTracks; iTrack++) {
660 const CbmTrack* track = static_cast<const CbmTrack*>(tracks->At(iTrack));
661 CbmTrackMatchNew* trackMatch = (editMode) ? static_cast<CbmTrackMatchNew*>(trackMatches->At(iTrack))
662 : new ((*trackMatches)[iTrack]) CbmTrackMatchNew();
663 Int_t nofHits = track->GetNofHits();
664 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
665 if ((track->GetHitType(iHit) == kMUCHPIXELHIT)
666 || (track->GetHitType(iHit) == kMUCHSTRAWHIT && hitMatches == fMuchPixelHitMatches))
667 continue; //AZ
668 const CbmMatch* hitMatch = static_cast<CbmMatch*>(hitMatches->At(track->GetHitIndex(iHit)));
669 Int_t nofLinks = hitMatch->GetNofLinks();
670 for (Int_t iLink = 0; iLink < nofLinks; iLink++) {
671 const FairMCPoint* point = static_cast<const FairMCPoint*>(points->Get(hitMatch->GetLink(iLink)));
672 if (nullptr == point) continue;
674 if (CbmStsAddress::GetSystemId(point->GetDetectorID()) == ECbmModuleId::kSts) {
675 Int_t mcTrackId = point->GetTrackID();
676 CbmMCTrack* mcTrack = (CbmMCTrack*) fMCTracks->Get(hitMatch->GetLink(iLink).GetFile(),
677 hitMatch->GetLink(iLink).GetEntry(), mcTrackId);
678 if (mcTrack->GetNPoints(ECbmModuleId::kSts) < 2) continue;
679 }
681 trackMatch->AddLink(CbmLink(1., point->GetTrackID(), hitMatch->GetLink(iLink).GetEntry()));
682 }
683 }
684 if (!trackMatch->GetNofLinks()) continue;
685 // Calculate number of true and wrong hits
686 Int_t trueCounter = trackMatch->GetNofTrueHits();
687 Int_t wrongCounter = trackMatch->GetNofWrongHits();
688 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
689 if ((track->GetHitType(iHit) == kMUCHPIXELHIT)
690 || (track->GetHitType(iHit) == kMUCHSTRAWHIT && hitMatches == fMuchPixelHitMatches))
691 continue; //AZ
692 const CbmMatch* hitMatch = static_cast<CbmMatch*>(hitMatches->At(track->GetHitIndex(iHit)));
693 Int_t nofLinks = hitMatch->GetNofLinks();
694 Bool_t hasTrue = false;
695 for (Int_t iLink = 0; iLink < nofLinks; iLink++) {
696 const FairMCPoint* point = static_cast<const FairMCPoint*>(points->Get(hitMatch->GetLink(iLink)));
697 if (nullptr == point) continue;
698 if (
699 /*point->GetEventID() == trackMatch->GetMatchedLink().GetEntry() && */
700 point->GetTrackID() == trackMatch->GetMatchedLink().GetIndex()) {
701 hasTrue = true;
702 break;
703 }
704 }
705 if (hasTrue)
706 trueCounter++;
707 else
708 wrongCounter++;
709 }
710 trackMatch->SetNofTrueHits(trueCounter);
711 trackMatch->SetNofWrongHits(wrongCounter);
712 // LOG(debug) << iTrack << " "; track->Print(); LOG(debug) << " " << trackMatch->ToString();
713 }
714}
715
716void CbmMatchRecoToMC::MatchStsTracks(const TClonesArray* mvdHitMatches, const TClonesArray* stsHitMatches,
717 CbmMCDataArray* mvdPoints, CbmMCDataArray* stsPoints, const TClonesArray* tracks,
718 TClonesArray* trackMatches)
719{
720 if (!tracks) {
721 return;
722 }
723
724 if (!(stsHitMatches && stsPoints && tracks && trackMatches)) {
725 LOG(error) << "CbmMatchRecoToMC: missing the necessary STS info for the "
726 "STS track match";
727 return;
728 }
729
730 if (fIsMvdActive && !(mvdHitMatches && mvdPoints)) {
731 LOG(error) << "CbmMatchRecoToMC: missing the necessary MVD info for the "
732 "STS track match";
733 return;
734 }
735
736 // make sure that the Mvd hits are not present in the reconstruced tracks
737 // when we assume that the Mvd reconstruction was disabled
738 if (!fIsMvdActive) {
739 for (Int_t iTrack = 0; iTrack < trackMatches->GetEntriesFast(); iTrack++) {
740 const CbmStsTrack* track = static_cast<const CbmStsTrack*>(tracks->At(iTrack));
741 if (track->GetNofMvdHits() > 0) {
742 LOG(error) << "CbmMatchRecoToMC: Sts track contains Mvd hits, but "
743 "there is no MVD data available";
744 return;
745 }
746 }
747 }
748
749 Bool_t editMode = (trackMatches->GetEntriesFast() != 0);
750
751 Int_t nofTracks = tracks->GetEntriesFast();
752
753 for (Int_t iTrack = 0; iTrack < nofTracks; iTrack++) {
754 const CbmStsTrack* track = static_cast<const CbmStsTrack*>(tracks->At(iTrack));
755 CbmTrackMatchNew* trackMatch = (editMode) ? static_cast<CbmTrackMatchNew*>(trackMatches->At(iTrack))
756 : new ((*trackMatches)[iTrack]) CbmTrackMatchNew();
757
758 Int_t nofStsHits = track->GetNofStsHits();
759 for (Int_t iHit = 0; iHit < nofStsHits; iHit++) {
760 const CbmMatch* hitMatch = static_cast<CbmMatch*>(stsHitMatches->At(track->GetStsHitIndex(iHit)));
761 Int_t nofLinks = hitMatch->GetNofLinks();
762 for (Int_t iLink = 0; iLink < nofLinks; iLink++) {
763 const CbmLink& link = hitMatch->GetLink(iLink);
764 if (link.IsNoise()) continue;
765 const FairMCPoint* point = static_cast<const FairMCPoint*>(stsPoints->Get(link));
766 assert(point);
768 if (CbmStsAddress::GetSystemId(point->GetDetectorID()) == ECbmModuleId::kSts) {
769 Int_t mcTrackId = point->GetTrackID();
770 CbmMCTrack* mcTrack = (CbmMCTrack*) fMCTracks->Get(link.GetFile(), link.GetEntry(), mcTrackId);
771 assert(mcTrack);
772 if (mcTrack->GetNPoints(ECbmModuleId::kSts) < 2) {
773 continue;
774 }
775 }
777 trackMatch->AddLink(CbmLink(1., point->GetTrackID(), link.GetEntry(), link.GetFile()));
778 }
779 }
780
781 // Include MVD hits in STS matching if MVD is active
782 Int_t nofMvdHits = track->GetNofMvdHits();
783 for (Int_t iHit = 0; iHit < nofMvdHits; iHit++) {
784 const CbmMatch* hitMatch = static_cast<CbmMatch*>(mvdHitMatches->At(track->GetMvdHitIndex(iHit)));
785 Int_t nofLinks = hitMatch->GetNofLinks();
786 for (Int_t iLink = 0; iLink < nofLinks; iLink++) {
787 const CbmLink& link = hitMatch->GetLink(iLink);
788 if (link.IsNoise()) continue;
789 const FairMCPoint* point = static_cast<const FairMCPoint*>(mvdPoints->Get(link));
790 assert(point);
791 trackMatch->AddLink(CbmLink(1., point->GetTrackID(), link.GetEntry(), link.GetFile()));
792 }
793 }
794
795 if (!trackMatch->GetNofLinks()) {
796 continue;
797 }
798
799 // Calculate number of true and wrong hits
800 Int_t trueCounter = trackMatch->GetNofTrueHits();
801 Int_t wrongCounter = trackMatch->GetNofWrongHits();
802 for (Int_t iHit = 0; iHit < nofStsHits; iHit++) {
803 const CbmMatch* hitMatch = static_cast<CbmMatch*>(stsHitMatches->At(track->GetStsHitIndex(iHit)));
804 Int_t nofLinks = hitMatch->GetNofLinks();
805 Bool_t hasTrue = false;
806 for (Int_t iLink = 0; iLink < nofLinks; iLink++) {
807 const FairMCPoint* point = static_cast<const FairMCPoint*>(stsPoints->Get(hitMatch->GetLink(iLink)));
808 if (nullptr == point) {
809 continue;
810 }
811 if (point->GetTrackID() == trackMatch->GetMatchedLink().GetIndex()) {
812 hasTrue = true;
813 break;
814 }
815 }
816 if (hasTrue) {
817 trueCounter++;
818 }
819 else {
820 wrongCounter++;
821 }
822 }
823
824 // Include MVD hits in STS matching if MVD is active
825 for (Int_t iHit = 0; iHit < nofMvdHits; iHit++) {
826 const CbmMatch* hitMatch = static_cast<CbmMatch*>(mvdHitMatches->At(track->GetMvdHitIndex(iHit)));
827 Int_t nofLinks = hitMatch->GetNofLinks();
828 Bool_t hasTrue = false;
829 for (Int_t iLink = 0; iLink < nofLinks; iLink++) {
830 const FairMCPoint* point = static_cast<const FairMCPoint*>(mvdPoints->Get(hitMatch->GetLink(iLink)));
831 if (nullptr == point) {
832 continue;
833 }
834 if (
835 /*point->GetEventID() == trackMatch->GetMatchedLink().GetEntry() && */
836 point->GetTrackID() == trackMatch->GetMatchedLink().GetIndex()) {
837 hasTrue = true;
838 break;
839 }
840 }
841 if (hasTrue) {
842 trueCounter++;
843 }
844 else {
845 wrongCounter++;
846 }
847 }
848
849 trackMatch->SetNofTrueHits(trueCounter);
850 trackMatch->SetNofWrongHits(wrongCounter);
851 // LOG(debug) << iTrack << " "; track->Print(); LOG(debug) << " " << trackMatch->ToString();
852 }
853}
854
855
856void CbmMatchRecoToMC::MatchRichRings(const TClonesArray* richRings, const TClonesArray* richHits,
857 CbmMCDataArray* richMcPoints, CbmMCDataArray* mcTracks, TClonesArray* ringMatches)
858{
859 // Loop over RichRings
860 Int_t nRings = richRings->GetEntriesFast();
861 for (Int_t iRing = 0; iRing < nRings; iRing++) {
862 const CbmRichRing* ring = static_cast<const CbmRichRing*>(richRings->At(iRing));
863 if (!ring) continue;
864 CbmTrackMatchNew* ringMatch = new ((*ringMatches)[iRing]) CbmTrackMatchNew();
865 Int_t nHits = ring->GetNofHits();
866 for (Int_t iHit = 0; iHit < nHits; iHit++) {
867 const CbmRichHit* hit = static_cast<const CbmRichHit*>(richHits->At(static_cast<Int_t>(ring->GetHit(iHit))));
868 if (!hit) continue;
869 vector<CbmLink> motherIds = GetMcTrackMotherIdsForRichHit(fDigiManager, hit, richMcPoints, mcTracks);
870 for (const auto& motherId : motherIds) {
871 ringMatch->AddLink(motherId);
872 }
873 }
874
875 if (ringMatch->GetNofLinks() != 0) {
876 CbmLink bestTrackLink = ringMatch->GetMatchedLink();
877 if (bestTrackLink.IsNoise()) continue;
878 Int_t trueCounter = 0;
879 Int_t wrongCounter = 0;
880 for (Int_t iHit = 0; iHit < nHits; iHit++) {
881 const CbmRichHit* hit = static_cast<const CbmRichHit*>(richHits->At(static_cast<Int_t>(ring->GetHit(iHit))));
882 if (!hit) continue;
883 vector<CbmLink> motherIds = GetMcTrackMotherIdsForRichHit(fDigiManager, hit, richMcPoints, mcTracks);
884 if (std::find(motherIds.begin(), motherIds.end(), bestTrackLink) != motherIds.end()) {
885 trueCounter++;
886 }
887 else {
888 wrongCounter++;
889 }
890 }
891 ringMatch->SetNofTrueHits(trueCounter);
892 ringMatch->SetNofWrongHits(wrongCounter);
893 }
894 else {
895 ringMatch->SetNofTrueHits(0);
896 ringMatch->SetNofWrongHits(0);
897 }
898 } // Ring loop
899}
900
903{
904 vector<CbmLink> result;
905 if (!hit) return result;
906 Int_t digiIndex = hit->GetRefId();
907 if (digiIndex < 0) return result;
908 const CbmRichDigi* digi = digiMan->Get<CbmRichDigi>(digiIndex);
909 if (!digi) return result;
910 const CbmMatch* digiMatch = digiMan->GetMatch(ECbmModuleId::kRich, digiIndex);
911 if (!digiMatch) return result;
912
913 vector<CbmLink> links = digiMatch->GetLinks();
914 for (const auto& link : links) {
915 if (link.IsNoise()) continue;
916 const CbmRichPoint* point = static_cast<const CbmRichPoint*>(richPoints->Get(link));
917 if (!point) continue;
918 Int_t mcTrackId = point->GetTrackID();
919 if (mcTrackId < 0) continue;
920 const CbmMCTrack* mcTrack =
921 static_cast<const CbmMCTrack*>(mcTracks->Get(link.GetFile(), link.GetEntry(), mcTrackId));
922 if (!mcTrack) continue;
923 if (mcTrack->GetPdgCode() != 50000050) continue; // select only Cherenkov photons
924 Int_t motherId = mcTrack->GetMotherId();
925 // several photons can have same mother track
926 // count only unique motherIds
927 CbmLink val(1., motherId, link.GetEntry(), link.GetFile());
928 if (std::find(result.begin(), result.end(), val) == result.end()) result.push_back(val);
929 }
930
931 return result;
932}
933
935 const CbmRichHit* hit,
936 CbmMCDataArray* richPoints,
938 Int_t /*eventNumber*/)
939{
940 vector<pair<Int_t, Int_t>> result;
941 if (nullptr == hit) return result;
942 Int_t digiIndex = hit->GetRefId();
943 if (digiIndex < 0) return result;
944 const CbmRichDigi* digi = digiMan->Get<CbmRichDigi>(digiIndex);
945 if (nullptr == digi) return result;
946 const CbmMatch* digiMatch = digiMan->GetMatch(ECbmModuleId::kRich, digiIndex);
947 if (digiMatch == nullptr) return result;
948
949 vector<CbmLink> links = digiMatch->GetLinks();
950 for (UInt_t i = 0; i < links.size(); i++) {
951 if (links[i].IsNoise()) continue;
952 Int_t pointId = links[i].GetIndex();
953 Int_t eventId = links[i].GetEntry();
954 const CbmRichPoint* pMCpt = static_cast<const CbmRichPoint*>(richPoints->Get(0, eventId, pointId));
955
956 if (nullptr == pMCpt) continue;
957 Int_t mcTrackIndex = pMCpt->GetTrackID();
958 if (mcTrackIndex < 0) continue;
959 //TODO: Currently we support only legacy mode of CbmMCDataArray
960 const CbmMCTrack* mcTrack = static_cast<const CbmMCTrack*>(mcTracks->Get(0, eventId, mcTrackIndex));
961 // CbmMCTrack* pMCtr = (CbmMCTrack*) mcTracks->At(mcTrackIndex);
962 if (nullptr == mcTrack) continue;
963 if (mcTrack->GetPdgCode() != 50000050) continue; // select only Cherenkov photons
964 Int_t motherId = mcTrack->GetMotherId();
965 // several photons can have same mother track
966 // count only unique motherIds
967 pair<Int_t, Int_t> val = std::make_pair(eventId, motherId);
968 if (std::find(result.begin(), result.end(), val) == result.end()) {
969 result.push_back(val);
970 }
971 }
972
973 return result;
974}
975
977 const TClonesArray* richPoints,
978 const TClonesArray* mcTracks)
979{
980 vector<Int_t> result;
981 if (nullptr == hit) return result;
982 Int_t digiIndex = hit->GetRefId();
983 if (digiIndex < 0) return result;
984 const CbmRichDigi* digi = digiMan->Get<CbmRichDigi>(digiIndex);
985 if (nullptr == digi) return result;
986 const CbmMatch* digiMatch = digiMan->GetMatch(ECbmModuleId::kRich, digiIndex);
987 if (digiMatch == nullptr) return result;
988
989 vector<CbmLink> links = digiMatch->GetLinks();
990 for (UInt_t i = 0; i < links.size(); i++) {
991 if (links[i].IsNoise()) continue;
992 Int_t pointId = links[i].GetIndex();
993 const CbmRichPoint* pMCpt = static_cast<const CbmRichPoint*>(richPoints->At(pointId));
994 if (nullptr == pMCpt) continue;
995 Int_t mcTrackIndex = pMCpt->GetTrackID();
996 if (mcTrackIndex < 0) continue;
997 //TODO: Currently we support only legacy mode of CbmMCDataArray
998 const CbmMCTrack* mcTrack = static_cast<const CbmMCTrack*>(mcTracks->At(mcTrackIndex));
999 // CbmMCTrack* pMCtr = (CbmMCTrack*) mcTracks->At(mcTrackIndex);
1000 if (nullptr == mcTrack) continue;
1001 if (mcTrack->GetPdgCode() != 50000050) continue; // select only Cherenkov photons
1002 Int_t motherId = mcTrack->GetMotherId();
1003 // several photons can have same mother track
1004 // count only unique motherIds
1005 if (std::find(result.begin(), result.end(), motherId) == result.end()) {
1006 result.push_back(motherId);
1007 }
1008 }
1009
1010 return result;
1011}
1012
1013
TClonesArray * tracks
TClonesArray * points
Base class for cluster objects.
ECbmModuleId
Definition CbmDefs.h:39
@ kMvd
Micro-Vertex Detector.
@ kTrd
Transition Radiation Detector.
@ kSts
Silicon Tracking System.
@ kMuch
Muon detection system.
@ kFsd
Forward spectator detector.
@ kRich
Ring-Imaging Cherenkov Detector.
@ kMUCHSTRAWHIT
Definition CbmHit.h:29
@ kMUCHPIXELHIT
Definition CbmHit.h:28
ClassImp(CbmMatchRecoToMC)
FairTask for matching RECO data to MC.
Data class for a reconstructed hit in the STS.
Data class for STS tracks.
static vector< vector< QAMCTrack > > mcTracks
static vector< vector< QAHit > > hits
Base class for cluster objects.
Definition CbmCluster.h:30
int32_t GetDigi(int32_t index) const
Get digi at position index.
Definition CbmCluster.h:76
int32_t GetNofDigis() const
Number of digis in cluster.
Definition CbmCluster.h:69
CbmDigiManager.
static Int_t GetNofDigis(ECbmModuleId systemId)
static Bool_t IsMatchPresent(ECbmModuleId systemId)
Presence of a digi match branch.
InitStatus Init()
Initialisation.
const Digi * Get(Int_t index) const
Get a digi object.
static CbmDigiManager * Instance()
Static instance.
const CbmMatch * GetMatch(ECbmModuleId systemId, UInt_t index) const
Get a match object.
Data class for FSD digital information.
Definition CbmFsdDigi.h:36
uint32_t GetAddress() const
Address.
Definition CbmFsdDigi.h:80
int32_t GetAddress() const
Definition CbmHit.h:74
double GetZ() const
Definition CbmHit.h:71
int32_t GetRefId() const
Definition CbmHit.h:73
Access to a MC data branch for time-based analysis.
TObject * Get(const CbmLink *lnk)
Task class creating and managing CbmMCDataArray objects.
CbmMCDataObject * GetObject(const char *name)
CbmMCDataArray * InitBranch(const char *name)
int32_t GetMotherId() const
Definition CbmMCTrack.h:69
int32_t GetPdgCode() const
Definition CbmMCTrack.h:68
int32_t GetNPoints(ECbmModuleId detId) const
const std::vector< CbmTofDigi > * fTofDigis
CbmTofPoint array.
TClonesArray * fMuchPixelHits
Clusters [in].
TClonesArray * fTrdTrackMatches
Hit matches [out].
CbmMCDataArray * fRichMcPoints
Track matches [out].
void MatchRichRings(const TClonesArray *richRings, const TClonesArray *richHits, CbmMCDataArray *richMcPoints, CbmMCDataArray *mcTracks, TClonesArray *ringMatches)
TClonesArray * fMvdCluster
MC points [in].
CbmMCDataArray * fMvdPoints
Interface to digi branches.
TClonesArray * fStsHitMatches
Cluster matches [out].
void MatchHitsTof(const TClonesArray *HitDigiMatches, const TClonesArray *hits, TClonesArray *hitMatches)
CbmMCDataArray * fTrdPoints
Track matches [out].
void MatchHits(const TClonesArray *matches, const TClonesArray *hits, TClonesArray *hitMatches)
CbmMCDataArray * fStsPoints
Hit matches [out].
void MatchHitsToPoints(CbmMCDataArray *points, const TClonesArray *hits, TClonesArray *hitMatches)
CbmMCDataArray * fTofPoints
Track matches [out].
TClonesArray * fTrdClusters
MC points [in].
TClonesArray * fMuchPixelHitMatches
Cluster matches [out].
TClonesArray * fTrdHits
Clusters [in].
TClonesArray * fMuchTracks
Hits [in].
void ReadAndCreateDataBranches()
Read and create data branches.
void MatchHitsSts(const TClonesArray *clusterMmatches, const TClonesArray *hits, TClonesArray *hitMatches)
Match STS hits, using cluster match objects.
void MatchTracks(const TClonesArray *hitMatches, CbmMCDataArray *points, const TClonesArray *tracks, TClonesArray *trackMatches)
TClonesArray * fRichTrackMatches
Rings [in].
void MatchClusters(const TClonesArray *digiMatches, const TClonesArray *clusters, TClonesArray *clusterMatches)
Generic creation of cluster match objects.
TClonesArray * fStsTrackMatches
Hit matches [out].
void MatchHitsFsd(const TClonesArray *hits, TClonesArray *hitMatches)
virtual InitStatus Init()
Derived from FairTask.
TClonesArray * fTofHitMatches
Match Hit -> Digi [out].
TClonesArray * fFsdHitMatches
Hits [in].
TClonesArray * fMvdHits
Clusters [in].
CbmDigiManager * fDigiManager
Monte-Carlo tracks.
TClonesArray * fRichHits
MC points [in].
TClonesArray * fMuchTrackMatches
Hit matches [out].
TClonesArray * fFsdHits
MC points [in].
TClonesArray * fMvdClusterMatches
Hits [in].
TClonesArray * fTrdTracks
Hits [in].
TClonesArray * fTofHits
virtual void Exec(Option_t *opt)
Derived from FairTask.
void MatchHitsMvd(const TClonesArray *hits, TClonesArray *hitMatches)
CbmMCDataArray * fMCTracks
const std::vector< CbmMatch > * fTofDigiMatch
TClonesArray * fStsHits
Clusters [in].
static std::vector< CbmLink > GetMcTrackMotherIdsForRichHit(CbmDigiManager *digiMan, const CbmRichHit *hit, CbmMCDataArray *richPoints, CbmMCDataArray *mcTracks)
Get CbmLinks of Mother MC Tracks for RICH hit.
TClonesArray * fStsTracks
Hits [in].
TClonesArray * fTrdHitMatches
Cluster matches [out].
TClonesArray * fTofHitDigiMatches
CbmTofHit array.
TClonesArray * fRichRings
Hits [in].
TClonesArray * fStsClusterMatches
Tracks [in].
TClonesArray * fStsClusters
MC points [in].
void MatchStsTracks(const TClonesArray *mvdHitMatches, const TClonesArray *stsHitMatches, CbmMCDataArray *mvdPoints, CbmMCDataArray *stsPoints, const TClonesArray *tracks, TClonesArray *trackMatches)
virtual void Finish()
Derived from FairTask.
static Int_t fEventNumber
CbmMatchRecoToMC()
Constructor.
TClonesArray * fMuchClusters
MC points [in].
TClonesArray * fTrdClusterMatches
Tracks [in].
TClonesArray * fMuchClusterMatches
Tracks [in].
CbmMCDataArray * fFsdPoints
Match Hit -> MC point [out].
TClonesArray * fMvdHitMatches
Cluster matches [out].
virtual ~CbmMatchRecoToMC()
Destructor.
CbmMCDataArray * fMuchPoints
Match Ring -> MC track [out].
void AddLinks(const CbmMatch &match)
Definition CbmMatch.cxx:39
const CbmLink & GetLink(int32_t i) const
Definition CbmMatch.h:39
int32_t GetNofLinks() const
Definition CbmMatch.h:42
void AddLink(const CbmLink &newLink)
Definition CbmMatch.cxx:47
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
const std::vector< CbmLink > & GetLinks() const
Definition CbmMatch.h:40
double GetY() const
Definition CbmPixelHit.h:74
double GetX() const
Definition CbmPixelHit.h:73
uint32_t GetHit(int32_t i) const
Definition CbmRichRing.h:39
int32_t GetNofHits() const
Definition CbmRichRing.h:37
data class for a reconstructed 3-d hit in the STS
Definition CbmStsHit.h:35
int32_t GetFrontClusterId() const
Definition CbmStsHit.h:105
int32_t GetBackClusterId() const
Definition CbmStsHit.h:65
int32_t GetNofMvdHits() const
Definition CbmStsTrack.h:87
int32_t GetMvdHitIndex(int32_t iHit) const
Definition CbmStsTrack.h:75
int32_t GetStsHitIndex(int32_t iHit) const
int32_t GetNofStsHits() const
Definition CbmStsTrack.h:93
Data class for expanded digital TOF information.
Definition CbmTofDigi.h:47
int32_t GetAddress() const
Inherited from CbmDigi.
Definition CbmTofDigi.h:112
double GetTot() const
Alias for GetCharge.
Definition CbmTofDigi.h:140
void SetNofWrongHits(int32_t nofWrongHits)
int32_t GetNofWrongHits() const
void SetNofTrueHits(int32_t nofTrueHits)
int32_t GetNofTrueHits() const
virtual int32_t GetNofHits() const
Definition CbmTrack.h:58
int32_t GetHitIndex(int32_t iHit) const
Definition CbmTrack.h:59
HitType GetHitType(int32_t iHit) const
Definition CbmTrack.h:60
ECbmModuleId GetSystemId(int32_t address)
Get system Id (should be ECbmModuleId::kSts)