1 /****************************************************************************
2 * Very important, delicate and rather obscure macro. *
4 * Creates list of "findable" V0s, *
5 * calculates efficiency, resolutions etc. *
7 * Origin: I.Belikov, IReS, Strasbourg, Jouri.Belikov@cern.ch *
8 ****************************************************************************/
10 #if !defined(__CINT__) || defined(__MAKECINT__)
13 #include <Riostream.h>
16 #include <TParticle.h>
20 #include <TBenchmark.h>
26 #include "AliHeader.h"
27 #include "AliTrackReference.h"
28 #include "AliRunLoader.h"
33 Int_t GoodV0s(const Char_t *dir=".");
35 extern AliRun *gAlice;
36 extern TBenchmark *gBenchmark;
39 static Int_t allgood=0;
40 static Int_t allfound=0;
42 Int_t AliV0Comparison(Int_t code=310, const Char_t *dir=".") {
43 //Lambda=3122, LambdaBar=-3122
44 gBenchmark->Start("AliV0Comparison");
46 ::Info("AliV0Comparison.C","Doing comparison...");
49 const Double_t V0window=0.05;
50 Double_t ptncut=0.13, ptpcut=0.13, kinecut=0.03;
51 Double_t V0mass=0.497672, V0width=0.020;
56 V0mass=1.115683; V0width=0.015; ptpcut=0.50; kinecut=0.002;
59 V0mass=1.115683; V0width=0.015; ptncut=0.50; kinecut=0.002;
61 default: cerr<<"Invalid PDG code !\n"; return 1;
64 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
65 if (!hp) hp=new TH1F("hp","PHI Resolution",30,-30.,30.);
66 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
68 TH1F *hl=(TH1F*)gROOT->FindObject("hl");
69 if (!hl) hl=new TH1F("hl","LAMBDA Resolution",30,-30,30);
70 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
72 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
73 if (!hpt) hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
74 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
76 TH1F *hx=(TH1F*)gROOT->FindObject("hx");
77 if (!hx) hx=new TH1F("hx","Position Resolution (X)",30,-3.,3.);
78 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
80 TH1F *hy=(TH1F*)gROOT->FindObject("hy");
81 if (!hy) hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.);
82 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
84 TH1F *hz=(TH1F*)gROOT->FindObject("hz");
85 if (!hz) hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.);
86 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
89 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
90 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
91 if (!hgood) hgood=new TH1F("hgood","Good Vertices",nchan,pmin,pmax);
93 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
94 if (!hfound) hfound=new TH1F("hfound","Found Vertices",nchan,pmin,pmax);
96 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
97 if (!hfake) hfake=new TH1F("hfake","Fake Vertices",nchan,pmin,pmax);
99 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
100 if (!hg) hg=new TH1F("hg","Efficiency for Good Vertices",nchan,pmin,pmax);
101 hg->SetLineColor(4); hg->SetLineWidth(2);
103 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
104 if (!hf) hf=new TH1F("hf","Probability of Fake Vertices",nchan,pmin,pmax);
105 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
108 Double_t mmin=V0mass-V0window, mmax=V0mass+V0window;
109 TH1F *v0s=(TH1F*)gROOT->FindObject("v0s");
110 if (!v0s) v0s=new TH1F("v0s","V0s Effective Mass",40, mmin, mmax);
111 v0s->SetXTitle("(GeV)"); v0s->SetLineColor(4); v0s->SetLineWidth(4);
113 TH1F *v0sf=(TH1F*)gROOT->FindObject("v0sf");
114 if (!v0sf) v0sf=new TH1F("v0sf","Fake V0s Effective Mass",40, mmin, mmax);
115 v0sf->SetXTitle("(GeV)"); v0sf->SetFillColor(6);
119 sprintf(fname,"%s/GoodV0s.root",dir);
121 TFile *refFile=TFile::Open(fname,"old");
122 if (!refFile || !refFile->IsOpen()) {
123 ::Info("AliV0Comparison.C","Marking good V0s (will take a while)...");
125 ::Error("AliV0Comparison.C","Can't generate the reference file !");
129 refFile=TFile::Open(fname,"old");
130 if (!refFile || !refFile->IsOpen()) {
131 ::Error("AliV0Comparison.C","Can't open the reference file !");
135 TTree *v0Tree=(TTree*)refFile->Get("v0Tree");
137 ::Error("AliV0Comparison.C","Can't get the reference tree !");
140 TBranch *pbranch=v0Tree->GetBranch("positive");
142 ::Error("AliV0Comparison.C","Can't get the positive daughter branch !");
145 TClonesArray dummy("AliTrackReference",1000), *prefs=&dummy;
146 pbranch->SetAddress(&prefs);
148 TBranch *nbranch=v0Tree->GetBranch("negative");
150 ::Error("AliV0Comparison.C","Can't get the negative daughter branch !");
153 TClonesArray dumm("AliTrackReference",1000), *nrefs=&dumm;
154 nbranch->SetAddress(&nrefs);
157 sprintf(fname,"%s/AliESDs.root",dir);
158 TFile *ef=TFile::Open(fname);
159 if ((!ef)||(!ef->IsOpen())) {
160 sprintf(fname,"%s/AliESDv0.root",dir);
161 ef=TFile::Open(fname);
162 if ((!ef)||(!ef->IsOpen())) {
163 ::Error("AliV0Comparison.C","Can't open AliESDv0.root !");
167 AliESD* event = new AliESD;
168 TTree* esdTree = (TTree*) ef->Get("esdTree");
170 ::Error("AliV0Comparison.C", "no ESD tree found");
173 esdTree->SetBranchAddress("ESD", &event);
176 //******* Loop over events *********
178 while (esdTree->GetEvent(e)) {
179 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
181 Int_t nentr=event->GetNumberOfV0s();
184 if (v0Tree->GetEvent(e++)==0) {
185 cerr<<"No reconstructable V0s !\n";
189 Int_t ngood=prefs->GetEntriesFast(),ng=0;
191 Double_t pxg=0.,pyg=0.,pzg=0.,ptg=0.;
192 Int_t nlab=-1, plab=-1;
194 for (i=0; i<nentr; i++) {
195 AliESDv0 *vertex=event->GetV0(i);
197 Int_t nidx=TMath::Abs(vertex->GetNindex());
198 Int_t pidx=TMath::Abs(vertex->GetPindex());
200 AliESDtrack *ntrack=event->GetTrack(nidx);
201 AliESDtrack *ptrack=event->GetTrack(pidx);
203 nlab=TMath::Abs(ntrack->GetLabel());
204 plab=TMath::Abs(ptrack->GetLabel());
206 /** Kinematical cuts **/
207 Double_t pxn,pyn,pzn; vertex->GetNPxPyPz(pxn,pyn,pzn);
208 Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
209 if (ptn < ptncut) continue;
210 Double_t pxp,pyp,pzp; vertex->GetPPxPyPz(pxp,pyp,pzp);
211 Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
212 if (ptp < ptpcut) continue;
213 Double_t kine=vertex->ChangeMassHypothesis(code);
214 //if (TMath::Abs(kine)>kinecut) continue;
216 Double_t mass=vertex->GetEffMass();
220 AliTrackReference *nref=0, *pref=0;
222 for (j=0; j<ngood; j++) {
223 nref=(AliTrackReference*)nrefs->UncheckedAt(j);
224 pref=(AliTrackReference*)prefs->UncheckedAt(j);
225 if (nref->Label() == nlab)
226 if (pref->Label() == plab) break;
229 if (TMath::Abs(mass-V0mass)>V0width) continue;
231 Double_t px,py,pz; vertex->GetPxPyPz(px,py,pz);
232 Double_t pt=TMath::Sqrt(px*px+py*py);
236 cout<<"Fake vertex: ("<<nlab<<","<<plab<<")\n";
241 pxg=nref->Px()+pref->Px(); pyg=nref->Py()+pref->Py();
242 pzg=nref->Pz()+pref->Pz();
243 ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
244 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
245 Double_t lamg=TMath::ATan2(pzg,ptg), lam=TMath::ATan2(pz,pt);
246 hp->Fill((phi - phig)*1000.);
247 hl->Fill((lam - lamg)*1000.);
248 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
250 Double_t x,y,z; vertex->GetXYZ(x,y,z);
251 hx->Fill((x - nref->X())*10);
252 hy->Fill((y - nref->Y())*10);
253 hz->Fill((z - nref->Z())*10);
259 for (i=0; i<ngood; i++) {
260 AliTrackReference *nref=(AliTrackReference*)nrefs->UncheckedAt(i);
261 AliTrackReference *pref=(AliTrackReference*)prefs->UncheckedAt(i);
262 Int_t pdg=(Int_t)nref->GetLength();//this is the mother's PDG !
263 if (code!=pdg) continue;
265 pxg=nref->Px()+pref->Px(); pyg=nref->Py()+pref->Py();
266 ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
268 nlab=nref->Label(); plab=pref->Label();
269 if (nlab < 0) continue;
270 cout<<"Vertex ("<<nlab<<','<<plab<<") has not been found !\n";
274 cout<<"Number of found V0s : "<<nentr<<endl;
275 cout<<"Number of \"good\" V0s : "<<ng<<endl;
280 } //**** End of the loop over events
288 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
289 if (ng!=0) cout<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
290 cout<<"Total number of found V0s: "<<allfound<<" ("<<nf<<" in the peak)\n";
291 cout<<"Total number of \"good\" V0s: "<<allgood<<endl;
293 gStyle->SetOptStat(111110);
294 gStyle->SetOptFit(1);
296 TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
298 c1=new TCanvas("c1","",0,0,580,610);
305 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
306 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus");
307 hl->Draw("same"); c1->cd();
310 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
311 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus");
315 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
316 if (hx->GetEntries()<minc) hx->Draw(); else hx->Fit("gaus");
317 hy->Draw("same"); c1->cd();
320 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
321 if (hz->GetEntries()<minc) hz->Draw(); else hz->Fit("gaus");
325 TCanvas *c2=(TCanvas*)gROOT->FindObject("c2");
327 c2=new TCanvas("c2","",600,0,580,610);
332 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
333 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
334 hg->Divide(hfound,hgood,1,1.,"b");
335 hf->Divide(hfake,hgood,1,1.,"b");
337 hg->SetYTitle("V0 reconstruction efficiency");
338 hg->SetXTitle("Pt (GeV/c)");
341 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
343 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
347 hf->SetFillStyle(3013);
350 hf->Draw("histsame");
351 TText *text = new TText(0.461176,0.248448,"Fake vertices");
352 text->SetTextSize(0.05);
354 text = new TText(0.453919,1.11408,"Good vertices");
355 text->SetTextSize(0.05);
360 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
361 if (v0s->GetEntries()<minc) v0s->Draw();
362 else v0s->Fit("gaus","","",V0mass-V0width,V0mass+V0width);
364 Double_t max=v0s->GetMaximum();
365 TLine *line3 = new TLine(V0mass-V0width,0.,V0mass-V0width,max);
367 TLine *line4 = new TLine(V0mass+V0width,0.,V0mass+V0width,max);
372 TFile fc("AliV0Comparison.root","RECREATE");
377 gBenchmark->Stop("AliV0Comparison");
378 gBenchmark->Show("AliV0Comparison");
384 Int_t GoodV0s(const Char_t *dir) {
386 delete gAlice->GetRunLoader();
387 delete gAlice;//if everything was OK here it is already NULL
392 sprintf(fname,"%s/galice.root",dir);
394 AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
396 ::Error("GoodTracksITS","Can't start session !");
402 rl->LoadKinematics();
405 Int_t nev=rl->GetNumberOfEvents();
406 ::Info("GoodV0s","Number of events : %d\n",nev);
408 sprintf(fname,"%s/GoodTracksITS.root",dir);
409 TFile *itsFile=TFile::Open(fname);
410 if ((!itsFile)||(!itsFile->IsOpen())) {
411 ::Error("GoodV0s","Can't open the GoodTracksITS.root !");
415 TClonesArray dum("AliTrackReference",1000), *itsRefs=&dum;
416 TTree *itsTree=(TTree*)itsFile->Get("itsTree");
418 ::Error("GoodV0s","Can't get the ITS reference tree !");
422 TBranch *itsBranch=itsTree->GetBranch("ITS");
424 ::Error("GoodV0s","Can't get the ITS reference branch !");
428 itsBranch->SetAddress(&itsRefs);
431 sprintf(fname,"%s/GoodV0s.root",dir);
432 TFile *v0File=TFile::Open(fname,"recreate");
433 TClonesArray dummy("AliTrackReference",1000), *nrefs=&dummy;
434 TClonesArray dumm("AliTrackReference",1000), *prefs=&dumm;
435 TTree v0Tree("v0Tree","Tree with info about the reconstructable V0s");
436 v0Tree.Branch("negative",&nrefs);
437 v0Tree.Branch("positive",&prefs);
440 /*** Some information about the cuts ***/
441 Double_t r2min=0.9*0.9;
442 Double_t r2max=2.9*2.9;
446 //******** Loop over generated events
447 for (Int_t e=0; e<nev; e++) {
448 rl->GetEvent(e); v0File->cd();
450 Int_t np = rl->GetHeader()->GetNtrack();
451 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
453 itsTree->GetEvent(e);
454 Int_t nk=itsRefs->GetEntriesFast();
456 AliStack *stack=rl->Stack();
458 AliTrackReference *nref=0, *pref=0;
463 TParticle *p0=stack->Particle(np);
465 /*** only these V0s are "good" ***/
466 Int_t code=p0->GetPdgCode();
469 if (code!=kLambda0Bar) continue;
471 /*** daughter tracks should be "good" ***/
472 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
473 if (nlab==plab) continue;
474 if (nlab<0) continue;
475 if (plab<0) continue;
477 if (stack->Particle(plab)->GetPDG()->Charge() < 0.) {
478 i=plab; plab=nlab; nlab=i;
481 for (i=0; i<nk; i++) {
482 nref=(AliTrackReference*)itsRefs->UncheckedAt(i);
483 if (nref->Label()==nlab) break;
486 for (i=0; i<nk; i++) {
487 pref=(AliTrackReference*)itsRefs->UncheckedAt(i);
488 if (pref->Label()==plab) break;
492 /*** fiducial volume ***/
493 TParticle *p=stack->Particle(nlab);
494 Double_t x=p->Vx(), y=p->Vy(), z=p->Vz(), r2=x*x+y*y;
495 if (r2<r2min) continue;
496 if (r2>r2max) continue;
498 Int_t pdg=p0->GetPdgCode();
499 nref->SetLength(pdg); //This will the V0's PDG !
501 new((*nrefs)[nv]) AliTrackReference(*nref);
502 new((*prefs)[nv]) AliTrackReference(*pref);
509 nrefs->Clear(); prefs->Clear();
511 } //**** end of the loop over generated events