CbmRoot
Loading...
Searching...
No Matches
CbmTrdTrackFinderIdeal.cxx
Go to the documentation of this file.
1/* Copyright (C) 2005-2019 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese, Semen Lebedev, Denis Bertini [committer], Florian Uhlig */
4
5// -------------------------------------------------------------------------
6// ----- CbmTrdTrackFinderIdeal source file -----
7// ----- Created 28/11/05 by V. Friese -----
8// ----- according to the CbmStsTrackFinderIdeal -----
9// -------------------------------------------------------------------------
10// CBM includes
12
13#include "CbmDigiManager.h"
14#include "CbmMCTrack.h"
15#include "CbmMatch.h"
16#include "CbmTrdHit.h"
17#include "CbmTrdTrack.h"
18#include "FairBaseParSet.h"
19#include "FairDetector.h"
20#include "FairMCPoint.h"
21#include "FairRootManager.h"
22#include "FairRunAna.h"
23#include "FairRuntimeDb.h"
24
25// ROOT includes
26#include "TClonesArray.h"
27#include "TObjArray.h"
28
29// C++ includes
30#include <iostream>
31#include <map>
32using std::cout;
33using std::endl;
34using std::map;
35
36
39 , fMcTracks(NULL)
40 , fTrdPoints(NULL)
41 , fTrdHitMatches(NULL)
42 , fTrdHitProducerType("")
43 , fEventNum(-1)
44{
45}
46
47
49
50
52{
53 FairRootManager* ioman = FairRootManager::Instance();
54 if (NULL == ioman) LOG(fatal) << GetName() << " Init: RootManager not instantised!";
55
56 fMcTracks = (TClonesArray*) ioman->GetObject("MCTrack");
57 if (NULL == fMcTracks) LOG(fatal) << GetName() << " Init: No MCTrack array!";
58
59 fTrdPoints = (TClonesArray*) ioman->GetObject("TrdPoint");
60 if (NULL == fTrdPoints) LOG(fatal) << GetName() << " Init: No TrdPoint array!";
61
62 TClonesArray* digis = (TClonesArray*) ioman->GetObject("TrdDigi");
63 if (NULL == digis) {
64 fTrdHitProducerType = "smearing";
65 }
66 else {
67 fTrdHitProducerType = "digi";
68 fTrdHitMatches = (TClonesArray*) ioman->GetObject("TrdHitMatch");
69 if (NULL == fTrdHitMatches) LOG(fatal) << GetName() << "Init: No TrdHitMatch array!";
70 }
71
72 fEventNum = 1;
73}
74
75Int_t CbmTrdTrackFinderIdeal::DoFind(TClonesArray* trdHits, TClonesArray* trdTracks)
76{
77 cout << "CbmTrdTrackFinderIdeal, event no. " << fEventNum << " " << fTrdHitProducerType << endl;
78 fEventNum++;
79 if (NULL == trdHits) LOG(fatal) << GetName() << " DoFind: No TrdHitArray!";
80 if (NULL == trdTracks) LOG(fatal) << GetName() << " DoFind: No TrdTrackArray!";
81
82 Int_t nofTrdHits = trdHits->GetEntriesFast();
83
84 // STL map from MCTrack index to CbmTrdTrack
85 map<Int_t, CbmTrdTrack*> mcTrackToTrdTrack;
86 for (Int_t iHit = 0; iHit < nofTrdHits; iHit++) {
87 CbmTrdHit* trdHit = static_cast<CbmTrdHit*>(trdHits->At(iHit));
88 if (NULL == trdHit) continue;
89
90 int trdPointInd = -1;
91 if (fTrdHitProducerType == "smearing") {
92 trdPointInd = trdHit->GetRefId();
93 }
94 else if (fTrdHitProducerType == "digi") {
95 const CbmMatch* trdHitMatch = static_cast<const CbmMatch*>(fTrdHitMatches->At(iHit));
96 trdPointInd = trdHitMatch->GetMatchedLink().GetIndex();
97 }
98
99 if (trdPointInd < 0) continue; // fake or background hit
100 FairMCPoint* trdPoint = static_cast<FairMCPoint*>(fTrdPoints->At(trdPointInd));
101 if (NULL == trdPoint) continue;
102
103 int mcTrackInd = trdPoint->GetTrackID();
104 if (NULL == mcTrackToTrdTrack[mcTrackInd]) mcTrackToTrdTrack[mcTrackInd] = new CbmTrdTrack();
105 mcTrackToTrdTrack[mcTrackInd]->AddHit(iHit, kTRDHIT);
106 }
107
108 int trackCount = 0;
109 for (map<int, CbmTrdTrack*>::iterator it = mcTrackToTrdTrack.begin(); it != mcTrackToTrdTrack.end(); ++it) {
110 new ((*trdTracks)[trackCount]) CbmTrdTrack(*it->second);
111 trackCount++;
112 }
113
114 return mcTrackToTrdTrack.size();
115}
116
ClassImp(CbmConverterManager)
@ kTRDHIT
Definition CbmHit.h:30
Class for hits in TRD detector.
int32_t GetRefId() const
Definition CbmHit.h:73
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
data class for a reconstructed Energy-4D measurement in the TRD
Definition CbmTrdHit.h:40
virtual Int_t DoFind(TClonesArray *hitArray, TClonesArray *trackArray)