1 /// \file AliTPCComparison.C
2 /// \brief Very important, delicate and rather obscure macro
4 /// Creates list of "trackable" tracks, calculates efficiency, resolutions etc.
6 /// There is a possibility to run this macro over several events.
8 /// \author I.Belikov, CERN, Jouri.Belikov@cern.ch, M.Ivanov, GSI, m.ivanov@gsi.de
10 #if !defined(__CINT__) || defined(__MAKECINT__)
13 #include <Riostream.h>
17 #include <TParticle.h>
21 #include <TBenchmark.h>
27 #include "AliHeader.h"
28 #include "AliTrackReference.h"
29 #include "AliRunLoader.h"
31 #include "AliESDEvent.h"
32 #include "AliESDtrack.h"
34 #include "AliSimDigits.h"
36 #include "AliTPCParamSR.h"
37 #include "AliTPCClustersArray.h"
38 #include "AliTPCClustersRow.h"
39 #include "AliTPCcluster.h"
40 #include "AliTPCLoader.h"
42 #include "AliCDBManager.h"
43 #include "AliTPCcalibDB.h"
46 Int_t GoodTracksTPC(const Char_t *dir=".");
48 extern TBenchmark *gBenchmark;
51 static Int_t allgood=0;
52 static Int_t allselected=0;
53 static Int_t allfound=0;
55 Int_t AliTPCComparison
56 (Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
57 gBenchmark->Start("AliTPCComparison");
59 ::Info("AliTPCComparison.C","Doing comparison...");
63 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
64 if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.);
67 TH1F *hl=(TH1F*)gROOT->FindObject("hl");
68 if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
71 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
72 if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
75 TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
77 hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
78 hmpt->SetFillColor(6);
82 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
83 if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
85 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
86 if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);
88 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
89 if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);
91 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
92 if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
93 hg->SetLineColor(4); hg->SetLineWidth(2);
95 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
96 if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
97 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
99 TH1F *he=(TH1F*)gROOT->FindObject("he");
101 he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
103 TH2F *hep=(TH2F*)gROOT->FindObject("hep");
104 if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
105 hep->SetMarkerStyle(8);
106 hep->SetMarkerSize(0.4);
110 sprintf(fname,"%s/GoodTracksTPC.root",dir);
112 TFile *refFile=TFile::Open(fname,"old");
113 if (!refFile || !refFile->IsOpen()) {
114 ::Info("AliTPCComparison.C","Marking good tracks (will take a while)...");
115 if (GoodTracksTPC(dir)) {
116 ::Error("AliTPCComparison.C","Can't generate the reference file !");
120 refFile=TFile::Open(fname,"old");
121 if (!refFile || !refFile->IsOpen()) {
122 ::Error("AliTPCComparison.C","Can't open the reference file !");
126 TTree *tpcTree=(TTree*)refFile->Get("tpcTree");
128 ::Error("AliTPCComparison.C","Can't get the reference tree !");
131 TBranch *branch=tpcTree->GetBranch("TPC");
133 ::Error("AliTPCComparison.C","Can't get the TPC branch !");
136 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
137 branch->SetAddress(&refs);
140 sprintf(fname,"%s/AliESDs.root",dir);
141 TFile *ef=TFile::Open(fname);
142 if ((!ef)||(!ef->IsOpen())) {
143 sprintf(fname,"%s/AliESDtpc.root",dir);
144 ef=TFile::Open(fname);
145 if ((!ef)||(!ef->IsOpen())) {
146 ::Error("AliTPCComparison.C","Can't open AliESDtpc.root !");
150 AliESDEvent* event = new AliESDEvent();
151 TTree* esdTree = (TTree*) ef->Get("esdTree");
153 ::Error("AliTPCComparison.C", "no ESD tree found");
156 event->ReadFromTree(esdTree);
159 // ******* Loop over events *********
162 while (esdTree->GetEvent(e)) {
163 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
165 Int_t nentr=event->GetNumberOfTracks();
168 if (tpcTree->GetEvent(e++)==0) {
169 cerr<<"No reconstructable tracks !\n";
173 Int_t ngood=refs->GetEntriesFast();
176 const Int_t MAX=15000;
178 Int_t track_notfound[MAX], itrack_notfound=0;
179 Int_t track_fake[MAX], itrack_fake=0;
180 Int_t track_multifound[MAX],track_multifound_n[MAX],itrack_multifound=0;
183 for (Int_t k=0; k<ngood; k++) {
184 AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k);
185 Int_t lab=ref->Label(), tlab=-1;
186 Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());
188 if (ptg<1e-33) continue; // for those not crossing 0 pad row
190 if (ptg<ptcutl) continue;
191 if (ptg>ptcuth) continue;
197 AliESDtrack *track=0;
198 for (i=0; i<nentr; i++) {
199 track=event->GetTrack(i);
200 tlab=track->GetTPCLabel();
201 if (lab==TMath::Abs(tlab)) break;
204 track_notfound[itrack_notfound++]=lab;
208 //MI change - addition
211 AliESDtrack * mitrack;
212 for (mi=0; mi<nentr; mi++) {
213 mitrack=event->GetTrack(mi);
214 if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++;
217 track_multifound[itrack_multifound]=lab;
218 track_multifound_n[itrack_multifound]=micount;
221 if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
222 if (lab==tlab) hfound->Fill(ptg);
224 track_fake[itrack_fake++]=lab;
228 Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz);
229 Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]);
230 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
231 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
232 Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]);
233 Float_t lam=TMath::ATan2(pxpypz[2],pt);
236 Int_t pdg=(Int_t)ref->GetLength(); //this is particle's PDG !
238 if (TMath::Abs(pdg)==11 && ptg>4.) {//high momentum electrons
239 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
241 Float_t phig=TMath::ATan2(ref->Py(),ref->Px());
242 hp->Fill((phi - phig)*1000.);
244 Float_t lamg=TMath::ATan2(ref->Pz(),ptg);
245 hl->Fill((lam - lamg)*1000.);
247 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
250 Float_t mom=pt/TMath::Cos(lam);
251 Float_t dedx=track->GetTPCsignal();
252 hep->Fill(mom,dedx,1.);
253 if (TMath::Abs(pdg)==211) //pions
254 if (mom>0.4 && mom<0.5) {
259 cout<<"\nList of Not found tracks :\n";
260 for ( i = 0; i< itrack_notfound; i++){
261 cout<<track_notfound[i]<<"\t";
262 if ((i%5)==4) cout<<"\n";
264 cout<<"\nList of fake tracks :\n";
265 for ( i = 0; i< itrack_fake; i++){
266 cout<<track_fake[i]<<"\t";
267 if ((i%5)==4) cout<<"\n";
269 cout<<"\nList of multiple found tracks :\n";
270 for ( i=0; i<itrack_multifound; i++) {
271 cout<<"id. "<<track_multifound[i]
272 <<" found - "<<track_multifound_n[i]<<"times\n";
275 cout<<"Number of found tracks ="<<nentr<<endl;
276 cout<<"Number of \"good\" tracks ="<<ngood<<endl;
279 }// ***** End of the loop over events
288 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
289 if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
290 cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl;
291 cout<<"Total number of found tracks ="<<allfound<<endl;
292 cout<<"Total number of \"good\" tracks ="<<allgood<<endl;
295 gStyle->SetOptStat(111110);
296 gStyle->SetOptFit(1);
298 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
302 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
303 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
304 hp->SetFillColor(4); hp->SetXTitle("(mrad)");
305 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();
307 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
308 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
309 hl->SetXTitle("(mrad)");
310 if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();
312 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
313 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
314 hpt->SetXTitle("(%)");
315 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();
317 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
318 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
319 hmpt->SetXTitle("(%)");
320 if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd();
322 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
323 p5->SetFillColor(41); p5->SetFrameFillColor(10);
324 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
325 hg->Divide(hfound,hgood,1,1.,"b");
326 hf->Divide(hfake,hgood,1,1.,"b");
328 hg->SetYTitle("Tracking efficiency");
329 hg->SetXTitle("Pt (GeV/c)");
332 TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
334 TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
338 hf->SetFillStyle(3013);
341 hf->Draw("histsame");
342 TText *text = new TText(0.461176,0.248448,"Fake tracks");
343 text->SetTextSize(0.05);
345 text = new TText(0.453919,1.11408,"Good tracks");
346 text->SetTextSize(0.05);
351 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
353 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
354 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
355 he->SetFillColor(2); he->SetFillStyle(3005);
356 he->SetXTitle("Arbitrary Units");
357 if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();
359 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
360 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
361 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
362 hep->Draw(); c1->cd();
364 TFile fc("AliTPCComparison.root","RECREATE");
369 gBenchmark->Stop("AliTPCComparison");
370 gBenchmark->Show("AliTPCComparison");
376 Int_t GoodTracksTPC(const Char_t *dir) {
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 AliCDBManager *man=AliCDBManager::Instance();
411 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
414 (AliTPCParamSR*)(AliTPCcalibDB::Instance()->GetParameters());
416 ::Error("AliTPCComparison.C","TPC parameters have not been found !");
421 Int_t nrow_up=digp->GetNRowUp();
422 Int_t nrows=digp->GetNRowLow()+nrow_up;
423 Int_t zero=digp->GetZeroSup();
424 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
425 Int_t good_number=Int_t(0.4*nrows);
427 Int_t nev=rl->GetNumberOfEvents();
428 ::Info("GoodTracksTPC","Number of events : %d\n",nev);
431 sprintf(fname,"%s/GoodTracksTPC.root",dir);
432 TFile *file=TFile::Open(fname,"recreate");
434 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
435 TTree tpcTree("tpcTree","Tree with info about the reconstructable TPC tracks");
436 tpcTree.Branch("TPC",&refs);
438 //******** Loop over generated events
439 for (Int_t e=0; e<nev; e++) {
445 rl->GetEvent(e); file->cd();
447 Int_t np = rl->GetHeader()->GetNtrack();
448 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
450 //******** Fill the "good" masks
451 Int_t *good=new Int_t[np]; for (i=0; i<np; i++) good[i]=0;
455 AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
457 ca.SetClusterType("AliTPCcluster");
458 ca.ConnectTree(tpcl->TreeR());
459 Int_t nrows=Int_t(ca.GetTree()->GetEntries());
460 for (Int_t n=0; n<nrows; n++) {
461 AliSegmentID *s=ca.LoadEntry(n);
463 digp->AdjustSectorRow(s->GetID(),sec,row);
464 AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
465 Int_t ncl=clrow.GetArray()->GetEntriesFast();
467 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
468 Int_t lab=c->GetLabel(0);
469 if (lab<0) continue; //noise cluster
472 if (sec>=digp->GetNInnerSector())
473 if (row==nrow_up-1) good[lab]|=0x4000;
474 if (sec>=digp->GetNInnerSector())
475 if (row==nrow_up-1-gap) good[lab]|=0x1000;
477 if (sec>=digp->GetNInnerSector())
478 if (row==nrow_up-1-shift) good[lab]|=0x2000;
479 if (sec>=digp->GetNInnerSector())
480 if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
484 ca.ClearRow(sec,row);
491 TTree *TD=tpcl->TreeD();
493 AliSimDigits da, *digits=&da;
494 TD->GetBranch("Segment")->SetAddress(&digits);
496 Int_t *count = new Int_t[np];
498 for (i=0; i<np; i++) count[i]=0;
500 Int_t sectors_by_rows=(Int_t)TD->GetEntries();
501 for (i=0; i<sectors_by_rows; i++) {
502 if (!TD->GetEvent(i)) continue;
504 digp->AdjustSectorRow(digits->GetID(),sec,row);
505 //cerr<<sec<<' '<<row<<'\r';
507 do { //Many thanks to J.Chudoba who noticed this
508 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
509 Short_t dig = digits->GetDigit(it,ip);
510 Int_t idx0=digits->GetTrackID(it,ip,0);
511 Int_t idx1=digits->GetTrackID(it,ip,1);
512 Int_t idx2=digits->GetTrackID(it,ip,2);
513 if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
514 if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
515 if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
516 } while (digits->Next());
517 for (Int_t j=0; j<np; j++) {
519 if (sec>=digp->GetNInnerSector())
520 if (row==nrow_up-1 ) good[j]|=0x4000;
521 if (sec>=digp->GetNInnerSector())
522 if (row==nrow_up-1-gap) good[j]|=0x1000;
524 if (sec>=digp->GetNInnerSector())
525 if (row==nrow_up-1-shift) good[j]|=0x2000;
526 if (sec>=digp->GetNInnerSector())
527 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
539 //****** select tracks which are "good" enough
540 AliStack* stack = rl->Stack();
541 for (i=0; i<np; i++) {
542 if ((good[i]&0x5000) != 0x5000)
543 if ((good[i]&0x2800) != 0x2800) continue;
544 if ((good[i]&0x7FF ) < good_number) continue;
546 TParticle *p = (TParticle*)stack->Particle(i);
548 cerr<<"Can not get particle "<<i<<endl;
551 if (p->Pt()<0.100) continue;
552 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
554 Double_t vx=p->Vx(),vy=p->Vy(),vz=p->Vz();
555 if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue;
556 if (TMath::Abs(vz) > 50.) continue;
558 AliTrackReference *ref=new((*refs)[nt]) AliTrackReference();
561 Int_t pdg=p->GetPdgCode();
562 ref->SetLength(pdg); //This will the particle's PDG !
563 ref->SetMomentum(0.,0.,0.);
564 ref->SetPosition(0.,0.,0.);
569 //**** check if there is also information at the entrance of the TPC
570 TTree *TR=rl->TreeTR();
571 TBranch *branch=TR->GetBranch("TrackReferences");
573 ::Error("GoodTracksTPC","No track references !");
577 TClonesArray tpcdummy("AliTrackReference",1000), *tpcRefs=&tpcdummy;
578 branch->SetAddress(&tpcRefs);
579 Int_t nr=(Int_t)TR->GetEntries();
580 for (Int_t r=0; r<nr; r++) {
581 //cerr<<r<<' '<<nr<<'\r';
584 Int_t nref = tpcRefs->GetEntriesFast();
586 AliTrackReference *tpcRef= 0x0;
587 for (Int_t iref=0; iref<nref; ++iref) {
588 tpcRef = (AliTrackReference*)tpcRefs->UncheckedAt(iref);
589 if (tpcRef->DetectorId() == AliTrackReference::kTPC) break;
593 if (!tpcRef) continue;
596 AliTrackReference *ref=0;
597 for (j=0; j<nt; j++) {
598 ref=(AliTrackReference *)refs->UncheckedAt(j);
599 if (ref->Label()==tpcRef->Label()) break;
602 if (ref==0) continue;
604 ref->SetMomentum(tpcRef->Px(),tpcRef->Py(),tpcRef->Pz());
605 ref->SetPosition(tpcRef->LocalX(),tpcRef->LocalY(),tpcRef->Z());
613 } //****** end of the loop over generated events