1 /****************************************************************************
2 * Very important, delicate and rather obscure macro. *
4 * Creates list of "findable" cascades, *
5 * calculates efficiency, resolutions etc. *
7 * Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr *
8 ****************************************************************************/
17 #include "TObjArray.h"
22 #include "TParticle.h"
23 #include "TStopwatch.h"
27 #include "AliCascadeVertex.h"
31 Int_t nlab,plab; // V0's daughter labels
32 Int_t blab; // Bachelor label
37 Int_t good_cascades(GoodCascade *gt, Int_t max);
39 Int_t AliCascadeComparison(Int_t code=3312) {
40 //code= 3312; //kXiMinus
41 //code=-3312; //kXiPlusBar
42 //code= 3334; //kOmegaMinus
43 //code=-3334; //kOmegaPlusBar
45 cerr<<"Doing comparison...\n";
47 const Double_t cascadeWindow=0.05, cascadeWidth=0.015;
48 Double_t cascadeMass=0.5;
51 case kXiPlusBar: cascadeMass=1.32131; break;
53 case kOmegaPlusBar: cascadeMass=1.67245; break;
54 default: cerr<<"Invalid PDG code !\n"; return 1;
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);
70 carray.AddLast(iovertex);
75 /*** Check if the file with the "good" cascades exists ***/
78 ifstream in("good_cascades");
80 cerr<<"Reading good cascades...\n";
81 while (in>>gc[ngood].nlab>>gc[ngood].plab>>
82 gc[ngood].blab>>gc[ngood].code>>
83 gc[ngood].px>>gc[ngood].py>>gc[ngood].pz>>
84 gc[ngood].x >>gc[ngood].y >>gc[ngood].z) {
88 cerr<<"Too many good cascades !\n";
92 if (!in.eof()) cerr<<"Read error (good_cascades) !\n";
94 /*** generate a file with the "good" cascades ***/
95 cerr<<"Marking good cascades (this will take a while)...\n";
96 ngood=good_cascades(gc,100);
97 ofstream out("good_cascades");
99 for (Int_t ngd=0; ngd<ngood; ngd++)
100 out<<gc[ngd].nlab<<' '<<gc[ngd].plab<<' '<<
101 gc[ngd].blab<<' '<<gc[ngd].code<<' '<<
102 gc[ngd].px<<' '<<gc[ngd].py<<' '<<gc[ngd].pz<<' '<<
103 gc[ngd].x <<' '<<gc[ngd].y <<' '<<gc[ngd].z <<endl;
104 } else cerr<<"Can not open file (good_cascades) !\n";
108 /*** create some histograms ***/
109 TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution
110 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
111 TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
112 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
113 TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
114 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
116 TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res.
117 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
118 TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.); //y res
119 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
120 TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.); //z res.
121 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
123 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
124 TH1F *hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);
125 TH1F *hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
126 TH1F *hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
127 TH1F *hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
128 hg->SetLineColor(4); hg->SetLineWidth(2);
129 TH1F *hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
130 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
132 Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
133 TH1F *cs =new TH1F("cs","Cascade Effective Mass",40, mmin, mmax);
134 cs->SetXTitle("(GeV)"); cs->SetFillColor(6);
136 Double_t pxg=0.,pyg=0.,ptg=0.;
137 Int_t nlab=-1, plab=-1, blab=-1;
139 for (i=0; i<nentr; i++) {
140 AliCascadeVertex *cascade=(AliCascadeVertex*)carray.UncheckedAt(i);
141 cascade->GetV0labels(nlab,plab);
142 nlab=TMath::Abs(nlab); plab=TMath::Abs(plab);
143 blab=TMath::Abs(cascade->GetBachelorLabel());
145 cascade->ChangeMassHypothesis(code);
147 Double_t mass=cascade->GetEffMass();
150 if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
153 for (j=0; j<ngood; j++) {
154 if (gc[j].code != cascade->GetPdgCode()) continue;
155 if (gc[j].nlab == nlab)
156 if (gc[j].plab == plab)
157 if (gc[j].blab == blab) break;
160 Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
161 Double_t pt=TMath::Sqrt(px*px+py*py);
165 cerr<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
169 pxg=gc[j].px; pyg=gc[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
170 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
171 Double_t lamg=TMath::ATan2(gc[j].pz,ptg), lam=TMath::ATan2(pz,pt);
172 hp->Fill((phi - phig)*1000.);
173 hl->Fill((lam - lamg)*1000.);
174 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
176 Double_t x,y,z; cascade->GetXYZ(x,y,z);
177 hx->Fill((x-gc[j].x)*10);
178 hy->Fill((y-gc[j].y)*10);
179 hz->Fill((z-gc[j].z)*10);
185 for (i=0; i<ngood; i++) {
186 if (gc[i].code != code) continue;
187 pxg=gc[i].px; pyg=gc[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
189 nlab=gc[i].nlab; plab=gc[i].plab; blab=gc[i].blab;
190 if (nlab < 0) continue;
191 cerr<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
196 Stat_t ng=hgood->GetEntries();
197 Stat_t nf=hfound->GetEntries();
199 cerr<<"Number of found cascades: "<<nentr<<" ("<<nf<<" in the peak)\n";
200 cerr<<"Number of good cascades: "<<ng<<endl;
203 cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
205 gStyle->SetOptStat(111110);
206 gStyle->SetOptFit(1);
208 TCanvas *c1=new TCanvas("c1","",0,0,580,610);
212 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
215 hl->Draw("same"); c1->cd();
218 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
219 //hpt->Fit("gaus"); c1->cd();
220 hpt->Draw(); c1->cd();
223 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
226 hy->Draw("same"); c1->cd();
229 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
234 TCanvas *c2=new TCanvas("c2","",600,0,580,610);
238 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
239 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
240 hg->Divide(hfound,hgood,1,1.,"b");
241 hf->Divide(hfake,hgood,1,1.,"b");
243 hg->SetYTitle("Cascade reconstruction efficiency");
244 hg->SetXTitle("Pt (GeV/c)");
247 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
249 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
253 hf->SetFillStyle(3013);
256 hf->Draw("histsame");
257 TText *text = new TText(0.461176,0.248448,"Fake cascades");
258 text->SetTextSize(0.05);
260 text = new TText(0.453919,1.11408,"Good cascades");
261 text->SetTextSize(0.05);
266 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
267 cs->SetXTitle("(GeV/c)");
268 //cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
270 Double_t max=cs->GetMaximum();
272 new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
275 new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
278 timer.Stop(); timer.Print();
285 Int_t good_cascades(GoodCascade *gc, Int_t max) {
287 /*** Get information about the cuts ***/
288 Double_t r2min=0.9*0.9;
289 Double_t r2max=2.9*2.9;
291 /*** Get labels of the "good" tracks ***/
292 Double_t dd; Int_t id, label[15000], ngt=0;
293 ifstream in("good_tracks_its");
295 cerr<<"Can't open the file good_tracks_its \n";
298 while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
301 cerr<<"Too many good ITS tracks !\n";
306 cerr<<"Read error (good_tracks_its) !\n";
310 /*** Get an access to the kinematics ***/
311 if (gAlice) {delete gAlice; gAlice=0;}
313 TFile *file=TFile::Open("galice.root");
314 if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
315 if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
316 cerr<<"gAlice has not been found on galice.root !\n";
320 Int_t np=gAlice->GetEvent(0);
323 TParticle *cp=gAlice->Particle(np);
325 /*** only these cascades are "good" ***/
326 Int_t code=cp->GetPdgCode();
327 if (code!=kXiMinus) if (code!=kXiPlusBar)
328 if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue;
330 /*** daughter tracks must be "good" ***/
331 Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
332 if (v0lab==blab) continue;
333 if (v0lab<0) continue;
334 if (blab<0) continue;
336 TParticle *p0=gAlice->Particle(v0lab);
337 TParticle *bp=gAlice->Particle(blab);
338 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
339 TParticle *p=p0; p0=bp; bp=p;
340 Int_t i=v0lab; v0lab=blab; blab=i;
341 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar))
345 /** is the bachelor "good" ? **/
347 for (i=0; i<ngt; i++) if (label[i]==blab) break;
348 if (i==ngt) continue;
350 /** is the V0 "good" ? **/
351 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
352 if (nlab==plab) continue;
353 if (nlab<0) continue;
354 if (plab<0) continue;
356 for (i=0; i<ngt; i++) if (label[i]==nlab) break;
357 if (i==ngt) continue;
358 for (i=0; i<ngt; i++) if (label[i]==plab) break;
359 if (i==ngt) continue;
361 /*** fiducial volume ***/
362 Double_t x=bp->Vx(), y=bp->Vy(), r2=x*x+y*y; //bachelor
363 if (r2<r2min) continue;
364 if (r2>r2max) continue;
365 TParticle *pp=gAlice->Particle(plab);
366 x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y; //V0
367 if (r2<r2min) continue;
368 if (r2>r2max) continue;
370 if (gAlice->Particle(plab)->GetPDG()->Charge() < 0.) {
371 i=plab; plab=nlab; nlab=i;
375 gc[nc].plab=plab; gc[nc].nlab=nlab; gc[nc].blab=blab;
376 gc[nc].px=cp->Px(); gc[nc].py=cp->Py(); gc[nc].pz=cp->Pz();
377 gc[nc].x=bp->Vx(); gc[nc].y=bp->Vy(); gc[nc].z=bp->Vz();