]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFComparison.C
Adapting script of README to new tracking option
[u/mrichter/AliRoot.git] / TOF / AliTOFComparison.C
1 /****************************************************************************
2  *      This macro estimates efficiency of matching with the TOF.           *
3  *      TOF "Good" tracks are those originating from the primary vertex,    *
4  *      being "good" in the ITS and having at least one digit in the TOF.   * 
5  *         (To get the list of "good" tracks one should first run           *
6  *          AliTPCComparison.C and AliITSComparisonV2.C macros)             *
7  *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
8  ****************************************************************************/
9
10 #if !defined(__CINT__) || defined(__MAKECINT__)
11   #include <TMath.h>
12   #include <TError.h>
13   #include <Riostream.h>
14   #include <TH1F.h>
15   #include <TTree.h>
16   #include <TParticle.h>
17   #include <TCanvas.h>
18   #include <TLine.h>
19   #include <TText.h>
20   #include <TBenchmark.h>
21   #include <TStyle.h>
22   #include <TFile.h>
23   #include <TROOT.h>
24
25   #include "AliStack.h"
26   #include "AliHeader.h"
27   #include "AliTrackReference.h"
28   #include "AliRunLoader.h"
29   #include "AliRun.h"
30   #include "AliESD.h"
31
32   #include "AliTOFdigit.h"
33   #include "AliLoader.h"
34
35   #include "TClonesArray.h"
36 #endif
37
38 Int_t GoodTracksTOF(const Char_t *dir=".");
39
40 extern AliRun *gAlice;
41 extern TBenchmark *gBenchmark;
42 extern TROOT *gROOT;
43
44 static Int_t allgood=0;
45 static Int_t allmatched=0;
46 static Int_t allmismatched=0;
47
48 Int_t AliTOFComparison(const Char_t *dir=".") {
49    gBenchmark->Start("AliTOFComparison");
50
51    ::Info("AliTOFComparison.C","Doing comparison...");
52    
53
54    Double_t pmin=0.2;
55    Double_t pmax=3.0;
56
57    TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");    
58    if (!hgood) hgood=new TH1F("hgood","Good tracks",30,pmin,pmax);
59     
60    TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
61    if (!hfound) hfound=new TH1F("hfound","Matched tracks",30,pmin,pmax);
62
63    TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
64    if (!hfake) hfake=new TH1F("hfake","Mismatched tracks",30,pmin,pmax);
65
66    TH1F *hgp=(TH1F*)gROOT->FindObject("hgp");
67    if (!hgp) hgp=new TH1F("hgp","",30,pmin,pmax);
68    hgp->SetLineColor(4); hgp->SetLineWidth(2);
69
70    TH1F *hfp=(TH1F*)gROOT->FindObject("hfp");
71    if (!hfp) hfp=new TH1F("hfp","Probability of mismatching",30,pmin,pmax);
72    hfp->SetFillColor(1); hfp->SetFillStyle(3013); hfp->SetLineWidth(2);
73
74    TH1F *hgoo=(TH1F*)gROOT->FindObject("hgoo");    
75    if (!hgoo) hgoo=new TH1F("hgoo","Good tracks",30,-1,1);
76     
77    TH1F *hfoun=(TH1F*)gROOT->FindObject("hfoun");
78    if (!hfoun) hfoun=new TH1F("hfoun","Matched tracks",30,-1,1);
79
80    TH1F *hfak=(TH1F*)gROOT->FindObject("hfak");
81    if (!hfak) hfak=new TH1F("hfak","Mismatched tracks",30,-1,1);
82
83    TH1F *hgl=(TH1F*)gROOT->FindObject("hgl");
84    if (!hgl) hgl=new TH1F("hgl","",30,-1,1);
85    hgl->SetLineColor(4); hgl->SetLineWidth(2);
86
87    TH1F *hfl=(TH1F*)gROOT->FindObject("hfl");
88    if (!hfl) hfl=new TH1F("hfl","Probability of mismatching",30,-1,1);
89    hfl->SetFillColor(1); hfl->SetFillStyle(3013); hfl->SetLineWidth(2);
90
91
92
93    Char_t fname[100];
94    sprintf(fname,"%s/GoodTracksTOF.root",dir);
95
96    TFile *refFile=TFile::Open(fname,"old");
97    if (!refFile || !refFile->IsOpen()) {
98    ::Info("AliTOFComparison.C","Marking good tracks (will take a while)...");
99      if (GoodTracksTOF(dir)) {
100         ::Error("AliTOFComparison.C","Can't generate the reference file !");
101         return 1;
102      }
103    }
104    refFile=TFile::Open(fname,"old");
105    if (!refFile || !refFile->IsOpen()) {
106      ::Error("AliTOFComparison.C","Can't open the reference file !");
107      return 1;
108    }   
109
110    TTree *tofTree=(TTree*)refFile->Get("tofTree");
111    if (!tofTree) {
112      ::Error("AliTOFComparison.C","Can't get the reference tree !");
113      return 2;
114    }
115    TBranch *branch=tofTree->GetBranch("TOF");
116    if (!branch) {
117      ::Error("AliTOFComparison.C","Can't get the TOF branch !");
118      return 3;
119    }
120    TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
121    branch->SetAddress(&refs);
122
123
124    if (gAlice) { 
125      delete gAlice->GetRunLoader();
126      delete gAlice;//if everything was OK here it is already NULL
127      gAlice = 0x0;
128    }
129    sprintf(fname,"%s/galice.root",dir);
130    AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
131    if (rl == 0x0) {
132       cerr<<"Can not open session"<<endl;
133       return 1;
134    }
135    AliLoader* tofl = rl->GetLoader("TOFLoader");
136    if (tofl == 0x0) {
137       cerr<<"Can not get the TOF loader"<<endl;
138       return 2;
139    }
140    tofl->LoadDigits("read");
141
142
143    sprintf(fname,"%s/AliESDs.root",dir);
144    TFile *ef=TFile::Open(fname);
145    if ((!ef)||(!ef->IsOpen())) {
146       ::Error("AliTOFComparison.C","Can't open AliESDs.root !");
147       delete rl;
148       return 4;
149    }
150    AliESD* event = new AliESD;
151    TTree* esdTree = (TTree*) ef->Get("esdTree");
152    if (!esdTree) {
153       ::Error("AliTOFComparison.C", "no ESD tree found");
154       return 5;
155    }
156    esdTree->SetBranchAddress("ESD", &event);
157
158
159
160    //******* Loop over events *********
161    Int_t e=0;
162    while (esdTree->GetEvent(e)) {
163      cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
164
165      rl->GetEvent(e);
166
167      TTree *digTree=tofl->TreeD();
168      if (!digTree) {
169         cerr<<"Can't get the TOF cluster tree !\n";
170         return 3;
171      } 
172      TBranch *branch=digTree->GetBranch("TOF");
173      if (!branch) { 
174         cerr<<"Can't get the branch with the TOF digits !\n";
175         return 4;
176      }
177      TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;   
178      branch->SetAddress(&digits);
179
180      digTree->GetEvent(0);
181      Int_t nd=digits->GetEntriesFast();
182      cerr<<"Number of the TOF digits: "<<nd<<endl;
183
184
185
186      Int_t ntrk=event->GetNumberOfTracks();
187      cerr<<"Number of ESD tracks : "<<ntrk<<endl; 
188
189
190      if (tofTree->GetEvent(e++)==0) {
191         cerr<<"No reconstructable tracks !\n";
192         continue;
193      }
194
195      Int_t ngood=refs->GetEntriesFast(); 
196
197      Int_t matched=0;
198      Int_t mismatched=0;
199      for (Int_t i=0; i<ngood; i++) {
200         AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(i); 
201         Int_t lab=ref->Label();
202         Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());
203
204         Double_t tgl=ref->Pz()/ptg; //tan(lambda)
205
206          if (ptg>pmin) { hgood->Fill(ptg); hgoo->Fill(tgl); }
207
208          Int_t j;
209          AliESDtrack *t=0;
210          for (j=0; j<ntrk; j++) {
211              AliESDtrack *tt=event->GetTrack(j);
212              if (lab!=TMath::Abs(tt->GetLabel())) continue;
213              t=tt;
214              //if ((tt->GetStatus()&AliESDtrack::kTOFpid) == 0) continue;
215              if (tt->GetTOFsignal() < 0) continue;
216              UInt_t idx=tt->GetTOFcluster();
217              if ((Int_t)idx>=nd) {
218                cerr<<"Wrong digit index ! "<<idx<<endl;
219                return 5;
220              }
221              AliTOFdigit *dig=(AliTOFdigit*)digits->UncheckedAt(idx);
222              Int_t *label=dig->GetTracks();
223              if (label[0]!=lab)
224              if (label[1]!=lab)
225                if (label[2]!=lab) {
226                  mismatched++; 
227                  if (ptg>pmin) { hfake->Fill(ptg); hfak->Fill(tgl); } 
228                  break;
229                }
230              if (ptg>pmin) { hfound->Fill(ptg); hfoun->Fill(tgl); }
231              matched++;
232              break;
233          }
234          if (j==ntrk) {
235             cerr<<"Not matched: "<<lab<<"   ";
236             if (t) {
237                cerr<<(t->GetStatus()&AliESDtrack::kITSout)<<' '
238                    <<(t->GetStatus()&AliESDtrack::kTPCout)<<' '
239                    <<(t->GetStatus()&AliESDtrack::kTRDout)<<' '
240                    <<(t->GetStatus()&AliESDtrack::kTIME);
241             } else cerr<<"No ESD track !";
242             cerr<<endl;
243          }
244      }
245      cout<<"Number of good tracks: "<<ngood<<endl;
246      cout<<"Number of matched tracks: "<<matched<<endl;
247      cout<<"Number of mismatched tracks: "<<mismatched<<endl;
248
249      allgood+=ngood; allmatched+=matched; allmismatched+=mismatched;
250
251      refs->Clear();
252    } //***** End of the loop over events
253
254    delete event;
255    ef->Close();
256    
257    delete tofTree;
258    refFile->Close();
259
260    if (allgood!=0) cerr<<"\n\nEfficiency: "<<Float_t(allmatched)/allgood<<endl;
261    cout<<"Total number of good tracks: "<<allgood<<endl;
262    cout<<"Total number of matched tracks: "<<allmatched<<endl;
263    cout<<"Total number of mismatched tracks: "<<allmismatched<<endl;
264
265    TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
266    if (!c1) {
267       c1=new TCanvas("c1","",0,0,600,900);
268       c1->Divide(1,2);
269    }
270    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
271    hgp->Divide(hfound,hgood,1,1.,"b");
272    hfp->Divide(hfake,hgood,1,1.,"b");
273    hgp->SetMaximum(1.4);
274    hgp->SetYTitle("Matching efficiency");
275    hgp->SetXTitle("Pt (GeV/c)");
276
277    hfoun->Sumw2(); hgoo->Sumw2(); hfak->Sumw2();
278    hgl->Divide(hfoun,hgoo,1,1.,"b");
279    hfl->Divide(hfak,hgoo,1,1.,"b");
280    hgl->SetMaximum(1.4);
281    hgl->SetYTitle("Matching efficiency");
282    hgl->SetXTitle("Tan(lambda)");
283
284    c1->cd(1);
285
286    hgp->Draw();
287    hfp->Draw("histsame");
288    TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
289    line1->Draw("same");
290    TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
291    line2->Draw("same");
292
293    c1->cd(2);
294
295    hgl->Draw();
296    hfl->Draw("histsame");
297    TLine *line3 = new TLine(-1,1.0,1,1.0); line3->SetLineStyle(4);
298    line3->Draw("same");
299    TLine *line4 = new TLine(-1,0.9,1,0.9); line4->SetLineStyle(4);
300    line4->Draw("same");
301
302    c1->Update();
303    
304    TFile fc("AliTOFComparison.root","RECREATE");
305    c1->Write();
306    fc.Close();
307
308    gBenchmark->Stop("AliTOFComparison");
309    gBenchmark->Show("AliTOFComparison");
310
311    delete rl;
312    return 0;
313 }
314
315 Int_t GoodTracksTOF(const Char_t *dir) {
316    if (gAlice) { 
317        delete gAlice->GetRunLoader();
318        delete gAlice;//if everything was OK here it is already NULL
319        gAlice = 0x0;
320    }
321
322    Char_t fname[100];
323    sprintf(fname,"%s/galice.root",dir);
324
325    AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
326    if (!rl) {
327       ::Error("GoodTracksTOF","Can't start session !");
328       return 1;
329    }
330
331    rl->LoadgAlice();
332    rl->LoadHeader();
333    rl->LoadKinematics();
334
335
336    AliLoader* tofl = rl->GetLoader("TOFLoader");
337    if (tofl == 0x0) {
338       ::Error("GoodTracksTOF","Can not get the TOF loader !");
339       delete rl;
340       return 2;
341    }
342    tofl->LoadDigits("read");
343
344    Int_t nev=rl->GetNumberOfEvents();
345    ::Info("GoodTracksTOF","Number of events : %d\n",nev);  
346
347    sprintf(fname,"%s/GoodTracksITS.root",dir);
348    TFile *itsFile=TFile::Open(fname);
349    if ((!itsFile)||(!itsFile->IsOpen())) {
350        ::Error("GoodTracksTOF","Can't open the GoodTracksITS.root !");
351        delete rl;
352        return 5; 
353    }
354    TClonesArray dum("AliTrackReference",1000), *itsRefs=&dum;
355    TTree *itsTree=(TTree*)itsFile->Get("itsTree");
356    if (!itsTree) {
357        ::Error("GoodTracksTOF","Can't get the ITS reference tree !");
358        delete rl;
359        return 6;
360    }
361    TBranch *itsBranch=itsTree->GetBranch("ITS");
362    if (!itsBranch) {
363       ::Error("GoodTracksTOF","Can't get the ITS reference branch !");
364       delete rl;
365       return 7;
366    }
367    itsBranch->SetAddress(&itsRefs);
368
369
370    sprintf(fname,"%s/GoodTracksTOF.root",dir);
371    TFile *tofFile=TFile::Open(fname,"recreate");
372    TClonesArray dummy("AliTrackReference",1000), *tofRefs=&dummy;
373    TTree tofTree("tofTree","Tree with info about the reconstructable TOF tracks");
374    tofTree.Branch("TOF",&tofRefs);
375
376
377    //********  Loop over generated events 
378    for (Int_t e=0; e<nev; e++) {
379
380      rl->GetEvent(e); tofFile->cd();
381
382      Int_t np = rl->GetHeader()->GetNtrack();
383      cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
384
385      //******** Fill the "good" masks
386      Int_t *good=new Int_t[np]; Int_t k; for (k=0; k<np; k++) good[k]=0;
387
388      TTree *dTree=tofl->TreeD();
389      if (!dTree) {
390         ::Error("GoodTracksTOF","Can't get the TOF cluster tree !");
391         delete rl;
392         return 8;
393      } 
394      TBranch *branch=dTree->GetBranch("TOF");
395      if (!branch) { 
396         ::Error("GoodTracksTOF","Can't get the branch with the TOF digits !");
397         return 9;
398      }
399      TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
400      branch->SetAddress(&digits);
401    
402      dTree->GetEvent(0);
403      Int_t nd=digits->GetEntriesFast();
404
405      for (Int_t i=0; i<nd; i++) {
406        AliTOFdigit *d=(AliTOFdigit*)digits->UncheckedAt(i);
407        Int_t l0=d->GetTrack(0);
408           if (l0>=np) {cerr<<"Wrong label: "<<l0<<endl; continue;}
409        Int_t l1=d->GetTrack(1);
410           if (l1>=np) {cerr<<"Wrong label: "<<l1<<endl; continue;}
411        Int_t l2=d->GetTrack(2);
412           if (l2>=np) {cerr<<"Wrong label: "<<l2<<endl; continue;}
413        if (l0>=0) good[l0]++; 
414        if (l1>=0) good[l1]++; 
415        if (l2>=0) good[l2]++;
416      }
417      digits->Clear();
418
419
420      //****** select tracks which are "good" enough
421      AliStack* stack = rl->Stack();
422
423      itsTree->GetEvent(e);
424      Int_t nk=itsRefs->GetEntriesFast();
425
426      Int_t nt=0;
427      for (k=0; k<nk; k++) {
428         AliTrackReference *itsRef=(AliTrackReference *)itsRefs->UncheckedAt(k);
429         Int_t lab=itsRef->Label();
430         if (good[lab] == 0) continue;
431         TParticle *p = (TParticle*)stack->Particle(lab);
432         if (p == 0x0) {
433            cerr<<"Can not get particle "<<lab<<endl;
434            continue;
435         }
436
437         if (TMath::Abs(p->Vx())>0.1) continue;
438         if (TMath::Abs(p->Vy())>0.1) continue;
439         //if (TMath::Abs(p->Vz())>0.1) continue;
440
441         new((*tofRefs)[nt]) AliTrackReference(*itsRef);
442         nt++;
443      }
444      itsRefs->Clear();
445
446      tofTree.Fill();
447      tofRefs->Clear();
448
449      delete[] good;
450
451    } //*** end of the loop over generated events
452
453    tofTree.Write();
454    tofFile->Close();
455
456    delete itsTree;
457    itsFile->Close();
458
459    delete rl;
460    return 0;
461 }