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