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__)
11 #include "Riostream.h"
16 #include "AliHeader.h"
17 #include "AliRunLoader.h"
18 #include "AliITSLoader.h"
24 #include "TObjArray.h"
29 #include "TParticle.h"
30 #include "TStopwatch.h"
34 #include "AliV0vertex.h"
43 Int_t good_vertices(GoodVertex *gt, Int_t max);
45 extern AliRun *gAlice;
47 Int_t AliV0Comparison(Int_t code=310) { //Lambda=3122, LambdaBar=-3122
48 cerr<<"Doing comparison...\n";
51 /*** Check if the file with the "good" vertices exists ***/
54 ifstream in("good_vertices");
56 cerr<<"Reading good vertices...\n";
57 while (in>>gv[ngood].nlab>>gv[ngood].plab>>gv[ngood].code>>
58 gv[ngood].px>>gv[ngood].py>>gv[ngood].pz>>
59 gv[ngood].x >>gv[ngood].y >>gv[ngood].z) {
63 cerr<<"Too many good vertices !\n";
67 if (!in.eof()) cerr<<"Read error (good_vertices) !\n";
69 /*** generate a file with the "good" vertices ***/
70 cerr<<"Marking good vertices (this will take a while)...\n";
71 ngood=good_vertices(gv,1000);
72 ofstream out("good_vertices");
74 for (Int_t ngd=0; ngd<ngood; ngd++)
75 out<<gv[ngd].nlab<<' '<<gv[ngd].plab<<' '<<gv[ngd].code<<' '<<
76 gv[ngd].px<<' '<<gv[ngd].py<<' '<<gv[ngd].pz<<' '<<
77 gv[ngd].x <<' '<<gv[ngd].y <<' '<<gv[ngd].z <<endl;
78 } else cerr<<"Can not open file (good_vertices) !\n";
84 { /*** Load reconstructed vertices ***/
85 TFile *ef=TFile::Open("AliESDv0.root");
86 if ((!ef)||(!ef->IsOpen())) {
87 cerr<<"AliV0Comparison.C: Can't open AliESDv0.root !\n";
91 TIter next(ef->GetListOfKeys());
92 if ((key=(TKey*)next())!=0) event=(AliESD*)key->ReadObj();
95 Int_t nentr=event->GetNumberOfV0s();
98 const Double_t V0window=0.05;
99 Double_t ptncut=0.13, ptpcut=0.13, kinecut=0.03;
100 Double_t V0mass=0.497672, V0width=0.020;
105 V0mass=1.115683; V0width=0.015; ptpcut=0.50; kinecut=0.002;
108 V0mass=1.115683; V0width=0.015; ptncut=0.50; kinecut=0.002;
110 default: cerr<<"Invalid PDG code !\n"; return 1;
113 TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution
114 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
115 TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
116 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
117 TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
118 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
120 TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res.
121 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
122 TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.); //y res
123 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
124 TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.); //z res.
125 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
128 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
129 TH1F *hgood=new TH1F("hgood","Good Vertices",nchan,pmin,pmax);
130 TH1F *hfound=new TH1F("hfound","Found Vertices",nchan,pmin,pmax);
131 TH1F *hfake=new TH1F("hfake","Fake Vertices",nchan,pmin,pmax);
132 TH1F *hg=new TH1F("hg","Efficiency for Good Vertices",nchan,pmin,pmax);
133 hg->SetLineColor(4); hg->SetLineWidth(2);
134 TH1F *hf=new TH1F("hf","Probability of Fake Vertices",nchan,pmin,pmax);
135 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
137 Double_t mmin=V0mass-V0window, mmax=V0mass+V0window;
138 TH1F *v0s =new TH1F("v0s","V0s Effective Mass",40, mmin, mmax);
139 v0s->SetXTitle("(GeV)");
140 v0s->SetLineColor(4); v0s->SetLineWidth(4);
141 TH1F *v0sf =new TH1F("v0sf","Fake V0s Effective Mass",40, mmin, mmax);
142 v0sf->SetXTitle("(GeV)"); v0sf->SetFillColor(6);
145 Double_t pxg=0.,pyg=0.,ptg=0.;
146 Int_t nlab=-1, plab=-1;
148 for (i=0; i<nentr; i++) {
149 AliESDv0 *vertex=event->GetV0(i);
151 Int_t nidx=TMath::Abs(vertex->GetNindex());
152 Int_t pidx=TMath::Abs(vertex->GetPindex());
154 AliESDtrack *ntrack=event->GetTrack(nidx);
155 AliESDtrack *ptrack=event->GetTrack(pidx);
157 nlab=TMath::Abs(ntrack->GetLabel());
158 plab=TMath::Abs(ptrack->GetLabel());
160 /** Kinematical cuts **/
161 Double_t pxn,pyn,pzn; vertex->GetNPxPyPz(pxn,pyn,pzn);
162 Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
163 if (ptn < ptncut) continue;
164 Double_t pxp,pyp,pzp; vertex->GetPPxPyPz(pxp,pyp,pzp);
165 Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
166 if (ptp < ptpcut) continue;
167 Double_t kine=vertex->ChangeMassHypothesis(code);
168 //if (TMath::Abs(kine)>kinecut) continue;
170 Double_t mass=vertex->GetEffMass();
175 for (j=0; j<ngood; j++) {
176 if (gv[j].code != vertex->GetPdgCode()) continue;
177 if (gv[j].nlab == nlab)
178 if (gv[j].plab == plab) break;
181 if (TMath::Abs(mass-V0mass)>V0width) continue;
183 Double_t px,py,pz; vertex->GetPxPyPz(px,py,pz);
184 Double_t pt=TMath::Sqrt(px*px+py*py);
188 cerr<<"Fake vertex: ("<<nlab<<","<<plab<<")\n";
193 pxg=gv[j].px; pyg=gv[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
194 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
195 Double_t lamg=TMath::ATan2(gv[j].pz,ptg), lam=TMath::ATan2(pz,pt);
196 hp->Fill((phi - phig)*1000.);
197 hl->Fill((lam - lamg)*1000.);
198 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
200 Double_t x,y,z; vertex->GetXYZ(x,y,z);
201 hx->Fill((x-gv[j].x)*10);
202 hy->Fill((y-gv[j].y)*10);
203 hz->Fill((z-gv[j].z)*10);
209 for (i=0; i<ngood; i++) {
210 if (gv[i].code != code) continue;
211 pxg=gv[i].px; pyg=gv[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
213 nlab=gv[i].nlab; plab=gv[i].plab;
214 if (nlab < 0) continue;
215 cerr<<"Vertex ("<<nlab<<','<<plab<<") has not been found !\n";
220 Stat_t ng=hgood->GetEntries();
221 Stat_t nf=hfound->GetEntries();
223 cerr<<"Number of found vertices: "<<nentr<<" ("<<nf<<" in the peak)\n";
224 cerr<<"Number of good vertices: "<<ng<<endl;
227 cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
229 gStyle->SetOptStat(111110);
230 gStyle->SetOptFit(1);
232 TCanvas *c1=new TCanvas("c1","",0,0,580,610);
236 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
239 hl->Draw("same"); c1->cd();
242 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
243 //hpt->Fit("gaus"); c1->cd();
244 hpt->Draw(); c1->cd();
247 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
250 hy->Draw("same"); c1->cd();
253 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
258 TCanvas *c2=new TCanvas("c2","",600,0,580,610);
262 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
263 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
264 hg->Divide(hfound,hgood,1,1.,"b");
265 hf->Divide(hfake,hgood,1,1.,"b");
267 hg->SetYTitle("V0 reconstruction efficiency");
268 hg->SetXTitle("Pt (GeV/c)");
271 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
273 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
277 hf->SetFillStyle(3013);
280 hf->Draw("histsame");
281 TText *text = new TText(0.461176,0.248448,"Fake vertices");
282 text->SetTextSize(0.05);
284 text = new TText(0.453919,1.11408,"Good vertices");
285 text->SetTextSize(0.05);
290 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
291 //v0s->Fit("gaus","","",V0mass-V0width,V0mass+V0width);
294 Double_t max=v0s->GetMaximum();
295 TLine *line3 = new TLine(V0mass-V0width,0.,V0mass-V0width,max);
297 TLine *line4 = new TLine(V0mass+V0width,0.,V0mass+V0width,max);
300 timer.Stop(); timer.Print();
305 Int_t good_vertices(GoodVertex *gv, Int_t max) {
307 /*** Get information about the cuts ***/
308 Double_t r2min=0.9*0.9;
309 Double_t r2max=2.9*2.9;
311 /*** Get labels of the "good" tracks ***/
312 Double_t dd; Int_t id, label[15000], ngt=0;
313 ifstream in("good_tracks_its");
315 cerr<<"Can't open the file good_tracks_its \n";
318 while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
321 cerr<<"Too many good ITS tracks !\n";
326 cerr<<"Read error (good_tracks_its) !\n";
330 /*** Get an access to the kinematics ***/
332 delete gAlice->GetRunLoader();
336 AliRunLoader *rl = AliRunLoader::Open("galice.root");
338 cerr<<"AliV0Comparison.C, good_vertices :Can't start session !\n";
344 rl->LoadKinematics();
345 Int_t np = rl->GetHeader()->GetNtrack();
349 TParticle *p0=gAlice->GetMCApp()->Particle(np);
351 /*** only these V0s are "good" ***/
352 Int_t code=p0->GetPdgCode();
355 if (code!=kLambda0Bar) continue;
357 /*** daughter tracks should be "good" ***/
358 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
359 if (nlab==plab) continue;
360 if (nlab<0) continue;
361 if (plab<0) continue;
363 for (i=0; i<ngt; i++) if (label[i]==nlab) break;
364 if (i==ngt) continue;
365 for (i=0; i<ngt; i++) if (label[i]==plab) break;
366 if (i==ngt) continue;
368 /*** fiducial volume ***/
369 TParticle *p=gAlice->GetMCApp()->Particle(nlab);
370 Double_t x=p->Vx(), y=p->Vy(), z=p->Vz(), r2=x*x+y*y;
371 if (r2<r2min) continue;
372 if (r2>r2max) continue;
374 if (gAlice->GetMCApp()->Particle(plab)->GetPDG()->Charge() < 0.) {
375 i=plab; plab=nlab; nlab=i;
379 gv[nv].plab=plab; gv[nv].nlab=nlab;
380 gv[nv].px=p0->Px(); gv[nv].py=p0->Py(); gv[nv].pz=p0->Pz();
381 gv[nv].x=x; gv[nv].y=y; gv[nv].z=z;