]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliCascadeComparison.C
fixxed def. of GetList.
[u/mrichter/AliRoot.git] / ITS / AliCascadeComparison.C
CommitLineData
a9a2d814 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
10#ifndef __CINT__
11 #include <iostream.h>
12 #include <fstream.h>
13
14 #include "TH1.h"
15 #include "TFile.h"
16 #include "TTree.h"
17 #include "TObjArray.h"
18 #include "TStyle.h"
19 #include "TCanvas.h"
20 #include "TLine.h"
21 #include "TText.h"
22 #include "TParticle.h"
23 #include "TStopwatch.h"
24
25 #include "AliRun.h"
26 #include "AliPDG.h"
27 #include "AliCascadeVertex.h"
28#endif
29
30struct GoodCascade {
31 Int_t nlab,plab; // V0's daughter labels
32 Int_t blab; // Bachelor label
33 Int_t code;
34 Float_t px,py,pz;
35 Float_t x,y,z;
36};
37Int_t good_cascades(GoodCascade *gt, Int_t max);
38
39Int_t AliCascadeComparison(Int_t code=3312) {
40 //code= 3312; //kXiMinus
41 //code=-3312; //kXiPlusBar
42 //code= 3334; //kOmegaMinus
43 //code=-3334; //kOmegaPlusBar
44
45 cerr<<"Doing comparison...\n";
46
47 const Double_t cascadeWindow=0.05, cascadeWidth=0.015;
48 Double_t cascadeMass=0.5;
49 switch (code) {
50 case kXiMinus:
51 case kXiPlusBar: cascadeMass=1.32131; break;
52 case kOmegaMinus:
53 case kOmegaPlusBar: cascadeMass=1.67245; break;
54 default: cerr<<"Invalid PDG code !\n"; return 1;
55 }
56
57 TStopwatch timer;
58
59 /*** Load reconstructed cascades ***/
60 TFile *cf=TFile::Open("AliCascadeVertices.root");
61 if (!cf->IsOpen()){cerr<<"Can't open AliCascadeVertices.root !\n";return 2;}
62 TObjArray carray(1000);
63 TTree *cTree=(TTree*)cf->Get("TreeCasc");
64 TBranch *branch=cTree->GetBranch("cascades");
65 Int_t nentr=(Int_t)cTree->GetEntries();
66 for (Int_t i=0; i<nentr; i++) {
67 AliCascadeVertex *iovertex=new AliCascadeVertex;
68 branch->SetAddress(&iovertex);
69 cTree->GetEvent(i);
70 carray.AddLast(iovertex);
71 }
72 cf->Close();
73
74 /*** Check if the file with the "good" cascades exists ***/
75 GoodCascade gc[100];
76 Int_t ngood=0;
77 ifstream in("good_cascades");
78 if (in) {
79 cerr<<"Reading good cascades...\n";
80 while (in>>gc[ngood].nlab>>gc[ngood].plab>>
81 gc[ngood].blab>>gc[ngood].code>>
82 gc[ngood].px>>gc[ngood].py>>gc[ngood].pz>>
83 gc[ngood].x >>gc[ngood].y >>gc[ngood].z) {
84 ngood++;
85 cerr<<ngood<<'\r';
86 if (ngood==100) {
87 cerr<<"Too many good cascades !\n";
88 break;
89 }
90 }
91 if (!in.eof()) cerr<<"Read error (good_cascades) !\n";
92 } else {
93 /*** generate a file with the "good" cascades ***/
94 cerr<<"Marking good cascades (this will take a while)...\n";
95 ngood=good_cascades(gc,100);
96 ofstream out("good_cascades");
97 if (out) {
98 for (Int_t ngd=0; ngd<ngood; ngd++)
99 out<<gc[ngd].nlab<<' '<<gc[ngd].plab<<' '<<
100 gc[ngd].blab<<' '<<gc[ngd].code<<' '<<
101 gc[ngd].px<<' '<<gc[ngd].py<<' '<<gc[ngd].pz<<' '<<
102 gc[ngd].x <<' '<<gc[ngd].y <<' '<<gc[ngd].z <<endl;
103 } else cerr<<"Can not open file (good_cascades) !\n";
104 out.close();
105 }
106
107 /*** create some histograms ***/
108 TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution
109 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
110 TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
111 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
112 TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
113 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
114
115 TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res.
116 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
117 TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.); //y res
118 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
119 TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.); //z res.
120 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
121
122 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
123 TH1F *hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);
124 TH1F *hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
125 TH1F *hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
126 TH1F *hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
127 hg->SetLineColor(4); hg->SetLineWidth(2);
128 TH1F *hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
129 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
130
131 Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
132 TH1F *cs =new TH1F("v0s","Cascade Effective Mass",40, mmin, mmax);
133 cs->SetXTitle("(GeV)"); cs->SetFillColor(6);
134
135 Double_t pxg=0.,pyg=0.,ptg=0.;
136 Int_t nlab=-1, plab=-1, blab=-1;
137 Int_t i;
138 for (i=0; i<nentr; i++) {
139 AliCascadeVertex *cascade=(AliCascadeVertex*)carray.UncheckedAt(i);
140 cascade->GetV0labels(nlab,plab);
141 nlab=TMath::Abs(nlab); plab=TMath::Abs(plab);
142 blab=TMath::Abs(cascade->GetBachelorLabel());
143
144 cascade->ChangeMassHypothesis(code);
145
146 Double_t mass=cascade->GetEffMass();
147 cs->Fill(mass);
148
149 if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
150
151 Int_t j;
152 for (j=0; j<ngood; j++) {
153 if (gc[j].code != cascade->GetPdgCode()) continue;
154 if (gc[j].nlab == nlab)
155 if (gc[j].plab == plab)
156 if (gc[j].blab == blab) break;
157 }
158
159 Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
160 Double_t pt=TMath::Sqrt(px*px+py*py);
161
162 if (j==ngood) {
163 hfake->Fill(pt);
164 cerr<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
165 continue;
166 }
167
168 pxg=gc[j].px; pyg=gc[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
169 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
170 Double_t lamg=TMath::ATan2(gc[j].pz,ptg), lam=TMath::ATan2(pz,pt);
171 hp->Fill((phi - phig)*1000.);
172 hl->Fill((lam - lamg)*1000.);
173 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
174
175 Double_t x,y,z; cascade->GetXYZ(x,y,z);
176 hx->Fill((x-gc[j].x)*10);
177 hy->Fill((y-gc[j].y)*10);
178 hz->Fill((z-gc[j].z)*10);
179
180 hfound->Fill(ptg);
181 gc[j].nlab=-1;
182
183 }
184 for (i=0; i<ngood; i++) {
185 if (gc[i].code != code) continue;
186 pxg=gc[i].px; pyg=gc[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
187 hgood->Fill(ptg);
188 nlab=gc[i].nlab; plab=gc[i].plab; blab=gc[i].blab;
189 if (nlab < 0) continue;
190 cerr<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
191 }
192
193 carray.Delete();
194
195 Stat_t ng=hgood->GetEntries();
196 Stat_t nf=hfound->GetEntries();
197
198 cerr<<"Number of found cascades: "<<nentr<<" ("<<nf<<" in the peak)\n";
199 cerr<<"Number of good cascades: "<<ng<<endl;
200
201 if (ng!=0)
202 cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
203
204 gStyle->SetOptStat(111110);
205 gStyle->SetOptFit(1);
206
207 TCanvas *c1=new TCanvas("c1","",0,0,580,610);
208 c1->Divide(2,2);
209
210 c1->cd(1);
211 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
212 //hp->Fit("gaus");
213 hp->Draw();
214 hl->Draw("same"); c1->cd();
215
216 c1->cd(2);
217 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
218 //hpt->Fit("gaus"); c1->cd();
219 hpt->Draw(); c1->cd();
220
221 c1->cd(3);
222 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
223 //hx->Fit("gaus");
224 hx->Draw();
225 hy->Draw("same"); c1->cd();
226
227 c1->cd(4);
228 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
229 //hz->Fit("gaus");
230 hz->Draw();
231
232
233 TCanvas *c2=new TCanvas("c2","",600,0,580,610);
234 c2->Divide(1,2);
235
236 c2->cd(1);
237 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
238 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
239 hg->Divide(hfound,hgood,1,1.,"b");
240 hf->Divide(hfake,hgood,1,1.,"b");
241 hg->SetMaximum(1.4);
242 hg->SetYTitle("Cascade reconstruction efficiency");
243 hg->SetXTitle("Pt (GeV/c)");
244 hg->Draw();
245
246 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
247 line1->Draw("same");
248 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
249 line2->Draw("same");
250
251 hf->SetFillColor(1);
252 hf->SetFillStyle(3013);
253 hf->SetLineColor(2);
254 hf->SetLineWidth(2);
255 hf->Draw("histsame");
256 TText *text = new TText(0.461176,0.248448,"Fake cascades");
257 text->SetTextSize(0.05);
258 text->Draw();
259 text = new TText(0.453919,1.11408,"Good cascades");
260 text->SetTextSize(0.05);
261 text->Draw();
262
263
264 c2->cd(2);
265 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
266 cs->SetXTitle("(GeV/c)");
267 cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
268 Double_t max=cs->GetMaximum();
269 TLine *line3 =
270 new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
271 line3->Draw("same");
272 TLine *line4 =
273 new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
274 line4->Draw("same");
275
276 timer.Stop(); timer.Print();
277
278 return 0;
279}
280
281
282
283Int_t good_cascades(GoodCascade *gc, Int_t max) {
284 Int_t nc=0;
285 /*** Get information about the cuts ***/
286 Double_t r2min=0.9*0.9;
287 Double_t r2max=2.9*2.9;
288
289 /*** Get labels of the "good" tracks ***/
290 Double_t dd; Int_t id, label[15000], ngt=0;
291 ifstream in("good_tracks_its");
292 if (!in) {
293 cerr<<"Can't open the file good_tracks_its \n";
294 return nc;
295 }
296 while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
297 ngt++;
298 if (ngt>=15000) {
299 cerr<<"Too many good ITS tracks !\n";
300 return nc;
301 }
302 }
303 if (!in.eof()) {
304 cerr<<"Read error (good_tracks_its) !\n";
305 return nc;
306 }
307
308 /*** Get an access to the kinematics ***/
309 if (gAlice) {delete gAlice; gAlice=0;}
310
311 TFile *file=TFile::Open("galice.root");
312 if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
313 if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
314 cerr<<"gAlice has not been found on galice.root !\n";
315 exit(5);
316 }
317
318 Int_t np=gAlice->GetEvent(0);
319 while (np--) {
320 cerr<<np<<'\r';
321 TParticle *cp=gAlice->Particle(np);
322
323 /*** only these cascades are "good" ***/
324 Int_t code=cp->GetPdgCode();
325 if (code!=kXiMinus) if (code!=kXiPlusBar)
326 if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue;
327
328 /*** daughter tracks must be "good" ***/
329 Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
330 if (v0lab==blab) continue;
331 if (v0lab<0) continue;
332 if (blab<0) continue;
333
334 TParticle *p0=gAlice->Particle(v0lab);
335 TParticle *bp=gAlice->Particle(blab);
336 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
337 TParticle *p=p0; p0=bp; bp=p;
338 Int_t i=v0lab; v0lab=blab; blab=i;
339 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar))
340 continue;
341 }
342
343 /** is the bachelor "good" ? **/
344 Int_t i;
345 for (i=0; i<ngt; i++) if (label[i]==blab) break;
346 if (i==ngt) continue;
347
348 /** is the V0 "good" ? **/
349 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
350 if (nlab==plab) continue;
351 if (nlab<0) continue;
352 if (plab<0) continue;
353
354 for (i=0; i<ngt; i++) if (label[i]==nlab) break;
355 if (i==ngt) continue;
356 for (i=0; i<ngt; i++) if (label[i]==plab) break;
357 if (i==ngt) continue;
358
359 /*** fiducial volume ***/
360 Double_t x=bp->Vx(), y=bp->Vy(), z=bp->Vz(), r2=x*x+y*y; //bachelor
361 if (r2<r2min) continue;
362 if (r2>r2max) continue;
363 TParticle *pp=gAlice->Particle(plab);
364 x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y; //V0
365 if (r2<r2min) continue;
366 if (r2>r2max) continue;
367
368 if (gAlice->Particle(plab)->GetPDG()->Charge() < 0.) {
369 i=plab; plab=nlab; nlab=i;
370 }
371
372 gc[nc].code=code;
373 gc[nc].plab=plab; gc[nc].nlab=nlab; gc[nc].blab=blab;
374 gc[nc].px=cp->Px(); gc[nc].py=cp->Py(); gc[nc].pz=cp->Pz();
375 gc[nc].x=x; gc[nc].y=y; gc[nc].z=z;
376 nc++;
377
378 }
379
380
381 return nc;
382}