1 /****************************************************************************
2 * Very important, delicate and rather obscure macro. *
4 * Creates list of "trackable" tracks, *
5 * calculates efficiency, resolutions etc. *
6 * There is a possibility to run this macro over several events. *
8 * Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch *
9 * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de *
10 ****************************************************************************/
12 #if !defined(__CINT__) || defined(__MAKECINT__)
15 #include <Riostream.h>
19 #include <TParticle.h>
23 #include <TBenchmark.h>
29 #include "AliHeader.h"
30 #include "AliTrackReference.h"
31 #include "AliRunLoader.h"
35 #include "AliSimDigits.h"
37 #include "AliTPCParamSR.h"
38 #include "AliTPCClustersArray.h"
39 #include "AliTPCClustersRow.h"
40 #include "AliTPCcluster.h"
41 #include "AliTPCLoader.h"
44 Int_t GoodTracksTPC(const Char_t *dir=".");
46 extern AliRun *gAlice;
47 extern TBenchmark *gBenchmark;
50 static Int_t allgood=0;
51 static Int_t allselected=0;
52 static Int_t allfound=0;
54 Int_t AliTPCComparison
55 (Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
56 gBenchmark->Start("AliTPCComparison");
58 ::Info("AliTPCComparison.C","Doing comparison...");
62 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
63 if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.);
66 TH1F *hl=(TH1F*)gROOT->FindObject("hl");
67 if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
70 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
71 if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
74 TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
76 hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
77 hmpt->SetFillColor(6);
81 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
82 if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
84 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
85 if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);
87 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
88 if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);
90 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
91 if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
92 hg->SetLineColor(4); hg->SetLineWidth(2);
94 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
95 if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
96 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
98 TH1F *he=(TH1F*)gROOT->FindObject("he");
100 he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
102 TH2F *hep=(TH2F*)gROOT->FindObject("hep");
103 if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
104 hep->SetMarkerStyle(8);
105 hep->SetMarkerSize(0.4);
109 sprintf(fname,"%s/GoodTracksTPC.root",dir);
111 TFile *refFile=TFile::Open(fname,"old");
112 if (!refFile || !refFile->IsOpen()) {
113 ::Info("AliTPCComparison.C","Marking good tracks (will take a while)...");
114 if (GoodTracksTPC(dir)) {
115 ::Error("AliTPCComparison.C","Can't generate the reference file !");
119 refFile=TFile::Open(fname,"old");
120 if (!refFile || !refFile->IsOpen()) {
121 ::Error("AliTPCComparison.C","Can't open the reference file !");
125 TTree *tpcTree=(TTree*)refFile->Get("tpcTree");
127 ::Error("AliTPCComparison.C","Can't get the reference tree !");
130 TBranch *branch=tpcTree->GetBranch("TPC");
132 ::Error("AliTPCComparison.C","Can't get the TPC branch !");
135 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
136 branch->SetAddress(&refs);
139 sprintf(fname,"%s/AliESDs.root",dir);
140 TFile *ef=TFile::Open(fname);
141 if ((!ef)||(!ef->IsOpen())) {
142 sprintf(fname,"%s/AliESDtpc.root",dir);
143 ef=TFile::Open(fname);
144 if ((!ef)||(!ef->IsOpen())) {
145 ::Error("AliTPCComparison.C","Can't open AliESDtpc.root !");
150 TIter next(ef->GetListOfKeys());
153 //******* Loop over events *********
155 while ((key=(TKey*)next())!=0) {
156 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
158 AliESD *event=(AliESD*)key->ReadObj();
160 Int_t nentr=event->GetNumberOfTracks();
163 if (tpcTree->GetEvent(e++)==0) {
164 cerr<<"No reconstructable tracks !\n";
168 Int_t ngood=refs->GetEntriesFast();
171 const Int_t MAX=15000;
173 Int_t track_notfound[MAX], itrack_notfound=0;
174 Int_t track_fake[MAX], itrack_fake=0;
175 Int_t track_multifound[MAX],track_multifound_n[MAX],itrack_multifound=0;
178 for (Int_t k=0; k<ngood; k++) {
179 AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k);
180 Int_t lab=ref->Label(), tlab=-1;
181 Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());
183 if (ptg<1e-33) continue; // for those not crossing 0 pad row
185 if (ptg<ptcutl) continue;
186 if (ptg>ptcuth) continue;
192 AliESDtrack *track=0;
193 for (i=0; i<nentr; i++) {
194 track=event->GetTrack(i);
195 tlab=track->GetTPCLabel();
196 if (lab==TMath::Abs(tlab)) break;
199 track_notfound[itrack_notfound++]=lab;
203 //MI change - addition
206 AliESDtrack * mitrack;
207 for (mi=0; mi<nentr; mi++) {
208 mitrack=event->GetTrack(mi);
209 if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++;
212 track_multifound[itrack_multifound]=lab;
213 track_multifound_n[itrack_multifound]=micount;
216 if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
217 if (lab==tlab) hfound->Fill(ptg);
219 track_fake[itrack_fake++]=lab;
223 Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz);
224 Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]);
225 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
226 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
227 Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]);
228 Float_t lam=TMath::ATan2(pxpypz[2],pt);
231 Int_t pdg=(Int_t)ref->GetLength(); //this is particle's PDG !
233 if (TMath::Abs(pdg)==11 && ptg>4.) {//high momentum electrons
234 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
236 Float_t phig=TMath::ATan2(ref->Py(),ref->Px());
237 hp->Fill((phi - phig)*1000.);
239 Float_t lamg=TMath::ATan2(ref->Pz(),ptg);
240 hl->Fill((lam - lamg)*1000.);
242 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
245 Float_t mom=pt/TMath::Cos(lam);
246 Float_t dedx=track->GetTPCsignal();
247 hep->Fill(mom,dedx,1.);
248 if (TMath::Abs(pdg)==211) //pions
249 if (mom>0.4 && mom<0.5) {
254 cout<<"\nList of Not found tracks :\n";
255 for ( i = 0; i< itrack_notfound; i++){
256 cout<<track_notfound[i]<<"\t";
257 if ((i%5)==4) cout<<"\n";
259 cout<<"\nList of fake tracks :\n";
260 for ( i = 0; i< itrack_fake; i++){
261 cout<<track_fake[i]<<"\t";
262 if ((i%5)==4) cout<<"\n";
264 cout<<"\nList of multiple found tracks :\n";
265 for ( i=0; i<itrack_multifound; i++) {
266 cout<<"id. "<<track_multifound[i]
267 <<" found - "<<track_multifound_n[i]<<"times\n";
270 cout<<"Number of found tracks ="<<nentr<<endl;
271 cout<<"Number of \"good\" tracks ="<<ngood<<endl;
275 }// ***** End of the loop over events
282 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
283 if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
284 cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl;
285 cout<<"Total number of found tracks ="<<allfound<<endl;
286 cout<<"Total number of \"good\" tracks ="<<allgood<<endl;
289 gStyle->SetOptStat(111110);
290 gStyle->SetOptFit(1);
292 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
296 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
297 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
298 hp->SetFillColor(4); hp->SetXTitle("(mrad)");
299 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();
301 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
302 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
303 hl->SetXTitle("(mrad)");
304 if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();
306 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
307 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
308 hpt->SetXTitle("(%)");
309 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();
311 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
312 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
313 hmpt->SetXTitle("(%)");
314 if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd();
316 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
317 p5->SetFillColor(41); p5->SetFrameFillColor(10);
318 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
319 hg->Divide(hfound,hgood,1,1.,"b");
320 hf->Divide(hfake,hgood,1,1.,"b");
322 hg->SetYTitle("Tracking efficiency");
323 hg->SetXTitle("Pt (GeV/c)");
326 TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
328 TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
332 hf->SetFillStyle(3013);
335 hf->Draw("histsame");
336 TText *text = new TText(0.461176,0.248448,"Fake tracks");
337 text->SetTextSize(0.05);
339 text = new TText(0.453919,1.11408,"Good tracks");
340 text->SetTextSize(0.05);
345 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
347 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
348 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
349 he->SetFillColor(2); he->SetFillStyle(3005);
350 he->SetXTitle("Arbitrary Units");
351 if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();
353 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
354 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
355 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
356 hep->Draw(); c1->cd();
358 TFile fc("AliTPCComparison.root","RECREATE");
363 gBenchmark->Stop("AliTPCComparison");
364 gBenchmark->Show("AliTPCComparison");
370 Int_t GoodTracksTPC(const Char_t *dir) {
372 delete gAlice->GetRunLoader();
373 delete gAlice;//if everything was OK here it is already NULL
378 sprintf(fname,"%s/galice.root",dir);
380 AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
382 ::Error("GoodTracksTPC","Can't start session !");
388 rl->LoadKinematics();
391 AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
393 ::Error("GoodTracksTPC","Can not find TPCLoader !");
398 AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
399 Int_t ver = TPC->IsVersion();
400 cout<<"TPC version "<<ver<<" has been found !\n";
401 if (ver==1) tpcl->LoadRecPoints();
402 else if (ver==2) tpcl->LoadDigits();
404 ::Error("GoodTracksTPC","Wrong TPC version !");
410 AliTPCParamSR *digp=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
412 ::Error("AliTPCHits2Digits.C","TPC parameters have not been found !");
418 Int_t nrow_up=digp->GetNRowUp();
419 Int_t nrows=digp->GetNRowLow()+nrow_up;
420 Int_t zero=digp->GetZeroSup();
421 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
422 Int_t good_number=Int_t(0.4*nrows);
424 Int_t nev=rl->GetNumberOfEvents();
425 ::Info("GoodTracksTPC","Number of events : %d\n",nev);
428 sprintf(fname,"%s/GoodTracksTPC.root",dir);
429 TFile *file=TFile::Open(fname,"recreate");
431 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
432 TTree tpcTree("tpcTree","Tree with info about the reconstructable TPC tracks");
433 tpcTree.Branch("TPC",&refs);
435 //******** Loop over generated events
436 for (Int_t e=0; e<nev; e++) {
442 rl->GetEvent(e); file->cd();
444 Int_t np = rl->GetHeader()->GetNtrack();
445 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
447 //******** Fill the "good" masks
448 Int_t *good=new Int_t[np]; for (i=0; i<np; i++) good[i]=0;
452 AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
454 ca.SetClusterType("AliTPCcluster");
455 ca.ConnectTree(tpcl->TreeR());
456 Int_t nrows=Int_t(ca.GetTree()->GetEntries());
457 for (Int_t n=0; n<nrows; n++) {
458 AliSegmentID *s=ca.LoadEntry(n);
460 digp->AdjustSectorRow(s->GetID(),sec,row);
461 AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
462 Int_t ncl=clrow.GetArray()->GetEntriesFast();
464 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
465 Int_t lab=c->GetLabel(0);
466 if (lab<0) continue; //noise cluster
469 if (sec>=digp->GetNInnerSector())
470 if (row==nrow_up-1) good[lab]|=0x4000;
471 if (sec>=digp->GetNInnerSector())
472 if (row==nrow_up-1-gap) good[lab]|=0x1000;
474 if (sec>=digp->GetNInnerSector())
475 if (row==nrow_up-1-shift) good[lab]|=0x2000;
476 if (sec>=digp->GetNInnerSector())
477 if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
481 ca.ClearRow(sec,row);
488 TTree *TD=tpcl->TreeD();
490 AliSimDigits da, *digits=&da;
491 TD->GetBranch("Segment")->SetAddress(&digits);
493 Int_t *count = new Int_t[np];
495 for (i=0; i<np; i++) count[i]=0;
497 Int_t sectors_by_rows=(Int_t)TD->GetEntries();
498 for (i=0; i<sectors_by_rows; i++) {
499 if (!TD->GetEvent(i)) continue;
501 digp->AdjustSectorRow(digits->GetID(),sec,row);
502 //cerr<<sec<<' '<<row<<'\r';
504 do { //Many thanks to J.Chudoba who noticed this
505 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
506 Short_t dig = digits->GetDigit(it,ip);
507 Int_t idx0=digits->GetTrackID(it,ip,0);
508 Int_t idx1=digits->GetTrackID(it,ip,1);
509 Int_t idx2=digits->GetTrackID(it,ip,2);
510 if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
511 if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
512 if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
513 } while (digits->Next());
514 for (Int_t j=0; j<np; j++) {
516 if (sec>=digp->GetNInnerSector())
517 if (row==nrow_up-1 ) good[j]|=0x4000;
518 if (sec>=digp->GetNInnerSector())
519 if (row==nrow_up-1-gap) good[j]|=0x1000;
521 if (sec>=digp->GetNInnerSector())
522 if (row==nrow_up-1-shift) good[j]|=0x2000;
523 if (sec>=digp->GetNInnerSector())
524 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
536 //****** select tracks which are "good" enough
537 AliStack* stack = rl->Stack();
538 for (i=0; i<np; i++) {
539 if ((good[i]&0x5000) != 0x5000)
540 if ((good[i]&0x2800) != 0x2800) continue;
541 if ((good[i]&0x7FF ) < good_number) continue;
543 TParticle *p = (TParticle*)stack->Particle(i);
545 cerr<<"Can not get particle "<<i<<endl;
548 if (p->Pt()<0.100) continue;
549 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
551 Double_t vx=p->Vx(),vy=p->Vy();
552 if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue;
554 AliTrackReference *ref=new((*refs)[nt]) AliTrackReference();
557 Int_t pdg=p->GetPdgCode();
558 ref->SetLength(pdg); //This will the particle's PDG !
559 ref->SetMomentum(0.,0.,0.);
560 ref->SetPosition(0.,0.,0.);
565 //**** check if there is also information at the entrance of the TPC
566 TTree *TR=rl->TreeTR();
567 TBranch *branch=TR->GetBranch("TPC");
569 ::Error("GoodTracksTPC","No track references !");
573 TClonesArray tpcdummy("AliTrackReference",1000), *tpcRefs=&tpcdummy;
574 branch->SetAddress(&tpcRefs);
575 Int_t nr=(Int_t)TR->GetEntries();
576 for (Int_t r=0; r<nr; r++) {
577 //cerr<<r<<' '<<nr<<'\r';
579 if (tpcRefs->GetEntriesFast()==0) continue;
580 AliTrackReference *tpcRef=(AliTrackReference*)tpcRefs->UncheckedAt(0);
582 AliTrackReference *ref=0;
583 for (j=0; j<nt; j++) {
584 ref=(AliTrackReference *)refs->UncheckedAt(j);
585 if (ref->Label()==tpcRef->Label()) break;
588 if (ref==0) continue;
590 ref->SetMomentum(tpcRef->Px(),tpcRef->Py(),tpcRef->Pz());
591 ref->SetPosition(tpcRef->LocalX(),tpcRef->LocalY(),tpcRef->Z());
599 } //****** end of the loop over generated events