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