]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCComparison.C
Updated version (Yu.Belikov)
[u/mrichter/AliRoot.git] / TPC / AliTPCComparison.C
1 /****************************************************************************
2  *           Very important, delicate and rather obscure macro.             *
3  *                                                                          *
4  *               Creates list of "trackable" tracks,                        *
5  *             calculates efficiency, resolutions etc.                      *
6  *                                                                          *
7  *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
8  * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de        *
9  ****************************************************************************/
10
11 #if !defined(__CINT__) || defined(__MAKECINT__)
12   #include "alles.h"
13   #include "AliTPCtracker.h"
14   #include "AliRunLoader.h"
15   #include "AliTPCLoader.h"
16   #include "AliMagF.h"
17   #include "AliRun.h"
18 #endif
19
20 struct GoodTrackTPC {
21   Int_t lab;
22   Int_t code;
23   Float_t px,py,pz;
24   Float_t x,y,z;
25 };
26
27 Int_t 
28 good_tracks_tpc(GoodTrackTPC *gt, const Int_t max, const char* evfoldname);
29
30 extern AliRun *gAlice;
31
32 Int_t AliTPCComparison(Bool_t thcut=1.0) {
33
34    if (gAlice) { 
35        delete gAlice->GetRunLoader();
36        delete gAlice;//if everything was OK here it is already NULL
37        gAlice = 0x0;
38     }
39
40    cerr<<"Doing comparison...\n";
41
42    const Int_t MAX=20000;
43    Int_t good_tracks_tpc(
44       GoodTrackTPC *gt, const Int_t max, 
45       const char* evfoldname = AliConfig::fgkDefaultEventFolderName
46    );//declaration only
47
48    gBenchmark->Start("AliTPCComparison");
49
50    AliRunLoader *rl = AliRunLoader::Open("galice.root","COMPARISON");
51    if (!rl) {
52        cerr<<"Can't start sesion !\n";
53        return 1;
54    }
55
56    rl->LoadgAlice();
57    if (rl->GetAliRun())
58    AliKalmanTrack::SetConvConst(
59      1000/0.299792458/rl->GetAliRun()->Field()->SolenoidField()
60    );
61    else {
62        cerr<<"AliTPCComparison.C :Can't get AliRun !\n";
63        return 1;
64    }
65    //rl->UnloadgAlice();
66   
67    AliTPCLoader * tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
68    if (tpcl == 0x0) {
69        cerr<<"AliTPCComparison.C : Can not find TPCLoader\n";
70        delete rl;
71        return 3;
72    }
73
74    /* Generate a list of "good" tracks */
75
76    GoodTrackTPC gt[MAX];
77    Int_t ngood=0;
78    ifstream in("good_tracks_tpc");
79    if (in) {
80       cerr<<"Reading good tracks...\n";
81       while (in>>gt[ngood].lab>>gt[ngood].code>>
82                  gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
83                  gt[ngood].x >>gt[ngood].y >>gt[ngood].z) {
84         Double_t rin = TMath::Sqrt( gt[ngood].x*gt[ngood].x+gt[ngood].y*gt[ngood].y);
85         if (rin<1) continue;
86         Double_t theta = gt[ngood].z/rin;
87         //theta =0;
88         if (TMath::Abs(theta)<thcut){ 
89           ngood++;
90           cerr<<ngood<<'\r';
91           if (ngood==MAX) {
92             cerr<<"Too many good tracks !\n";
93             break;
94           }     
95         //ngood++;
96         // cerr<<ngood<<'\r';
97         // if (ngood==MAX) {
98         //    cerr<<"Too many good tracks !\n";
99         //    break;
100         }
101       }
102       if (!in.eof()) cerr<<"Read error (good_tracks_tpc) !\n";
103    } else {
104      cerr<<"Marking good tracks (this will take a while)...\n";
105       ngood=good_tracks_tpc(gt,MAX,"COMPARISON");
106       ofstream out("good_tracks_tpc");
107       if (out) {
108          for (Int_t ngd=0; ngd<ngood; ngd++)            
109             out<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
110                  gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
111                  gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
112       } else cerr<<"Can not open file (good_tracks_tpc) !\n";
113       out.close();
114    }
115
116    Int_t nentr=0,i=0; TObjArray tarray(MAX);
117    { /*Load tracks*/
118    
119      tpcl->LoadTracks();
120      
121      TTree *tracktree=tpcl->TreeT();
122      if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n"; return 4;}
123
124      TBranch *tbranch=tracktree->GetBranch("tracks");
125      nentr=(Int_t)tracktree->GetEntries();
126      AliTPCtrack *iotrack=0;
127      for (i=0; i<nentr; i++) {
128        iotrack=new AliTPCtrack;
129        tbranch->SetAddress(&iotrack);
130        tracktree->GetEvent(i);
131        tarray.AddLast(iotrack);
132      }   
133      tpcl->UnloadTracks();
134    }
135
136    TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4);
137    TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4);
138    TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
139    hpt->SetFillColor(2); 
140    TH1F *hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60); 
141    hmpt->SetFillColor(6);
142
143    TH1F *hgood=new TH1F("hgood","Good tracks",30,0.1,6.1);    
144    TH1F *hfound=new TH1F("hfound","Found tracks",30,0.1,6.1);
145    TH1F *hfake=new TH1F("hfake","Fake tracks",30,0.1,6.1);
146    TH1F *hg=new TH1F("hg","",30,0.1,6.1); //efficiency for good tracks
147    hg->SetLineColor(4); hg->SetLineWidth(2);
148    TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,0.1,6.1);
149    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
150
151    TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
152    TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
153    hep->SetMarkerStyle(8);
154    hep->SetMarkerSize(0.4);
155
156    //MI change
157    Int_t mingood=ngood;
158    Int_t track_notfound[MAX], itrack_notfound=0;
159    Int_t track_fake[MAX], itrack_fake=0;
160    Int_t track_multifound[MAX], track_multifound_n[MAX], itrack_multifound=0;
161
162    while (ngood--) {
163       Int_t lab=gt[ngood].lab,tlab=-1;
164       Float_t ptg=
165       TMath::Sqrt(gt[ngood].px*gt[ngood].px + gt[ngood].py*gt[ngood].py);
166
167       if (ptg<1e-33) continue; // for those not crossing 0 pad row
168
169       hgood->Fill(ptg);
170
171       AliTPCtrack *track=0;
172       for (i=0; i<nentr; i++) {
173           track=(AliTPCtrack*)tarray.UncheckedAt(i);
174           tlab=track->GetLabel();
175           if (lab==TMath::Abs(tlab)) break;
176       }
177       if (i==nentr) {
178         track_notfound[itrack_notfound++]=lab;
179         continue;
180       }
181       
182       //MI change  - addition
183       Int_t micount=0;
184       Int_t mi;
185       AliTPCtrack * mitrack;
186       for (mi=0; mi<nentr; mi++) {
187         mitrack=(AliTPCtrack*)tarray.UncheckedAt(mi);          
188         if (lab==TMath::Abs(mitrack->GetLabel())) micount++;
189       }
190       if (micount>1) {
191         track_multifound[itrack_multifound]=lab;
192         track_multifound_n[itrack_multifound]=micount;
193         itrack_multifound++;
194       }
195       
196       Double_t xk=gt[ngood].x;
197       track->PropagateTo(xk);
198
199       if (lab==tlab) hfound->Fill(ptg);
200       else { 
201         track_fake[itrack_fake++]=lab;
202         hfake->Fill(ptg); 
203       }
204
205       Double_t par[5]; track->GetExternalParameters(xk,par);
206       Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
207       if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
208       if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
209       Float_t lam=TMath::ATan(par[3]); 
210       Float_t pt_1=TMath::Abs(par[4]);
211
212       if (TMath::Abs(gt[ngood].code)==11 && ptg>4.) {
213          hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
214       } else {
215          Float_t phig=TMath::ATan2(gt[ngood].py,gt[ngood].px);
216          hp->Fill((phi - phig)*1000.);
217
218          Float_t lamg=TMath::ATan2(gt[ngood].pz,ptg);
219          hl->Fill((lam - lamg)*1000.);
220
221          hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
222       }
223
224       Float_t mom=1./(pt_1*TMath::Cos(lam));
225       Float_t dedx=track->GetdEdx();
226       hep->Fill(mom,dedx,1.);
227       if (TMath::Abs(gt[ngood].code)==211)
228          if (mom>0.4 && mom<0.5) {
229             he->Fill(dedx,1.);
230          }
231
232    }
233
234    cout<<"\nList of Not found tracks :\n";
235    for ( i = 0; i< itrack_notfound; i++){
236      cout<<track_notfound[i]<<"\t";
237      if ((i%5)==4) cout<<"\n";
238    }
239    cout<<"\nList of fake  tracks :\n";
240    for ( i = 0; i< itrack_fake; i++){
241      cout<<track_fake[i]<<"\t";
242      if ((i%5)==4) cout<<"\n";
243    }
244     cout<<"\nList of multiple found tracks :\n";
245    for ( i=0; i<itrack_multifound; i++) {
246        cout<<"id.   "<<track_multifound[i]
247            <<"     found - "<<track_multifound_n[i]<<"times\n";
248    }
249
250    Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
251    if (ng!=0) cerr<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
252    cout<<"Total number of found tracks ="<<nentr<<endl;
253    cout<<"Total number of \"good\" tracks ="
254        <<mingood<<"   (selected for comparison: "<<ng<<')'<<endl<<endl;
255
256
257    gStyle->SetOptStat(111110);
258    gStyle->SetOptFit(1);
259
260    TCanvas *c1=new TCanvas("c1","",0,0,700,850);
261
262    TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
263    p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); 
264    hp->SetFillColor(4);  hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd();
265
266    TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); 
267    p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
268    hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd();
269
270    TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
271    p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); 
272    hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd();
273
274    TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
275    p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
276    hmpt->SetXTitle("(%)"); hmpt->Fit("gaus"); c1->cd();
277
278    TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); 
279    p5->SetFillColor(41); p5->SetFrameFillColor(10);
280    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
281    hg->Divide(hfound,hgood,1,1.,"b");
282    hf->Divide(hfake,hgood,1,1.,"b");
283    hg->SetMaximum(1.4);
284    hg->SetYTitle("Tracking efficiency");
285    hg->SetXTitle("Pt (GeV/c)");
286    hg->Draw();
287
288    TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
289    line1->Draw("same");
290    TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
291    line2->Draw("same");
292
293    hf->SetFillColor(1);
294    hf->SetFillStyle(3013);
295    hf->SetLineColor(2);
296    hf->SetLineWidth(2);
297    hf->Draw("histsame");
298    TText *text = new TText(0.461176,0.248448,"Fake tracks");
299    text->SetTextSize(0.05);
300    text->Draw();
301    text = new TText(0.453919,1.11408,"Good tracks");
302    text->SetTextSize(0.05);
303    text->Draw();
304
305
306
307    TCanvas *c2=new TCanvas("c2","",320,32,530,590);
308
309    TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
310    p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); 
311    he->SetFillColor(2); he->SetFillStyle(3005);  
312    he->SetXTitle("Arbitrary Units"); 
313    he->Fit("gaus"); c2->cd();
314
315    TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); 
316    p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
317    hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); 
318    hep->Draw(); c1->cd();
319
320    gBenchmark->Stop("AliTPCComparison");
321    gBenchmark->Show("AliTPCComparison");
322
323    delete rl;
324    return 0;
325 }
326
327
328 Int_t 
329 good_tracks_tpc(GoodTrackTPC *gt, const Int_t max, const char* evfoldname) {
330    Int_t nt=0;
331
332    AliRunLoader *rl = AliRunLoader::GetRunLoader(evfoldname);
333    if (rl == 0x0) {
334       ::Fatal("AliTPCComparison.C::good_tracks_tpc",
335               "Can not find Run Loader in Folder Named %s",
336               evfoldname);
337    }
338    AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
339    if (tpcl == 0x0) {
340        cerr<<"AliTPCHits2Digits.C : Can not find TPCLoader\n";
341        delete rl;
342        return 0;
343    }
344    
345    rl->LoadgAlice();
346    
347    AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
348    Int_t ver = TPC->IsVersion(); 
349    cerr<<"TPC version "<<ver<<" has been found !\n";
350
351    rl->CdGAFile();
352    AliTPCParamSR *digp=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
353    if (!digp) { cerr<<"TPC parameters have not been found !\n"; exit(6); }
354    TPC->SetParam(digp);
355
356    rl->LoadHeader();
357
358    Int_t np = rl->GetHeader()->GetNtrack();
359       
360
361    Int_t nrow_up=digp->GetNRowUp();
362    Int_t nrows=digp->GetNRowLow()+nrow_up;
363    Int_t zero=digp->GetZeroSup();
364    Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
365    Int_t good_number=Int_t(0.4*nrows);
366
367    Int_t *good=new Int_t[np];
368    Int_t i;
369    for (i=0; i<np; i++) good[i]=0;
370
371
372    switch (ver) {
373    case 1:
374      {
375       tpcl->LoadRecPoints();      
376       AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
377       ca.Setup(digp);
378       ca.SetClusterType("AliTPCcluster");
379       ca.ConnectTree(tpcl->TreeR());
380       Int_t nrows=Int_t(ca.GetTree()->GetEntries());
381       for (Int_t n=0; n<nrows; n++) {
382           AliSegmentID *s=ca.LoadEntry(n);
383           Int_t sec,row;
384           digp->AdjustSectorRow(s->GetID(),sec,row);
385           AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
386           Int_t ncl=clrow.GetArray()->GetEntriesFast();
387           while (ncl--) {
388               AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
389               Int_t lab=c->GetLabel(0);
390               if (lab<0) continue; //noise cluster
391               lab=TMath::Abs(lab);
392
393               if (sec>=digp->GetNInnerSector())
394                  if (row==nrow_up-1) good[lab]|=0x4000;
395               if (sec>=digp->GetNInnerSector())
396                  if (row==nrow_up-1-gap) good[lab]|=0x1000;
397
398               if (sec>=digp->GetNInnerSector())
399                  if (row==nrow_up-1-shift) good[lab]|=0x2000;
400               if (sec>=digp->GetNInnerSector())
401                  if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
402
403               good[lab]++;
404           }
405           ca.ClearRow(sec,row);
406       }
407       delete pca;
408       tpcl->UnloadRecPoints();
409      }
410      break;
411    case 2:
412      {
413       tpcl->LoadDigits();
414       TTree *TD=tpcl->TreeD();
415       
416       AliSimDigits da, *digits=&da;
417       TD->GetBranch("Segment")->SetAddress(&digits);
418
419       Int_t *count = new Int_t[np];
420       Int_t i;
421       for (i=0; i<np; i++) count[i]=0;
422
423       Int_t sectors_by_rows=(Int_t)TD->GetEntries();
424       for (i=0; i<sectors_by_rows; i++) {
425           if (!TD->GetEvent(i)) continue;
426           Int_t sec,row;
427           digp->AdjustSectorRow(digits->GetID(),sec,row);
428           cerr<<sec<<' '<<row<<"                                     \r";
429           digits->First();
430           do { //Many thanks to J.Chudoba who noticed this
431               Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
432               Short_t dig = digits->GetDigit(it,ip);
433               Int_t idx0=digits->GetTrackID(it,ip,0); 
434               Int_t idx1=digits->GetTrackID(it,ip,1);
435               Int_t idx2=digits->GetTrackID(it,ip,2);
436               if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
437               if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
438               if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
439           } while (digits->Next());
440           for (Int_t j=0; j<np; j++) {
441               if (count[j]>1) {
442                  if (sec>=digp->GetNInnerSector())
443                    if (row==nrow_up-1    ) good[j]|=0x4000;
444                  if (sec>=digp->GetNInnerSector())
445                    if (row==nrow_up-1-gap) good[j]|=0x1000;
446
447                  if (sec>=digp->GetNInnerSector())
448                    if (row==nrow_up-1-shift) good[j]|=0x2000;
449                  if (sec>=digp->GetNInnerSector())
450                    if (row==nrow_up-1-gap-shift) good[j]|=0x800;
451                  good[j]++;
452               }
453               count[j]=0;
454           }
455       }
456       delete[] count;
457       tpcl->UnloadDigits();
458      }
459       break;
460    default:
461       cerr<<"Invalid TPC version !\n";
462       delete rl;
463       exit(7);
464    }
465
466    rl->LoadKinematics();
467    AliStack* stack = rl->Stack();
468    /** select tracks which are "good" enough **/
469    for (i=0; i<np; i++) {
470       if ((good[i]&0x5000) != 0x5000)
471       if ((good[i]&0x2800) != 0x2800) continue;
472       if ((good[i]&0x7FF ) < good_number) continue;
473
474       TParticle *p = (TParticle*)stack->Particle(i);
475       if (p == 0x0)
476        {
477          cerr<<"Can not get particle "<<i<<endl;
478          continue;
479        }
480       if (p->Pt()<0.100) continue;
481       if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
482
483       Int_t j=p->GetFirstMother();
484       if (j>=0) {
485         TParticle *pp = (TParticle*)stack->Particle(j);
486         if (pp == 0x0)
487          {
488            cerr<<"Can not get particle "<<j<<endl;
489            continue;
490          }
491         if (pp->GetFirstMother()>=0) continue;//only one decay is allowed
492         /*  for cascade hyperons only
493         Int_t jj=pp->GetFirstMother();
494         if (jj>=0) {
495           TParticle *ppp = (TParticle*)stack->Particle(jj);
496           if (ppp->GetFirstMother()>=0) continue;//two decays are allowed
497         }
498         */
499       }
500
501       gt[nt].lab=i;
502       gt[nt].code=p->GetPdgCode();
503       gt[nt].px=0.; gt[nt].py=0.; gt[nt].pz=0.;
504       gt[nt].x=0.; gt[nt].y=0.; gt[nt].z=0.;
505       nt++;
506       if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
507       cerr<<np-i<<"                \r";
508    }
509    rl->UnloadKinematics();
510
511    /** check if there is also information at the entrance of the TPC **/
512    
513    tpcl->LoadHits();
514    TTree *TH=tpcl->TreeH();
515    TPC->SetTreeAddress();
516    np=(Int_t)TH->GetEntries();
517    for (i=0; i<np; i++) {
518       TPC->ResetHits();
519       TH->GetEvent(i);
520       AliTPChit *phit = (AliTPChit*)TPC->FirstHit(-1);
521       for ( ; phit; phit=(AliTPChit*)TPC->NextHit() ) {
522          if (phit->fQ !=0. ) continue;
523
524          Double_t px=phit->X(), py=phit->Y(), pz=phit->Z();
525
526          if ((phit=(AliTPChit*)TPC->NextHit())==0) break;
527          if (phit->fQ != 0.) continue;
528
529          Double_t x=phit->X(), y=phit->Y(), z=phit->Z();
530          if (TMath::Sqrt(x*x+y*y)>90.) continue;
531
532          Int_t j, lab=phit->Track();
533          for (j=0; j<nt; j++) {if (gt[j].lab==lab) break;}
534          if (j==nt) continue;         
535
536          // (px,py,pz) - in global coordinate system, (x,y,z) - in local !
537          gt[j].px=px; gt[j].py=py; gt[j].pz=pz;
538          Float_t cs,sn; digp->AdjustCosSin(phit->fSector,cs,sn);
539          gt[j].x = x*cs + y*sn;
540          gt[j].y =-x*sn + y*cs;
541          gt[j].z = z;
542       }
543       cerr<<np-i<<"                \r";
544    }
545
546    delete[] good;
547    
548    tpcl->UnloadHits();
549    rl->UnloadgAlice();
550
551    gBenchmark->Stop("AliTPCComparison");
552    gBenchmark->Show("AliTPCComparison");   
553    return nt;
554 }