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