]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSComparisonV2.C
Missing copy constructor and assignment operator are added
[u/mrichter/AliRoot.git] / ITS / AliITSComparisonV2.C
1 /****************************************************************************
2  *           Very important, delicate and rather obscure macro.             *
3  *                                                                          *
4  *               Creates list of "trackable" tracks,                        *
5  *             calculates efficiency, resolutions etc.                      *
6  *  The ESD tracks must be in an appropriate state (kITSin or kITSrefit)    *
7  *                                                                          *
8  *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
9  ****************************************************************************/
10
11 #if !defined(__CINT__) || defined(__MAKECINT__)
12   #include <TMath.h>
13   #include <TError.h>
14   #include <Riostream.h>
15   #include <TH1F.h>
16   #include <TH2F.h>
17   #include <TTree.h>
18   #include <TParticle.h>
19   #include <TCanvas.h>
20   #include <TLine.h>
21   #include <TText.h>
22   #include <TBenchmark.h>
23   #include <TStyle.h>
24   #include <TFile.h>
25   #include <TROOT.h>
26
27   #include "AliStack.h"
28   #include "AliHeader.h"
29   #include "AliTrackReference.h"
30   #include "AliRunLoader.h"
31   #include "AliRun.h"
32   #include "AliESD.h"
33
34   #include "AliITS.h"
35   #include "AliITSgeom.h"
36   #include "AliITSRecPoint.h"
37   #include "AliITSLoader.h"
38 #endif
39
40 Int_t GoodTracksITS(const Char_t *dir=".");
41
42 extern AliRun *gAlice;
43 extern TBenchmark *gBenchmark;
44 extern TROOT *gROOT;
45
46 static Int_t allgood=0;
47 static Int_t allselected=0;
48 static Int_t allfound=0;
49
50 Int_t AliITSComparisonV2
51 (Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".", Float_t ratio=0.0) {
52    gBenchmark->Start("AliITSComparisonV2");
53
54    ::Info("AliITSComparisonV2.C","Doing comparison...");
55    
56
57    TH1F *hp=(TH1F*)gROOT->FindObject("hp");
58    if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.); 
59    hp->SetFillColor(4);
60
61    TH1F *hl=(TH1F*)gROOT->FindObject("hl");
62    if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
63    hl->SetFillColor(4);
64
65    TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
66    if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
67    hpt->SetFillColor(2);
68  
69    TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
70    if (!hmpt) 
71       hmpt=new TH1F("hmpt","Transverse impact parameter",30,-777,777); 
72    hmpt->SetFillColor(6);
73
74    TH1F *hz=(TH1F*)gROOT->FindObject("hz");
75    if (!hz) hz=new TH1F("hz","Longitudinal impact parameter",30,-777,777); 
76
77
78
79    TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
80    if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
81     
82    TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
83    if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);
84
85    TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
86    if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);
87
88    TH1F *hg=(TH1F*)gROOT->FindObject("hg");
89    if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
90    hg->SetLineColor(4); hg->SetLineWidth(2);
91
92    TH1F *hf=(TH1F*)gROOT->FindObject("hf");
93    if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
94    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
95
96    TH1F *he=(TH1F*)gROOT->FindObject("he");
97    if (!he) 
98       he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
99
100    TH2F *hep=(TH2F*)gROOT->FindObject("hep");
101    if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
102    hep->SetMarkerStyle(8);
103    hep->SetMarkerSize(0.4);
104
105
106    Char_t fname[100];
107    sprintf(fname,"%s/GoodTracksITS.root",dir);
108
109    TFile *refFile=TFile::Open(fname,"old");
110    if (!refFile || !refFile->IsOpen()) {
111    ::Info("AliITSComparisonV2.C","Marking good tracks (will take a while)...");
112      if (GoodTracksITS(dir)) {
113         ::Error("AliITSComparisonV2.C","Can't generate the reference file !");
114         return 1;
115      }
116    }
117    refFile=TFile::Open(fname,"old");
118    if (!refFile || !refFile->IsOpen()) {
119      ::Error("AliITSComparisonV2.C","Can't open the reference file !");
120      return 1;
121    }   
122   
123    TTree *itsTree=(TTree*)refFile->Get("itsTree");
124    if (!itsTree) {
125      ::Error("AliITSComparisonV2.C","Can't get the reference tree !");
126      return 2;
127    }
128    TBranch *branch=itsTree->GetBranch("ITS");
129    if (!branch) {
130      ::Error("AliITSComparisonV2.C","Can't get the ITS branch !");
131      return 3;
132    }
133    TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
134    branch->SetAddress(&refs);
135
136
137    sprintf(fname,"%s/AliESDs.root",dir);
138    TFile *ef=TFile::Open(fname);
139    if ((!ef)||(!ef->IsOpen())) {
140       sprintf(fname,"%s/AliESDits.root",dir);
141       ef=TFile::Open(fname);
142       if ((!ef)||(!ef->IsOpen())) {
143          ::Error("AliITSComparisonV2.C","Can't open AliESDits.root !");
144          return 4;
145       }
146    }
147    AliESD* event = new AliESD;
148    TTree* esdTree = (TTree*) ef->Get("esdTree");
149    if (!esdTree) {
150       ::Error("AliITSComparison.C", "no ESD tree found");
151       return 6;
152    }
153    esdTree->SetBranchAddress("ESD", &event);
154
155
156    //******* Loop over events *********
157    Int_t e=0;
158    while (esdTree->GetEvent(e)) {
159      cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
160  
161      Int_t nentr=event->GetNumberOfTracks();
162      allfound+=nentr;
163
164      if (itsTree->GetEvent(e++)==0) {
165         cerr<<"No reconstructable tracks !\n";
166         continue;
167      }
168
169      Int_t ngood=refs->GetEntriesFast(); 
170      allgood+=ngood;
171
172      const Int_t MAX=15000;
173      Int_t notf[MAX], nnotf=0;
174      Int_t fake[MAX], nfake=0;
175      Int_t mult[MAX], numb[MAX], nmult=0;
176      Int_t k;
177      for (k=0; k<ngood; k++) {
178         AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k); 
179         Int_t lab=ref->Label(), tlab=-1;
180         Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());
181
182         if (ptg<ptcutl) continue;
183         if (ptg>ptcuth) continue;
184
185         allselected++;
186
187         hgood->Fill(ptg);
188
189         AliESDtrack *esd=0;
190         Int_t cnt=0;
191         for (Int_t i=0; i<nentr; i++) {
192            AliESDtrack *t=event->GetTrack(i);
193            UInt_t status=t->GetStatus();
194
195            if ((status&AliESDtrack::kITSrefit)==0) continue;
196
197            Int_t lbl=t->GetLabel();
198            if (lab==TMath::Abs(lbl)) {
199               if (cnt==0) {esd=t; tlab=lbl;}
200               cnt++;
201            }
202         }
203         if (cnt==0) {
204            notf[nnotf++]=lab;
205            continue;
206         } else if (cnt>1){
207            mult[nmult]=lab;
208            numb[nmult]=cnt; nmult++;        
209         }
210
211         if (esd->GetITSFakeRatio()<ratio) tlab = TMath::Abs(tlab);
212
213         if (lab==tlab) hfound->Fill(ptg);
214         else {
215           fake[nfake++]=lab;
216           hfake->Fill(ptg); 
217         }
218
219         Double_t alpha=esd->GetAlpha(),xv,par[5]; 
220         esd->GetExternalParameters(xv,par);
221         Float_t phi=TMath::ASin(par[2]) + alpha;
222         if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
223         if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
224         Float_t lam=TMath::ATan(par[3]); 
225         Float_t pt_1=TMath::Abs(par[4]);
226
227         Float_t phig=TMath::ATan2(ref->Py(),ref->Px());
228         hp->Fill((phi - phig)*1000.);
229
230         Float_t lamg=TMath::ATan2(ref->Pz(),ptg);
231         hl->Fill((lam - lamg)*1000.);
232
233         Float_t d,z; esd->GetImpactParameters(d,z);
234         hmpt->Fill(10000*d);
235         hz->Fill(10000*z);
236
237         hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
238
239         Float_t mom=1./(pt_1*TMath::Cos(lam));
240         Float_t dedx=esd->GetITSsignal();
241         hep->Fill(mom,dedx,1.);
242
243         Int_t pdg=(Int_t)ref->GetLength();  //this is particle's PDG !
244
245         if (TMath::Abs(pdg)==211) //pions
246            if (mom>0.4 && mom<0.5) he->Fill(dedx,1.);
247
248      }
249
250      cout<<"\nList of Not found tracks :\n";
251      for (k=0; k<nnotf; k++){
252        cout<<notf[k]<<"\t";
253        if ((k%9)==8) cout<<"\n";
254      }
255      cout<<"\n\nList of fake  tracks :\n";
256      for (k=0; k<nfake; k++){
257        cout<<fake[k]<<"\t";
258        if ((k%9)==8) cout<<"\n";
259      }
260      cout<<"\n\nList of multiple found tracks :\n";
261      for (k=0; k<nmult; k++) {
262          cout<<"id.   "<<mult[k]
263              <<"     found - "<<numb[k]<<"times\n";
264      }
265      cout<<endl;
266
267      cout<<"Number of found tracks : "<<nentr<<endl;
268      cout<<"Number of \"good\" tracks : "<<ngood<<endl;
269
270      refs->Clear();
271    } //***** End of the loop over events
272
273    delete event;
274    ef->Close();
275    
276    delete itsTree;
277    refFile->Close();
278
279    Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
280    if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
281    cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl;
282    cout<<"Total number of found tracks ="<<allfound<<endl;
283    cout<<"Total number of \"good\" tracks ="<<allgood<<endl;
284    cout<<endl;
285
286    gStyle->SetOptStat(111110);
287    gStyle->SetOptFit(1);
288
289    TCanvas *c1=new TCanvas("c1","",0,0,700,850);
290
291    Int_t minc=33; 
292
293    TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
294    p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); 
295    hp->SetFillColor(4);  hp->SetXTitle("(mrad)"); 
296    if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();
297
298    TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); 
299    p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
300    hl->SetXTitle("(mrad)");
301    if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();
302
303    TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
304    p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); 
305    hpt->SetXTitle("(%)");
306    if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();
307
308    TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
309    p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
310    hmpt->SetXTitle("(micron)");
311    if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); 
312    hz->Draw("same"); c1->cd();
313    
314
315    TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); 
316    p5->SetFillColor(41); p5->SetFrameFillColor(10);
317    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
318    hg->Divide(hfound,hgood,1,1.,"b");
319    hf->Divide(hfake,hgood,1,1.,"b");
320    hg->SetMaximum(1.4);
321    hg->SetYTitle("Tracking efficiency");
322    hg->SetXTitle("Pt (GeV/c)");
323    hg->Draw();
324
325    TLine *line1 = new TLine(0.2,1.0,6.1,1.0); line1->SetLineStyle(4);
326    line1->Draw("same");
327    TLine *line2 = new TLine(0.2,0.9,6.1,0.9); line2->SetLineStyle(4);
328    line2->Draw("same");
329
330    hf->SetFillColor(1);
331    hf->SetFillStyle(3013);
332    hf->SetLineColor(2);
333    hf->SetLineWidth(2);
334    hf->Draw("histsame");
335    TText *text = new TText(0.461176,0.248448,"Fake tracks");
336    text->SetTextSize(0.05);
337    text->Draw();
338    text = new TText(0.453919,1.11408,"Good tracks");
339    text->SetTextSize(0.05);
340    text->Draw();
341
342    TCanvas *c2=new TCanvas("c2","",320,32,530,590);
343    TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
344    p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); 
345    he->SetFillColor(2); he->SetFillStyle(3005);  
346    he->SetXTitle("Arbitrary Units"); 
347    if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();
348
349    TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); 
350    p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
351    hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); 
352    hep->Draw(); c1->cd();
353
354    TFile fc("AliITSComparisonV2.root","RECREATE");
355    c1->Write();
356    c2->Write();
357    fc.Close();
358
359    gBenchmark->Stop("AliITSComparisonV2");
360    gBenchmark->Show("AliITSComparisonV2");
361
362    return 0;
363 }
364
365
366
367 Int_t GoodTracksITS(const Char_t *dir) {
368    if (gAlice) { 
369        delete gAlice->GetRunLoader();
370        delete gAlice;//if everything was OK here it is already NULL
371        gAlice = 0x0;
372    }
373
374    Char_t fname[100];
375    sprintf(fname,"%s/galice.root",dir);
376
377    AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
378    if (!rl) {
379       ::Error("GoodTracksITS","Can't start session !");
380       return 1;
381    }
382
383    rl->LoadgAlice();
384    rl->LoadHeader();
385    rl->LoadKinematics();
386
387    AliITS *ITS=(AliITS*)rl->GetAliRun()->GetDetector("ITS");
388    if (!ITS) {
389       ::Error("GoodTracksITS","Can't get the ITS !");
390       delete rl;
391       return 2;
392    }
393    AliITSgeom *geom=ITS->GetITSgeom();
394    if (!geom) {
395       ::Error("GoodTracksITS","Can't get the ITS geometry !"); 
396       delete rl;
397       return 3;
398    }
399
400    AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
401    if (itsl == 0x0) {
402        ::Error("GoodTracksITS","Can not find the ITSLoader");
403        delete rl;
404        return 4;
405    }
406    itsl->LoadRecPoints();
407   
408
409    Int_t nev=rl->GetNumberOfEvents();
410    ::Info("GoodTracksITS","Number of events : %d\n",nev);  
411
412    sprintf(fname,"%s/GoodTracksTPC.root",dir);
413    TFile *tpcFile=TFile::Open(fname);
414    if ((!tpcFile)||(!tpcFile->IsOpen())) {
415        ::Error("GoodTracksITS","Can't open the GoodTracksTPC.root !");
416        delete rl;
417        return 5; 
418    }
419    TClonesArray dum("AliTrackReference",1000), *tpcRefs=&dum;
420    TTree *tpcTree=(TTree*)tpcFile->Get("tpcTree");
421    if (!tpcTree) {
422        ::Error("GoodTracksITS","Can't get the TPC reference tree !");
423        delete rl;
424        return 6;
425    }
426    TBranch *tpcBranch=tpcTree->GetBranch("TPC");
427    if (!tpcBranch) {
428       ::Error("GoodTracksITS","Can't get the TPC reference branch !");
429       delete rl;
430       return 7;
431    }
432    tpcBranch->SetAddress(&tpcRefs);
433
434    sprintf(fname,"%s/GoodTracksITS.root",dir);
435    TFile *itsFile=TFile::Open(fname,"recreate");
436    TClonesArray dummy2("AliTrackReference",1000), *itsRefs=&dummy2;
437    TTree itsTree("itsTree","Tree with info about the reconstructable ITS tracks");
438    itsTree.Branch("ITS",&itsRefs);
439
440    //********  Loop over generated events 
441    for (Int_t e=0; e<nev; e++) {
442      Int_t k;
443
444      rl->GetEvent(e);  itsFile->cd();
445
446      Int_t np = rl->GetHeader()->GetNtrack();
447      cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
448
449      //******** Fill the "good" masks
450      Int_t *good=new Int_t[np]; for (k=0; k<np; k++) good[k]=0;
451
452      TTree *cTree=itsl->TreeR();
453      if (!cTree) {
454         ::Error("GoodTracksITS","Can't get the cluster tree !"); 
455         delete rl;
456         return 8;
457      }
458      TBranch *branch=cTree->GetBranch("ITSRecPoints");
459      if (!branch) {
460         ::Error("GoodTracksITS","Can't get the clusters branch !"); 
461         delete rl;
462         return 9;
463      }
464      TClonesArray dummy("AliITSRecPoint",10000), *clusters=&dummy;
465      branch->SetAddress(&clusters);
466
467      Int_t entr=(Int_t)cTree->GetEntries();
468      for (k=0; k<entr; k++) {
469          cTree->GetEvent(k);
470          Int_t ncl=clusters->GetEntriesFast(); if (ncl==0) continue;
471          Int_t lay,lad,det;  geom->GetModuleId(k,lay,lad,det);
472          if (lay<1 || lay>6) {
473             ::Error("GoodTracksITS","Wrong layer !"); 
474             delete rl;
475             return 10;
476          }
477          while (ncl--) {
478             AliITSRecPoint *pnt=(AliITSRecPoint*)clusters->UncheckedAt(ncl);
479             Int_t l0=pnt->GetLabel(0);
480                if (l0>=np) {
481 //               cerr<<"Wrong label: "<<l0<<endl;
482                  continue;
483                }
484             Int_t l1=pnt->GetLabel(1);
485                if (l1>=np) {
486 //               cerr<<"Wrong label: "<<l1<<endl;
487                  continue;
488                }
489             Int_t l2=pnt->GetLabel(2);
490                if (l2>=np) {
491 //               cerr<<"Wrong label: "<<l2<<endl;
492                  continue;
493                }
494             Int_t mask=1<<(lay-1);
495             if (l0>=0) good[l0]|=mask; 
496             if (l1>=0) good[l1]|=mask; 
497             if (l2>=0) good[l2]|=mask;
498          }
499          clusters->Clear();
500      }
501    
502
503      //****** select tracks which are "good" enough
504      AliStack* stack = rl->Stack();
505
506      tpcTree->GetEvent(e);
507      Int_t nk=tpcRefs->GetEntriesFast();
508      Int_t nt=0;
509      for (k=0; k<nk; k++) {
510         AliTrackReference *tpcRef=(AliTrackReference *)tpcRefs->UncheckedAt(k);
511         Int_t lab=tpcRef->Label();
512         if (good[lab] != 0x3F) continue;
513         TParticle *p = (TParticle*)stack->Particle(lab);
514         if (p == 0x0) {
515            cerr<<"Can not get particle "<<lab<<endl;
516            continue;
517         }
518
519         AliTrackReference *ref=new((*itsRefs)[nt]) AliTrackReference(*tpcRef);
520         ref->SetMomentum(p->Px(),p->Py(),p->Pz());
521         ref->SetPosition(p->Vx(),p->Vy(),p->Vz());
522         nt++;
523      }
524      tpcRefs->Clear();
525
526      itsTree.Fill();
527      itsRefs->Clear();
528
529      delete[] good;
530
531    } //*** end of the loop over generated events
532
533    itsTree.Write();
534    itsFile->Close();
535
536    delete tpcTree;
537    tpcFile->Close();
538
539    delete rl;
540    return 0;
541 }
542
543