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"
33 #include "AliESDEvent.h"
34 #include "AliESDtrack.h"
36 #include "AliSimDigits.h"
38 #include "AliTPCParamSR.h"
39 #include "AliTPCClustersArray.h"
40 #include "AliTPCClustersRow.h"
41 #include "AliTPCcluster.h"
42 #include "AliTPCLoader.h"
44 #include "AliCDBManager.h"
45 #include "AliTPCcalibDB.h"
48 Int_t GoodTracksTPC(const Char_t *dir=".");
50 extern TBenchmark *gBenchmark;
53 static Int_t allgood=0;
54 static Int_t allselected=0;
55 static Int_t allfound=0;
57 Int_t AliTPCComparison
58 (Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
59 gBenchmark->Start("AliTPCComparison");
61 ::Info("AliTPCComparison.C","Doing comparison...");
65 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
66 if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.);
69 TH1F *hl=(TH1F*)gROOT->FindObject("hl");
70 if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
73 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
74 if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
77 TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
79 hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
80 hmpt->SetFillColor(6);
84 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
85 if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
87 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
88 if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);
90 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
91 if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);
93 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
94 if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
95 hg->SetLineColor(4); hg->SetLineWidth(2);
97 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
98 if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
99 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
101 TH1F *he=(TH1F*)gROOT->FindObject("he");
103 he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
105 TH2F *hep=(TH2F*)gROOT->FindObject("hep");
106 if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
107 hep->SetMarkerStyle(8);
108 hep->SetMarkerSize(0.4);
112 sprintf(fname,"%s/GoodTracksTPC.root",dir);
114 TFile *refFile=TFile::Open(fname,"old");
115 if (!refFile || !refFile->IsOpen()) {
116 ::Info("AliTPCComparison.C","Marking good tracks (will take a while)...");
117 if (GoodTracksTPC(dir)) {
118 ::Error("AliTPCComparison.C","Can't generate the reference file !");
122 refFile=TFile::Open(fname,"old");
123 if (!refFile || !refFile->IsOpen()) {
124 ::Error("AliTPCComparison.C","Can't open the reference file !");
128 TTree *tpcTree=(TTree*)refFile->Get("tpcTree");
130 ::Error("AliTPCComparison.C","Can't get the reference tree !");
133 TBranch *branch=tpcTree->GetBranch("TPC");
135 ::Error("AliTPCComparison.C","Can't get the TPC branch !");
138 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
139 branch->SetAddress(&refs);
142 sprintf(fname,"%s/AliESDs.root",dir);
143 TFile *ef=TFile::Open(fname);
144 if ((!ef)||(!ef->IsOpen())) {
145 sprintf(fname,"%s/AliESDtpc.root",dir);
146 ef=TFile::Open(fname);
147 if ((!ef)||(!ef->IsOpen())) {
148 ::Error("AliTPCComparison.C","Can't open AliESDtpc.root !");
152 AliESDEvent* event = new AliESDEvent();
153 TTree* esdTree = (TTree*) ef->Get("esdTree");
155 ::Error("AliTPCComparison.C", "no ESD tree found");
158 event->ReadFromTree(esdTree);
161 //******* Loop over events *********
163 while (esdTree->GetEvent(e)) {
164 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
166 Int_t nentr=event->GetNumberOfTracks();
169 if (tpcTree->GetEvent(e++)==0) {
170 cerr<<"No reconstructable tracks !\n";
174 Int_t ngood=refs->GetEntriesFast();
177 const Int_t MAX=15000;
179 Int_t track_notfound[MAX], itrack_notfound=0;
180 Int_t track_fake[MAX], itrack_fake=0;
181 Int_t track_multifound[MAX],track_multifound_n[MAX],itrack_multifound=0;
184 for (Int_t k=0; k<ngood; k++) {
185 AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k);
186 Int_t lab=ref->Label(), tlab=-1;
187 Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());
189 if (ptg<1e-33) continue; // for those not crossing 0 pad row
191 if (ptg<ptcutl) continue;
192 if (ptg>ptcuth) continue;
198 AliESDtrack *track=0;
199 for (i=0; i<nentr; i++) {
200 track=event->GetTrack(i);
201 tlab=track->GetTPCLabel();
202 if (lab==TMath::Abs(tlab)) break;
205 track_notfound[itrack_notfound++]=lab;
209 //MI change - addition
212 AliESDtrack * mitrack;
213 for (mi=0; mi<nentr; mi++) {
214 mitrack=event->GetTrack(mi);
215 if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++;
218 track_multifound[itrack_multifound]=lab;
219 track_multifound_n[itrack_multifound]=micount;
222 if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
223 if (lab==tlab) hfound->Fill(ptg);
225 track_fake[itrack_fake++]=lab;
229 Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz);
230 Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]);
231 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
232 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
233 Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]);
234 Float_t lam=TMath::ATan2(pxpypz[2],pt);
237 Int_t pdg=(Int_t)ref->GetLength(); //this is particle's PDG !
239 if (TMath::Abs(pdg)==11 && ptg>4.) {//high momentum electrons
240 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
242 Float_t phig=TMath::ATan2(ref->Py(),ref->Px());
243 hp->Fill((phi - phig)*1000.);
245 Float_t lamg=TMath::ATan2(ref->Pz(),ptg);
246 hl->Fill((lam - lamg)*1000.);
248 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
251 Float_t mom=pt/TMath::Cos(lam);
252 Float_t dedx=track->GetTPCsignal();
253 hep->Fill(mom,dedx,1.);
254 if (TMath::Abs(pdg)==211) //pions
255 if (mom>0.4 && mom<0.5) {
260 cout<<"\nList of Not found tracks :\n";
261 for ( i = 0; i< itrack_notfound; i++){
262 cout<<track_notfound[i]<<"\t";
263 if ((i%5)==4) cout<<"\n";
265 cout<<"\nList of fake tracks :\n";
266 for ( i = 0; i< itrack_fake; i++){
267 cout<<track_fake[i]<<"\t";
268 if ((i%5)==4) cout<<"\n";
270 cout<<"\nList of multiple found tracks :\n";
271 for ( i=0; i<itrack_multifound; i++) {
272 cout<<"id. "<<track_multifound[i]
273 <<" found - "<<track_multifound_n[i]<<"times\n";
276 cout<<"Number of found tracks ="<<nentr<<endl;
277 cout<<"Number of \"good\" tracks ="<<ngood<<endl;
280 }// ***** End of the loop over events
289 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
290 if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
291 cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl;
292 cout<<"Total number of found tracks ="<<allfound<<endl;
293 cout<<"Total number of \"good\" tracks ="<<allgood<<endl;
296 gStyle->SetOptStat(111110);
297 gStyle->SetOptFit(1);
299 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
303 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
304 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
305 hp->SetFillColor(4); hp->SetXTitle("(mrad)");
306 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();
308 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
309 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
310 hl->SetXTitle("(mrad)");
311 if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();
313 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
314 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
315 hpt->SetXTitle("(%)");
316 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();
318 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
319 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
320 hmpt->SetXTitle("(%)");
321 if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd();
323 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
324 p5->SetFillColor(41); p5->SetFrameFillColor(10);
325 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
326 hg->Divide(hfound,hgood,1,1.,"b");
327 hf->Divide(hfake,hgood,1,1.,"b");
329 hg->SetYTitle("Tracking efficiency");
330 hg->SetXTitle("Pt (GeV/c)");
333 TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
335 TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
339 hf->SetFillStyle(3013);
342 hf->Draw("histsame");
343 TText *text = new TText(0.461176,0.248448,"Fake tracks");
344 text->SetTextSize(0.05);
346 text = new TText(0.453919,1.11408,"Good tracks");
347 text->SetTextSize(0.05);
352 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
354 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
355 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
356 he->SetFillColor(2); he->SetFillStyle(3005);
357 he->SetXTitle("Arbitrary Units");
358 if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();
360 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
361 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
362 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
363 hep->Draw(); c1->cd();
365 TFile fc("AliTPCComparison.root","RECREATE");
370 gBenchmark->Stop("AliTPCComparison");
371 gBenchmark->Show("AliTPCComparison");
377 Int_t GoodTracksTPC(const Char_t *dir) {
379 sprintf(fname,"%s/galice.root",dir);
381 AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
383 ::Error("GoodTracksTPC","Can't start session !");
389 rl->LoadKinematics();
392 AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
394 ::Error("GoodTracksTPC","Can not find TPCLoader !");
399 AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
400 Int_t ver = TPC->IsVersion();
401 cout<<"TPC version "<<ver<<" has been found !\n";
402 if (ver==1) tpcl->LoadRecPoints();
403 else if (ver==2) tpcl->LoadDigits();
405 ::Error("GoodTracksTPC","Wrong TPC version !");
411 AliCDBManager *man=AliCDBManager::Instance();
412 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
415 (AliTPCParamSR*)(AliTPCcalibDB::Instance()->GetParameters());
417 ::Error("AliTPCComparison.C","TPC parameters have not been found !");
422 Int_t nrow_up=digp->GetNRowUp();
423 Int_t nrows=digp->GetNRowLow()+nrow_up;
424 Int_t zero=digp->GetZeroSup();
425 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
426 Int_t good_number=Int_t(0.4*nrows);
428 Int_t nev=rl->GetNumberOfEvents();
429 ::Info("GoodTracksTPC","Number of events : %d\n",nev);
432 sprintf(fname,"%s/GoodTracksTPC.root",dir);
433 TFile *file=TFile::Open(fname,"recreate");
435 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
436 TTree tpcTree("tpcTree","Tree with info about the reconstructable TPC tracks");
437 tpcTree.Branch("TPC",&refs);
439 //******** Loop over generated events
440 for (Int_t e=0; e<nev; e++) {
446 rl->GetEvent(e); file->cd();
448 Int_t np = rl->GetHeader()->GetNtrack();
449 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
451 //******** Fill the "good" masks
452 Int_t *good=new Int_t[np]; for (i=0; i<np; i++) good[i]=0;
456 AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
458 ca.SetClusterType("AliTPCcluster");
459 ca.ConnectTree(tpcl->TreeR());
460 Int_t nrows=Int_t(ca.GetTree()->GetEntries());
461 for (Int_t n=0; n<nrows; n++) {
462 AliSegmentID *s=ca.LoadEntry(n);
464 digp->AdjustSectorRow(s->GetID(),sec,row);
465 AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
466 Int_t ncl=clrow.GetArray()->GetEntriesFast();
468 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
469 Int_t lab=c->GetLabel(0);
470 if (lab<0) continue; //noise cluster
473 if (sec>=digp->GetNInnerSector())
474 if (row==nrow_up-1) good[lab]|=0x4000;
475 if (sec>=digp->GetNInnerSector())
476 if (row==nrow_up-1-gap) good[lab]|=0x1000;
478 if (sec>=digp->GetNInnerSector())
479 if (row==nrow_up-1-shift) good[lab]|=0x2000;
480 if (sec>=digp->GetNInnerSector())
481 if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
485 ca.ClearRow(sec,row);
492 TTree *TD=tpcl->TreeD();
494 AliSimDigits da, *digits=&da;
495 TD->GetBranch("Segment")->SetAddress(&digits);
497 Int_t *count = new Int_t[np];
499 for (i=0; i<np; i++) count[i]=0;
501 Int_t sectors_by_rows=(Int_t)TD->GetEntries();
502 for (i=0; i<sectors_by_rows; i++) {
503 if (!TD->GetEvent(i)) continue;
505 digp->AdjustSectorRow(digits->GetID(),sec,row);
506 //cerr<<sec<<' '<<row<<'\r';
508 do { //Many thanks to J.Chudoba who noticed this
509 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
510 Short_t dig = digits->GetDigit(it,ip);
511 Int_t idx0=digits->GetTrackID(it,ip,0);
512 Int_t idx1=digits->GetTrackID(it,ip,1);
513 Int_t idx2=digits->GetTrackID(it,ip,2);
514 if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
515 if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
516 if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
517 } while (digits->Next());
518 for (Int_t j=0; j<np; j++) {
520 if (sec>=digp->GetNInnerSector())
521 if (row==nrow_up-1 ) good[j]|=0x4000;
522 if (sec>=digp->GetNInnerSector())
523 if (row==nrow_up-1-gap) good[j]|=0x1000;
525 if (sec>=digp->GetNInnerSector())
526 if (row==nrow_up-1-shift) good[j]|=0x2000;
527 if (sec>=digp->GetNInnerSector())
528 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
540 //****** select tracks which are "good" enough
541 AliStack* stack = rl->Stack();
542 for (i=0; i<np; i++) {
543 if ((good[i]&0x5000) != 0x5000)
544 if ((good[i]&0x2800) != 0x2800) continue;
545 if ((good[i]&0x7FF ) < good_number) continue;
547 TParticle *p = (TParticle*)stack->Particle(i);
549 cerr<<"Can not get particle "<<i<<endl;
552 if (p->Pt()<0.100) continue;
553 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
555 Double_t vx=p->Vx(),vy=p->Vy(),vz=p->Vz();
556 if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue;
557 if (TMath::Abs(vz) > 50.) continue;
559 AliTrackReference *ref=new((*refs)[nt]) AliTrackReference();
562 Int_t pdg=p->GetPdgCode();
563 ref->SetLength(pdg); //This will the particle's PDG !
564 ref->SetMomentum(0.,0.,0.);
565 ref->SetPosition(0.,0.,0.);
570 //**** check if there is also information at the entrance of the TPC
571 TTree *TR=rl->TreeTR();
572 TBranch *branch=TR->GetBranch("TrackReferences");
574 ::Error("GoodTracksTPC","No track references !");
578 TClonesArray tpcdummy("AliTrackReference",1000), *tpcRefs=&tpcdummy;
579 branch->SetAddress(&tpcRefs);
580 Int_t nr=(Int_t)TR->GetEntries();
581 for (Int_t r=0; r<nr; r++) {
582 //cerr<<r<<' '<<nr<<'\r';
585 Int_t nref = tpcRefs->GetEntriesFast();
587 AliTrackReference *tpcRef= 0x0;
588 for (Int_t iref=0; iref<nref; ++iref) {
589 tpcRef = (AliTrackReference*)tpcRefs->UncheckedAt(iref);
590 if (tpcRef->DetectorId() == AliTrackReference::kTPC) break;
594 if (!tpcRef) continue;
597 AliTrackReference *ref=0;
598 for (j=0; j<nt; j++) {
599 ref=(AliTrackReference *)refs->UncheckedAt(j);
600 if (ref->Label()==tpcRef->Label()) break;
603 if (ref==0) continue;
605 ref->SetMomentum(tpcRef->Px(),tpcRef->Py(),tpcRef->Pz());
606 ref->SetPosition(tpcRef->LocalX(),tpcRef->LocalY(),tpcRef->Z());
614 } //****** end of the loop over generated events