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