PandaRoot
PndTrackingQualityBarrelAnalysisNewLinks.h
Go to the documentation of this file.
1 /*
2  * PndTrackingQualityBarrelAnalysisNewLinks.h
3  *
4  * Created on: Aug 23, 2013
5  * Author: stockman
6  */
7 
8 #ifndef PndTrackingQualityBarrelAnalysisNewLinks_H_
9 #define PndTrackingQualityBarrelAnalysisNewLinks_H_
10 
11 #include "FairMultiLinkedData.h"
12 #include "FairRootManager.h"
13 
14 #include "PndTrackCand.h"
16 
17 #include <TObject.h>
18 #include <TString.h>
19 #include <TClonesArray.h>
20 
21 #include <vector>
22 #include <map>
23 
24 #include <functional>
25 #include "PndTrackingQualityAnalysis.h"
26 
27 /*class PossibleTrackFunctor : public std::binary_function<FairMultiLinkedData* , Bool_t, Bool_t>
28 {
29  public :
30  virtual Bool_t operator() (FairMultiLinkedData* a, Bool_t primary) {return Call(a, primary);};
31  virtual Bool_t Call(FairMultiLinkedData* a, Bool_t primary) = 0;
32  virtual void Print() = 0;
33 
34  virtual ~PossibleTrackFunctor() {};
35 
36 };
37 
38 class StandardTrackFunctor : public PossibleTrackFunctor
39 {
40  Bool_t Call(FairMultiLinkedData* a, Bool_t primary){
41  FairRootManager* ioman = FairRootManager::Instance();
42  Bool_t possibleTrack = kFALSE;
43  possibleTrack = possibleTrack | (a->GetLinksWithType(ioman->GetBranchId("MVDHitsPixel")).GetNLinks() +
44  a->GetLinksWithType(ioman->GetBranchId("MVDHitsStrip")).GetNLinks() > 3);
45 
46  possibleTrack = possibleTrack | (a->GetLinksWithType(ioman->GetBranchId("MVDHitsPixel")).GetNLinks() +
47  a->GetLinksWithType(ioman->GetBranchId("MVDHitsStrip")).GetNLinks() +
48  a->GetLinksWithType(ioman->GetBranchId("STTHit")).GetNLinks() ) > 5;
49 
50  return possibleTrack;
51  }
52 
53  void Print(){
54  std::cout << "StandardTrackFunctor: > 3 Hits in MVD or > 5 Hits in (MVD+Stt)" << std::endl;
55  }
56 };
57 
58 class OnlySttFunctor : public PossibleTrackFunctor
59 {
60  Bool_t Call(FairMultiLinkedData* a, Bool_t primary){
61  FairRootManager* ioman = FairRootManager::Instance();
62  Bool_t possibleTrack = kFALSE;
63 
64  possibleTrack = possibleTrack | a->GetLinksWithType(ioman->GetBranchId("STTHit")).GetNLinks() > 5;
65 
66  return possibleTrack;
67  }
68  void Print(){
69  std::cout << "OnlySttFunctor: > 5 Hits in Stt" << std::endl;
70  }
71 };
72 
73 class RiemannMvdSttGemFunctor : public PossibleTrackFunctor
74 {
75  Bool_t Call(FairMultiLinkedData* a, Bool_t primary){
76  FairRootManager* ioman = FairRootManager::Instance();
77  Bool_t possibleTrack = kFALSE;
78  Bool_t mvdHits = ((a->GetLinksWithType(ioman->GetBranchId("MVDHitsPixel")).GetNLinks() +
79  a->GetLinksWithType(ioman->GetBranchId("MVDHitsStrip")).GetNLinks()) > 2);
80 
81  if (mvdHits){
82  possibleTrack = a->GetLinksWithType(ioman->GetBranchId("STTHit")).GetNLinks() > 1 | a->GetLinksWithType(ioman->GetBranchId("GEMHit")).GetNLinks() > 1;
83  }
84  return possibleTrack;
85  }
86  void Print(){
87  std::cout << "RiemannMvdSttGemFunctor: > 2 Hits in MVD and >0 Hits in (Stt+Gem)" << std::endl;
88  }
89 
90 };
91 
92 class CircleHoughTrackFunctor : public PossibleTrackFunctor
93 {
94  Bool_t Call(FairMultiLinkedData* a, Bool_t primary){
95  if (primary == kFALSE) return kFALSE;
96  FairRootManager* ioman = FairRootManager::Instance();
97  Bool_t possibleTrack = kFALSE;
98 
99  Int_t nHitsMvdPixel = a->GetLinksWithType(ioman->GetBranchId("MVDHitsPixel")).GetNLinks();
100  Int_t nHitsMvdStrip = a->GetLinksWithType(ioman->GetBranchId("MVDHitsStrip")).GetNLinks();
101  if (nHitsMvdPixel + nHitsMvdStrip > 2) { // First requirement: more than two MVD hits
102  Int_t nHitsStt = a->GetLinksWithType(ioman->GetBranchId("STTHit")).GetNLinks();
103  Int_t nHitsGem = a->GetLinksWithType(ioman->GetBranchId("GEMHit")).GetNLinks();
104  possibleTrack = (nHitsMvdPixel + nHitsMvdStrip + nHitsStt + nHitsGem > 6); // Second requirement: More than six hits total
105  }
106  return possibleTrack;
107  }
108  void Print(){
109  std::cout << "CircleHoughTrackFunctor: > 3 Hits in MVD and primary track" << std::endl;
110  }
111 
112 };
113 */
114 
125 /*
126 struct qualityNumbers {
127  static const int
128  // Following: Default statuses.
129  // Are the 'not found' tracks in the quality histogram of PndTrackingQualityTask.
130  kPossibleSec = -1, // possible: As defined through the possibleFunctor; secondary: a non-primary particle
131  kPossiblePrim = -2, // possible: As defined through the possibleFunctor; primary: coming directly from particle generator (e.g. EvtGen)
132  kAtLeastThreeSec = -3, // atLeastThree: min. 3 hit points in central tracking detectors (MVD, STT, GEM); secondary: a non-primary particle
133  kAtLeastThreePrim = -4, // atLeastThree: min. 3 hit points in central tracking detectors (MVD, STT, GEM); primary: coming directly from particle generator (e.g. EvtGen)
134  kLessThanThreePrim = -5, // LessThanThree: fewer than 3 hit points in central tracking detectors (MVD, STT, GEM); primar: coming directly from particle generator (e.g. EvtGen)
135 
136  // Following: MC statuses of all (found+notfound) tracks
137  kMcPossibleSec = -7, // see above
138  kMcPossiblePrim = -8,
139  kMcAtLeastThreeSec = -9,
140  kMcAtLeastThreePrim = -10,
141  kMcLessThanThreePrim = -11,
142 
143  // Following: Status of reconstructed tracks (= created PndTracks)
144  kSpuriousFound = 1, // spuriousFound: at least 70% of hits of reco'd track come from one MC track ('mostProbableTrack')
145  kPartiallyFound = 2, // partiallyFound: all hits of reco'd track come from one single MC track; at least 70% of hits of MC track have been found in reco'd track
146  kFullyFound = 3, // fullyFound: all hits of reco'd track come from one single MC track; all hits of MC track have been found in reco'd track
147 
148  kGhost = 5, // ghost: less than 70% of hits of reco'd track come from one MC track ('mostProbableTrack')
149 
150  kNotFound = 7, // notFound: total number of not reco'd tracks
151  kFound = 8; // found: total number of reco'd tracks; the sum of fullyFound, partiallyFound, spuriousFound
152 };
153 
154 */
155 
157  public:
158  PndTrackingQualityBarrelAnalysisNewLinks(TString trackBranchName, TString idealTrackName, Bool_t pndTrackData = kTRUE);
159  PndTrackingQualityBarrelAnalysisNewLinks(TString trackBranchName, TString idealTrackName, PossibleTrackFunctor *posTrack, Bool_t pndTrackData = kTRUE);
161 
162  virtual void Init();
163  void SetVerbose(Int_t val) { fVerbose = val; }
164 
166  void AddHitsBranchName(TString name) { fBranchNames.push_back(name); }
167  void SetHitsBranchNames(std::vector<TString> names) { fBranchNames = names; }
168 
169  void AnalyseEvent(TClonesArray *recoTrackInfo);
170 
171  // Int_t GetNIdealHits(Int_t trackId, TString branchName);
172  Int_t GetNIdealHits(FairMultiLinkedData &track, TString branchName);
173  std::map<Int_t, Int_t> GetMCTrackFound() { return fMCTrackFound; }
174  std::map<Int_t, Int_t> GetTrackQualification() { return fMapTrackQualification; }
175  std::map<Int_t, Int_t> GetTrackMCStatus() { return fMapTrackMCStatus; }
176  std::map<Int_t, std::map<TString, std::pair<Double_t, Int_t>>> GetEfficiencies() { return fMapEfficiencies; }
177  std::map<Int_t, Double_t> GetPResolution() { return fMapPResolution; }
178  std::map<Int_t, TVector3> GetP() { return fMapP; }
179  std::map<Int_t, Double_t> GetPtResolution() { return fMapPtResolution; }
180  std::map<Int_t, Double_t> GetPt() { return fMapPt; }
181  std::map<Int_t, Double_t> GetPResolutionRel() { return fMapPResolutionRel; }
182  std::map<Int_t, Double_t> GetPtResolutionRel() { return fMapPtResolutionRel; }
183  std::map<Int_t, Int_t> GetTrackIdMCId() { return fTrackIdMCId; }
184  Int_t GetNGhosts() { return fNGhosts; }
185 
186  void PrintTrackDataSummary(FairMultiLinkedData &trackData, Bool_t detailedInfo = kFALSE);
187 
191  void PrintTrackQualityMap(Bool_t detailedInfo = kFALSE);
192  void PrintTrackMCStatusMap();
193  void PrintTrackInfo(std::map<TString, FairMultiLinkedData> info);
194 
195  Int_t GetIdealTrackIdFromMCTrackId(int mctrackid) { return fMCIdIdealTrackId[mctrackid]; }
197  {
198  int mctrackid = fTrackIdMCId[trackid];
199  return fMCIdIdealTrackId[mctrackid];
200  }
201 
202  PndTrackingQualityRecoInfo GetRecoInfoFromRecoTrack(Int_t trackId, Int_t mctrackId);
203 
204  static Bool_t IsBarrelMVD(FairMultiLinkedData &links, int iHit);
205  static Int_t NBarrelMVD(FairMultiLinkedData &links);
206 
207  private:
208  virtual void FillMapTrackQualifikation();
209  Bool_t IsBetterTrackExisting(Int_t &mcIndex, int quality);
210  // virtual Bool_t PossibleTrack(FairMultiLinkedData& mcForward);
211  Int_t GetSumOfAllValidMCHits(FairMultiLinkedData *trackData);
212  virtual Int_t AnalyseTrackInfo(std::map<TString, FairMultiLinkedData> &trackInfo, Int_t trackId);
213  virtual void CalcEfficiencies(Int_t mostProbableTrack, std::map<TString, FairMultiLinkedData> &trackInfo);
214  FairMultiLinkedData GetMCInfoForBranch(TString branchName, PndTrackCand *trackCand);
215  std::map<TString, FairMultiLinkedData> AnalyseTrackCand(PndTrackCand *trackCand);
216 
217  // virtual Bool_t IsCorrectGemHit(FairLink& gemLink);
218 
219  Int_t fVerbose;
220  FairRootManager *ioman;
221  Int_t fNGhosts;
222 
223  TString fTrackBranchName;
224  TString fIdealTrackName;
225  Bool_t fPndTrackOrTrackCand; // kTRUE if track and kFALSE if track cand
226  PossibleTrackFunctor *fPossibleTrack;
227 
228  Bool_t fUseCorrectedSkewedHits;
229 
230  std::vector<TString> fBranchNames;
231  std::map<Int_t, Int_t> fTrackIdMCId;
232  std::map<Int_t, Int_t> fMCIdTrackId;
233  std::map<Int_t, Int_t> fMCIdIdealTrackId;
234 
235  std::map<Int_t, Int_t> fMCTrackFound;
236 
237  std::map<Int_t, Int_t> fMapTrackMCStatus;
238  std::map<Int_t, Int_t> fMapTrackQualification;
239  std::map<Int_t, std::map<TString, std::pair<Double_t, Int_t>>> fMapEfficiencies;
240  std::map<Int_t, Double_t> fMapPResolution;
241  std::map<Int_t, TVector3> fMapP;
242  std::map<Int_t, Double_t> fMapPtResolution;
243  std::map<Int_t, Double_t> fMapPt;
244  std::map<Int_t, Double_t> fMapPResolutionRel;
245  std::map<Int_t, Double_t> fMapPtResolutionRel;
246 
247  TClonesArray *fTrack;
248  TClonesArray *fMCTrack;
249  TClonesArray *fIdealTrack;
250  TClonesArray *fIdealTrackCand;
251 
253 };
254 
255 #endif /* PNDTRACKINGQUALITY_H_ */