]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/macros/trackMatching/MatchComparison.C
AliAlignObjAngles becomes AliAlignObjParams (Raffaele)
[u/mrichter/AliRoot.git] / EMCAL / macros / trackMatching / MatchComparison.C
CommitLineData
e10611f0 1//
2// This macros compares the matches found by algorithm with the true matches
3// found with the "SaveTrueMatchesSimple.C" macro.
4// It saves 4 histogram, which contain the Pt distribution of:
5// - all found matches
6// - all correctly found matches
7// - all wrong (fake) found matches
8// - all true matches
9// which will then be availabel for computing efficiency and contamination.
10//
11
12class match_t
13{
14 public:
15
16 Int_t label; // GEANT label of particle
17 Int_t indexT; // index of track in ESD collection
18 Int_t indexC; // index of cluster in ESD collection
19 Double_t p[3]; // track momentum
20 Double_t v[3]; // track vertex
21};
22
23void MatchComparison()
24{
25 //
26 // Initialize AliRun manager
27 //
28 if (gAlice) {
29 delete gAlice;
30 gAlice = 0;
31 }
32
33 //
34 // Initialize run loader and load Kinematics
35 //
36 AliRunLoader *runLoader = AliRunLoader::Open("galice.root");
37 if (!runLoader) return;
38 runLoader->LoadgAlice();
39 gAlice = runLoader->GetAliRun();
40 runLoader->LoadKinematics();
41
42 //
43 // Initialize histograms with their error computation
44 //
45 TH1D *hgood = new TH1D("hgood", "Well matched tracks", 20, 0.0, 10.0);
46 TH1D *hfake = new TH1D("hfake", "Fake matched tracks", 20, 0.0, 10.0);
47 TH1D *htrue = new TH1D("htrue", "True matches", 20, 0.0, 10.0);
48 TH1D *hfound = new TH1D("hfound", "Found matches", 20, 0.0, 10.0);
49 hgood->Sumw2();
50 hfake->Sumw2();
51 htrue->Sumw2();
52 hfound->Sumw2();
53
54 //
55 // Open file containing true matches,
56 // retrieve the Tree and link to a cursor.
57 //
58 TFile *fileTrue = TFile::Open("true-matches.root");
59 match_t trueMatch;
60
61 //
62 // Open file of found matches,
63 // link the modified ESD container.
64 //
65 TFile *fileFound = TFile::Open("matchESD.root");
66 TTree *treeFound = (TTree*)fileFound->Get("esdTree");
67 AliESD *esd = 0;
68 treeFound->SetBranchAddress("ESD", &esd);
69 Long64_t nEvents = treeFound->GetEntries();
70
71 //
72 // Loop on all events
73 //
74 Int_t im, it, ic, nTrueMatches, nTracks;
75 Int_t label, trkLabel, cluLabel;
76 for (Long64_t iev = 0; iev < nEvents; iev++) {
77
78 // get true matches tree of given event
79 TTree *treeTrue = (TTree*)fileTrue->Get(Form("tm_%d", iev));
80 treeTrue->SetBranchAddress("matches", &trueMatch);
81 nTrueMatches = treeTrue->GetEntries();
82
83 // set TTree pointers to selected event
84 runLoader->GetEvent(iev);
85 treeFound->GetEntry(iev);
86 AliStack *stack = runLoader->Stack();
87 nTracks = esd->GetNumberOfTracks();
88
89 // read all true pairs
90 for (im = 0; im < nTrueMatches; im++) {
91 treeTrue->GetEntry(im);
92 AliESDtrack *track = esd->GetTrack(trueMatch.indexT);
93 label = TMath::Abs(track->GetLabel());
94 TParticle *p = stack->Particle(label);
95 htrue->Fill(p->Pt());
96 }
97
98 // compare found matches
99 for (Int_t it = 0; it < nTracks; it++) {
100 AliESDtrack *track = esd->GetTrack(it);
101 ic = track->GetEMCALcluster();
a0270705 102 if (ic == AliEMCALTracker::kUnmatched) continue;
e10611f0 103 ic = TMath::Abs(ic);
104 AliESDCaloCluster *cl = esd->GetCaloCluster(ic);
105 if (!cl) continue;
a0270705 106 trkLabel = TMath::Abs(track->GetLabel());
e10611f0 107 cluLabel = cl->GetPrimaryIndex();
a0270705 108 if (trkLabel == cluLabel && trkLabel >= 0) {
e10611f0 109 TParticle *p = stack->Particle(TMath::Abs(trkLabel));
110 hgood->Fill(p->Pt());
111 hfound->Fill(p->Pt());
112 }
113 else {
114 TParticle *p = stack->Particle(TMath::Abs(trkLabel));
115 hfake->Fill(p->Pt());
116 hfound->Fill(p->Pt());
117 }
118 }
119 }
120
121 cout << "True matches : " << htrue->GetEntries() << endl;
122 cout << "Found matches: " << hfound->GetEntries() << endl;
123 cout << "Good matches : " << hgood->GetEntries() << endl;
124 cout << "Fake matches : " << hfake->GetEntries() << endl;
125
126 TFile *fout = TFile::Open("match-comparison.root", "RECREATE");
127 hgood->Write();
128 hfake->Write();
129 htrue->Write();
130 hfound->Write();
131 fout->Close();
132}