]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliCascadeComparison.C
Comparison working with ESD (Yu.Belikov)
[u/mrichter/AliRoot.git] / ITS / AliCascadeComparison.C
CommitLineData
ca28c5f5 1/****************************************************************************
2 * Very important, delicate and rather obscure macro. *
3 * *
4 * Creates list of "findable" cascades, *
5 * calculates efficiency, resolutions etc. *
6 * *
7 * Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr *
8 ****************************************************************************/
9
566abf75 10#if !defined(__CINT__) || defined(__MAKECINT__)
ddae8318 11 #include <Riostream.h>
ca28c5f5 12 #include <fstream.h>
13
566abf75 14 #include "AliRun.h"
5102bab6 15 #include "AliMC.h"
566abf75 16 #include "AliHeader.h"
17 #include "AliRunLoader.h"
18 #include "AliITSLoader.h"
19
ca28c5f5 20 #include "TH1.h"
21 #include "TFile.h"
5102bab6 22 #include "TKey.h"
ca28c5f5 23 #include "TObjArray.h"
24 #include "TStyle.h"
25 #include "TCanvas.h"
26 #include "TLine.h"
27 #include "TText.h"
28 #include "TParticle.h"
29 #include "TStopwatch.h"
ad6e4d9a 30 #include "TPDGCode.h"
ca28c5f5 31
32 #include "AliRun.h"
5102bab6 33 #include "AliESD.h"
ca28c5f5 34 #include "AliCascadeVertex.h"
35#endif
36
37struct GoodCascade {
38 Int_t nlab,plab; // V0's daughter labels
39 Int_t blab; // Bachelor label
40 Int_t code;
41 Float_t px,py,pz;
42 Float_t x,y,z;
43};
44Int_t good_cascades(GoodCascade *gt, Int_t max);
45
566abf75 46extern AliRun *gAlice;
47
ca28c5f5 48Int_t AliCascadeComparison(Int_t code=3312) {
49 //code= 3312; //kXiMinus
50 //code=-3312; //kXiPlusBar
51 //code= 3334; //kOmegaMinus
52 //code=-3334; //kOmegaPlusBar
ca28c5f5 53 cerr<<"Doing comparison...\n";
566abf75 54 TStopwatch timer;
55
5102bab6 56 /*** Check if the file with the "good" cascades exists ***/
57 GoodCascade gc[100];
58 Int_t ngood=0;
59 ifstream in("good_cascades");
60 if (in) {
61 cerr<<"Reading good cascades...\n";
62 while (in>>gc[ngood].nlab>>gc[ngood].plab>>
63 gc[ngood].blab>>gc[ngood].code>>
64 gc[ngood].px>>gc[ngood].py>>gc[ngood].pz>>
65 gc[ngood].x >>gc[ngood].y >>gc[ngood].z) {
66 ngood++;
67 cerr<<ngood<<'\r';
68 if (ngood==100) {
69 cerr<<"Too many good cascades !\n";
70 break;
71 }
72 }
73 if (!in.eof()) cerr<<"Read error (good_cascades) !\n";
74 } else {
75 /*** generate a file with the "good" cascades ***/
76 cerr<<"Marking good cascades (this will take a while)...\n";
77 ngood=good_cascades(gc,100);
78 ofstream out("good_cascades");
79 if (out) {
80 for (Int_t ngd=0; ngd<ngood; ngd++)
81 out<<gc[ngd].nlab<<' '<<gc[ngd].plab<<' '<<
82 gc[ngd].blab<<' '<<gc[ngd].code<<' '<<
83 gc[ngd].px<<' '<<gc[ngd].py<<' '<<gc[ngd].pz<<' '<<
84 gc[ngd].x <<' '<<gc[ngd].y <<' '<<gc[ngd].z <<endl;
85 } else cerr<<"Can not open file (good_cascades) !\n";
86 out.close();
87 }
88
89 AliESD *event=0;
90 { /*** Load reconstructed vertices ***/
91 TFile *ef=TFile::Open("AliESDcas.root");
92 if ((!ef)||(!ef->IsOpen())) {
93 cerr<<"AliCascadeComparison.C: Can't open AliESDcas.root !\n";
94 return 1;
95 }
96 TKey *key=0;
97 TIter next(ef->GetListOfKeys());
98 if ((key=(TKey*)next())!=0) event=(AliESD*)key->ReadObj();
99 ef->Close();
100 }
101 Int_t nentr=event->GetNumberOfCascades();
102
103
ca28c5f5 104 const Double_t cascadeWindow=0.05, cascadeWidth=0.015;
04b2a5f1 105 Double_t ptncut=0.12, ptpcut=0.33, kine0cut=0.003;
106 Double_t ptbcut=0.11, kinecut=0.002;
107 Double_t cascadeMass=1.32131;
ca28c5f5 108 switch (code) {
109 case kXiMinus:
04b2a5f1 110 break;
111 case kXiPlusBar:
112 ptncut=0.33; ptpcut=0.12;
113 break;
ca28c5f5 114 case kOmegaMinus:
04b2a5f1 115 cascadeMass=1.67245;
116 kine0cut=0.001;
117 ptbcut=0.22; kinecut=0.006;
118 break;
119 case kOmegaPlusBar:
120 cascadeMass=1.67245;
121 kine0cut=0.001;
122 ptncut=0.33; ptpcut=0.12;
123 ptbcut=0.22; kinecut=0.006;
124 break;
ca28c5f5 125 default: cerr<<"Invalid PDG code !\n"; return 1;
126 }
127
566abf75 128 /*** create some histograms ***/
129 TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution
130 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
131 TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
132 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
133 TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
134 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
135
136 TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res.
137 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
138 TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.); //y res
139 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
140 TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.); //z res.
141 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
142
143 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
144 TH1F *hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);
145 TH1F *hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
146 TH1F *hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
147 TH1F *hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
148 hg->SetLineColor(4); hg->SetLineWidth(2);
149 TH1F *hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
150 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
151
152 Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
153 TH1F *cs =new TH1F("cs","Cascade Effective Mass",40, mmin, mmax);
154 cs->SetXTitle("(GeV)");
155 cs->SetLineColor(4); cs->SetLineWidth(4);
156 TH1F *csf =new TH1F("csf","Fake Cascade Effective Mass",40, mmin, mmax);
157 csf->SetXTitle("(GeV)"); csf->SetFillColor(6);
158
ca28c5f5 159 Double_t pxg=0.,pyg=0.,ptg=0.;
160 Int_t nlab=-1, plab=-1, blab=-1;
161 Int_t i;
162 for (i=0; i<nentr; i++) {
5102bab6 163 AliESDcascade *cascade=event->GetCascade(i);
164
165 Int_t nidx=TMath::Abs(cascade->GetNindex());
166 Int_t pidx=TMath::Abs(cascade->GetPindex());
167 Int_t bidx=TMath::Abs(cascade->GetBindex());
168
169 AliESDtrack *ntrack=event->GetTrack(nidx);
170 AliESDtrack *ptrack=event->GetTrack(pidx);
171 AliESDtrack *btrack=event->GetTrack(bidx);
172
173 nlab=TMath::Abs(ntrack->GetLabel());
174 plab=TMath::Abs(ptrack->GetLabel());
175 blab=TMath::Abs(btrack->GetLabel());
04b2a5f1 176
177 /** Kinematical cuts **/
178 Double_t pxn,pyn,pzn; cascade->GetNPxPyPz(pxn,pyn,pzn);
179 Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
180 if (ptn < ptncut) continue;
181 Double_t pxp,pyp,pzp; cascade->GetPPxPyPz(pxp,pyp,pzp);
182 Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
183 if (ptp < ptpcut) continue;
184 Double_t pxb,pyb,pzb; cascade->GetBPxPyPz(pxb,pyb,pzb);
185 Double_t ptb=TMath::Sqrt(pxb*pxb + pyb*pyb);
186 if (ptb < ptbcut) continue;
187 Double_t kine0;
188 Double_t kine=cascade->ChangeMassHypothesis(kine0,code);
189 if (TMath::Abs(kine0)>kine0cut) continue;
190 //if (TMath::Abs(kine)>kinecut) continue;
ca28c5f5 191
192 Double_t mass=cascade->GetEffMass();
193 cs->Fill(mass);
04b2a5f1 194 csf->Fill(mass);
ca28c5f5 195
196 if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
197
198 Int_t j;
199 for (j=0; j<ngood; j++) {
200 if (gc[j].code != cascade->GetPdgCode()) continue;
201 if (gc[j].nlab == nlab)
202 if (gc[j].plab == plab)
203 if (gc[j].blab == blab) break;
204 }
205
206 Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
207 Double_t pt=TMath::Sqrt(px*px+py*py);
208
209 if (j==ngood) {
210 hfake->Fill(pt);
211 cerr<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
212 continue;
213 }
04b2a5f1 214 csf->Fill(mass,-1);
ca28c5f5 215
216 pxg=gc[j].px; pyg=gc[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
217 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
218 Double_t lamg=TMath::ATan2(gc[j].pz,ptg), lam=TMath::ATan2(pz,pt);
219 hp->Fill((phi - phig)*1000.);
220 hl->Fill((lam - lamg)*1000.);
221 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
222
223 Double_t x,y,z; cascade->GetXYZ(x,y,z);
224 hx->Fill((x-gc[j].x)*10);
225 hy->Fill((y-gc[j].y)*10);
226 hz->Fill((z-gc[j].z)*10);
227
228 hfound->Fill(ptg);
229 gc[j].nlab=-1;
230
231 }
232 for (i=0; i<ngood; i++) {
233 if (gc[i].code != code) continue;
234 pxg=gc[i].px; pyg=gc[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
235 hgood->Fill(ptg);
236 nlab=gc[i].nlab; plab=gc[i].plab; blab=gc[i].blab;
237 if (nlab < 0) continue;
238 cerr<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
239 }
240
5102bab6 241 delete event;
ca28c5f5 242
243 Stat_t ng=hgood->GetEntries();
244 Stat_t nf=hfound->GetEntries();
245
246 cerr<<"Number of found cascades: "<<nentr<<" ("<<nf<<" in the peak)\n";
247 cerr<<"Number of good cascades: "<<ng<<endl;
248
249 if (ng!=0)
250 cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
251
252 gStyle->SetOptStat(111110);
253 gStyle->SetOptFit(1);
254
255 TCanvas *c1=new TCanvas("c1","",0,0,580,610);
256 c1->Divide(2,2);
257
258 c1->cd(1);
259 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
260 //hp->Fit("gaus");
261 hp->Draw();
262 hl->Draw("same"); c1->cd();
263
264 c1->cd(2);
265 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
266 //hpt->Fit("gaus"); c1->cd();
267 hpt->Draw(); c1->cd();
268
269 c1->cd(3);
270 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
271 //hx->Fit("gaus");
272 hx->Draw();
273 hy->Draw("same"); c1->cd();
274
275 c1->cd(4);
276 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
277 //hz->Fit("gaus");
278 hz->Draw();
279
280
281 TCanvas *c2=new TCanvas("c2","",600,0,580,610);
282 c2->Divide(1,2);
283
284 c2->cd(1);
285 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
286 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
287 hg->Divide(hfound,hgood,1,1.,"b");
288 hf->Divide(hfake,hgood,1,1.,"b");
289 hg->SetMaximum(1.4);
290 hg->SetYTitle("Cascade reconstruction efficiency");
291 hg->SetXTitle("Pt (GeV/c)");
292 hg->Draw();
293
294 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
295 line1->Draw("same");
296 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
297 line2->Draw("same");
298
299 hf->SetFillColor(1);
300 hf->SetFillStyle(3013);
301 hf->SetLineColor(2);
302 hf->SetLineWidth(2);
303 hf->Draw("histsame");
304 TText *text = new TText(0.461176,0.248448,"Fake cascades");
305 text->SetTextSize(0.05);
306 text->Draw();
307 text = new TText(0.453919,1.11408,"Good cascades");
308 text->SetTextSize(0.05);
309 text->Draw();
310
311
312 c2->cd(2);
313 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
ca28c5f5 314 //cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
315 cs->Draw();
04b2a5f1 316 csf->Draw("same");
ca28c5f5 317 Double_t max=cs->GetMaximum();
318 TLine *line3 =
319 new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
320 line3->Draw("same");
321 TLine *line4 =
322 new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
323 line4->Draw("same");
324
325 timer.Stop(); timer.Print();
326
327 return 0;
328}
329
330
ca28c5f5 331Int_t good_cascades(GoodCascade *gc, Int_t max) {
332 Int_t nc=0;
333 /*** Get information about the cuts ***/
334 Double_t r2min=0.9*0.9;
335 Double_t r2max=2.9*2.9;
336
337 /*** Get labels of the "good" tracks ***/
338 Double_t dd; Int_t id, label[15000], ngt=0;
339 ifstream in("good_tracks_its");
340 if (!in) {
341 cerr<<"Can't open the file good_tracks_its \n";
342 return nc;
343 }
344 while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
345 ngt++;
346 if (ngt>=15000) {
347 cerr<<"Too many good ITS tracks !\n";
348 return nc;
349 }
350 }
351 if (!in.eof()) {
352 cerr<<"Read error (good_tracks_its) !\n";
353 return nc;
354 }
355
356 /*** Get an access to the kinematics ***/
5102bab6 357 if (gAlice) {
358 delete gAlice->GetRunLoader();
359 delete gAlice;
360 gAlice=0;
361 }
362 AliRunLoader *rl = AliRunLoader::Open("galice.root");
363 if (!rl) {
364 cerr<<"AliCascadeComparison.C, good_cascades :Can't start session !\n";
566abf75 365 return 1;
ca28c5f5 366 }
5102bab6 367
566abf75 368 rl->LoadgAlice();
369 rl->LoadHeader();
370 rl->LoadKinematics();
371 Int_t np = rl->GetHeader()->GetNtrack();
ca28c5f5 372
ca28c5f5 373 while (np--) {
374 cerr<<np<<'\r';
0555741a 375 TParticle *cp=gAlice->GetMCApp()->Particle(np);
ca28c5f5 376
377 /*** only these cascades are "good" ***/
378 Int_t code=cp->GetPdgCode();
379 if (code!=kXiMinus) if (code!=kXiPlusBar)
380 if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue;
381
382 /*** daughter tracks must be "good" ***/
383 Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
384 if (v0lab==blab) continue;
385 if (v0lab<0) continue;
386 if (blab<0) continue;
387
0555741a 388 TParticle *p0=gAlice->GetMCApp()->Particle(v0lab);
389 TParticle *bp=gAlice->GetMCApp()->Particle(blab);
ca28c5f5 390 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
391 TParticle *p=p0; p0=bp; bp=p;
392 Int_t i=v0lab; v0lab=blab; blab=i;
393 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar))
394 continue;
395 }
396
397 /** is the bachelor "good" ? **/
398 Int_t i;
399 for (i=0; i<ngt; i++) if (label[i]==blab) break;
400 if (i==ngt) continue;
401
402 /** is the V0 "good" ? **/
403 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
404 if (nlab==plab) continue;
405 if (nlab<0) continue;
406 if (plab<0) continue;
407
408 for (i=0; i<ngt; i++) if (label[i]==nlab) break;
409 if (i==ngt) continue;
410 for (i=0; i<ngt; i++) if (label[i]==plab) break;
411 if (i==ngt) continue;
412
413 /*** fiducial volume ***/
414 Double_t x=bp->Vx(), y=bp->Vy(), r2=x*x+y*y; //bachelor
415 if (r2<r2min) continue;
416 if (r2>r2max) continue;
0555741a 417 TParticle *pp=gAlice->GetMCApp()->Particle(plab);
ca28c5f5 418 x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y; //V0
419 if (r2<r2min) continue;
420 if (r2>r2max) continue;
421
0555741a 422 if (gAlice->GetMCApp()->Particle(plab)->GetPDG()->Charge() < 0.) {
ca28c5f5 423 i=plab; plab=nlab; nlab=i;
424 }
425
426 gc[nc].code=code;
427 gc[nc].plab=plab; gc[nc].nlab=nlab; gc[nc].blab=blab;
428 gc[nc].px=cp->Px(); gc[nc].py=cp->Py(); gc[nc].pz=cp->Pz();
429 gc[nc].x=bp->Vx(); gc[nc].y=bp->Vy(); gc[nc].z=bp->Vz();
430 nc++;
431
432 }
433
566abf75 434 delete rl;
ca28c5f5 435
436 return nc;
437}