]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/macros/trackMatching/MatchComparison.C
Update tracking class and macros, and pid class to read correctly the ESDs
[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 //
e10611f0 28
29 //
30 // Initialize run loader and load Kinematics
31 //
32 AliRunLoader *runLoader = AliRunLoader::Open("galice.root");
33 if (!runLoader) return;
34 runLoader->LoadgAlice();
35 gAlice = runLoader->GetAliRun();
36 runLoader->LoadKinematics();
37
38 //
39 // Initialize histograms with their error computation
40 //
3a2a23e1 41 TH1D *hgood = new TH1D("hgood", "Well matched tracks", 40, 0.0, 40.0);
42 TH1D *hfake = new TH1D("hfake", "Fake matched tracks", 40, 0.0, 40.0);
43 TH1D *htrue = new TH1D("htrue", "True matches" , 40, 0.0, 40.0);
44 TH1D *hfound = new TH1D("hfound","Found matches" , 40, 0.0, 40.0);
e10611f0 45 hgood->Sumw2();
46 hfake->Sumw2();
47 htrue->Sumw2();
48 hfound->Sumw2();
49
50 //
51 // Open file containing true matches,
52 // retrieve the Tree and link to a cursor.
53 //
54 TFile *fileTrue = TFile::Open("true-matches.root");
55 match_t trueMatch;
56
57 //
58 // Open file of found matches,
59 // link the modified ESD container.
60 //
61 TFile *fileFound = TFile::Open("matchESD.root");
62 TTree *treeFound = (TTree*)fileFound->Get("esdTree");
3a2a23e1 63 AliESDEvent* esd = new AliESDEvent();
64 esd->ReadFromTree(treeFound);
e10611f0 65 Long64_t nEvents = treeFound->GetEntries();
66
67 //
68 // Loop on all events
69 //
70 Int_t im, it, ic, nTrueMatches, nTracks;
71 Int_t label, trkLabel, cluLabel;
72 for (Long64_t iev = 0; iev < nEvents; iev++) {
73
74 // get true matches tree of given event
75 TTree *treeTrue = (TTree*)fileTrue->Get(Form("tm_%d", iev));
76 treeTrue->SetBranchAddress("matches", &trueMatch);
77 nTrueMatches = treeTrue->GetEntries();
78
79 // set TTree pointers to selected event
80 runLoader->GetEvent(iev);
81 treeFound->GetEntry(iev);
82 AliStack *stack = runLoader->Stack();
83 nTracks = esd->GetNumberOfTracks();
84
85 // read all true pairs
86 for (im = 0; im < nTrueMatches; im++) {
87 treeTrue->GetEntry(im);
88 AliESDtrack *track = esd->GetTrack(trueMatch.indexT);
3a2a23e1 89 if (!track) continue;
90
e10611f0 91 label = TMath::Abs(track->GetLabel());
92 TParticle *p = stack->Particle(label);
93 htrue->Fill(p->Pt());
3a2a23e1 94 cout <<"filling true"<< endl;
e10611f0 95 }
96
97 // compare found matches
98 for (Int_t it = 0; it < nTracks; it++) {
99 AliESDtrack *track = esd->GetTrack(it);
100 ic = track->GetEMCALcluster();
a0270705 101 if (ic == AliEMCALTracker::kUnmatched) continue;
e10611f0 102 ic = TMath::Abs(ic);
103 AliESDCaloCluster *cl = esd->GetCaloCluster(ic);
104 if (!cl) continue;
3a2a23e1 105 if (!cl->IsEMCAL()) continue ;
a0270705 106 trkLabel = TMath::Abs(track->GetLabel());
3a2a23e1 107 cluLabel = cl->GetLabel();
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());
3a2a23e1 112 cout <<"filling GOOD, pt:" << p->Pt()<< endl;
e10611f0 113 }
114 else {
115 TParticle *p = stack->Particle(TMath::Abs(trkLabel));
116 hfake->Fill(p->Pt());
117 hfound->Fill(p->Pt());
3a2a23e1 118 cout <<"filling FAKE" << endl;
e10611f0 119 }
120 }
121 }
122
123 cout << "True matches : " << htrue->GetEntries() << endl;
124 cout << "Found matches: " << hfound->GetEntries() << endl;
125 cout << "Good matches : " << hgood->GetEntries() << endl;
126 cout << "Fake matches : " << hfake->GetEntries() << endl;
127
128 TFile *fout = TFile::Open("match-comparison.root", "RECREATE");
129 hgood->Write();
130 hfake->Write();
131 htrue->Write();
132 hfound->Write();
133 fout->Close();
134}