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