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