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 ****************************************************************************/
11 #include <Riostream.h>
17 #include "TObjArray.h"
22 #include "TParticle.h"
23 #include "TStopwatch.h"
28 #include "AliCascadeVertex.h"
32 Int_t nlab,plab; // V0's daughter labels
33 Int_t blab; // Bachelor label
38 Int_t good_cascades(GoodCascade *gt, Int_t max);
40 Int_t AliCascadeComparison(Int_t code=3312) {
41 //code= 3312; //kXiMinus
42 //code=-3312; //kXiPlusBar
43 //code= 3334; //kOmegaMinus
44 //code=-3334; //kOmegaPlusBar
46 cerr<<"Doing comparison...\n";
48 const Double_t cascadeWindow=0.05, cascadeWidth=0.015;
49 Double_t ptncut=0.12, ptpcut=0.33, kine0cut=0.003;
50 Double_t ptbcut=0.11, kinecut=0.002;
51 Double_t cascadeMass=1.32131;
56 ptncut=0.33; ptpcut=0.12;
61 ptbcut=0.22; kinecut=0.006;
66 ptncut=0.33; ptpcut=0.12;
67 ptbcut=0.22; kinecut=0.006;
69 default: cerr<<"Invalid PDG code !\n"; return 1;
74 /*** Load reconstructed cascades ***/
75 TFile *cf=TFile::Open("AliCascadeVertices.root");
76 if (!cf->IsOpen()){cerr<<"Can't open AliCascadeVertices.root !\n";return 2;}
77 TObjArray carray(1000);
78 TTree *cTree=(TTree*)cf->Get("TreeCasc");
79 TBranch *branch=cTree->GetBranch("cascades");
80 Int_t nentr=(Int_t)cTree->GetEntries();
81 for (Int_t i=0; i<nentr; i++) {
82 AliCascadeVertex *iovertex=new AliCascadeVertex;
83 branch->SetAddress(&iovertex);
85 carray.AddLast(iovertex);
90 /*** Check if the file with the "good" cascades exists ***/
93 ifstream in("good_cascades");
95 cerr<<"Reading good cascades...\n";
96 while (in>>gc[ngood].nlab>>gc[ngood].plab>>
97 gc[ngood].blab>>gc[ngood].code>>
98 gc[ngood].px>>gc[ngood].py>>gc[ngood].pz>>
99 gc[ngood].x >>gc[ngood].y >>gc[ngood].z) {
103 cerr<<"Too many good cascades !\n";
107 if (!in.eof()) cerr<<"Read error (good_cascades) !\n";
109 /*** generate a file with the "good" cascades ***/
110 cerr<<"Marking good cascades (this will take a while)...\n";
111 ngood=good_cascades(gc,100);
112 ofstream out("good_cascades");
114 for (Int_t ngd=0; ngd<ngood; ngd++)
115 out<<gc[ngd].nlab<<' '<<gc[ngd].plab<<' '<<
116 gc[ngd].blab<<' '<<gc[ngd].code<<' '<<
117 gc[ngd].px<<' '<<gc[ngd].py<<' '<<gc[ngd].pz<<' '<<
118 gc[ngd].x <<' '<<gc[ngd].y <<' '<<gc[ngd].z <<endl;
119 } else cerr<<"Can not open file (good_cascades) !\n";
123 /*** create some histograms ***/
124 TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution
125 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
126 TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
127 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
128 TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
129 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
131 TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res.
132 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
133 TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.); //y res
134 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
135 TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.); //z res.
136 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
138 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
139 TH1F *hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);
140 TH1F *hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
141 TH1F *hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
142 TH1F *hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
143 hg->SetLineColor(4); hg->SetLineWidth(2);
144 TH1F *hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
145 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
147 Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
148 TH1F *cs =new TH1F("cs","Cascade Effective Mass",40, mmin, mmax);
149 cs->SetXTitle("(GeV)");
150 cs->SetLineColor(4); cs->SetLineWidth(4);
151 TH1F *csf =new TH1F("csf","Fake Cascade Effective Mass",40, mmin, mmax);
152 csf->SetXTitle("(GeV)"); csf->SetFillColor(6);
154 Double_t pxg=0.,pyg=0.,ptg=0.;
155 Int_t nlab=-1, plab=-1, blab=-1;
157 for (i=0; i<nentr; i++) {
158 AliCascadeVertex *cascade=(AliCascadeVertex*)carray.UncheckedAt(i);
159 nlab=TMath::Abs(cascade->GetNlabel());
160 plab=TMath::Abs(cascade->GetPlabel());
161 blab=TMath::Abs(cascade->GetBlabel());
163 /** Kinematical cuts **/
164 Double_t pxn,pyn,pzn; cascade->GetNPxPyPz(pxn,pyn,pzn);
165 Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
166 if (ptn < ptncut) continue;
167 Double_t pxp,pyp,pzp; cascade->GetPPxPyPz(pxp,pyp,pzp);
168 Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
169 if (ptp < ptpcut) continue;
170 Double_t pxb,pyb,pzb; cascade->GetBPxPyPz(pxb,pyb,pzb);
171 Double_t ptb=TMath::Sqrt(pxb*pxb + pyb*pyb);
172 if (ptb < ptbcut) continue;
174 Double_t kine=cascade->ChangeMassHypothesis(kine0,code);
175 if (TMath::Abs(kine0)>kine0cut) continue;
176 //if (TMath::Abs(kine)>kinecut) continue;
178 Double_t mass=cascade->GetEffMass();
182 if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
185 for (j=0; j<ngood; j++) {
186 if (gc[j].code != cascade->GetPdgCode()) continue;
187 if (gc[j].nlab == nlab)
188 if (gc[j].plab == plab)
189 if (gc[j].blab == blab) break;
192 Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
193 Double_t pt=TMath::Sqrt(px*px+py*py);
197 cerr<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
202 pxg=gc[j].px; pyg=gc[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
203 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
204 Double_t lamg=TMath::ATan2(gc[j].pz,ptg), lam=TMath::ATan2(pz,pt);
205 hp->Fill((phi - phig)*1000.);
206 hl->Fill((lam - lamg)*1000.);
207 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
209 Double_t x,y,z; cascade->GetXYZ(x,y,z);
210 hx->Fill((x-gc[j].x)*10);
211 hy->Fill((y-gc[j].y)*10);
212 hz->Fill((z-gc[j].z)*10);
218 for (i=0; i<ngood; i++) {
219 if (gc[i].code != code) continue;
220 pxg=gc[i].px; pyg=gc[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
222 nlab=gc[i].nlab; plab=gc[i].plab; blab=gc[i].blab;
223 if (nlab < 0) continue;
224 cerr<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
229 Stat_t ng=hgood->GetEntries();
230 Stat_t nf=hfound->GetEntries();
232 cerr<<"Number of found cascades: "<<nentr<<" ("<<nf<<" in the peak)\n";
233 cerr<<"Number of good cascades: "<<ng<<endl;
236 cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
238 gStyle->SetOptStat(111110);
239 gStyle->SetOptFit(1);
241 TCanvas *c1=new TCanvas("c1","",0,0,580,610);
245 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
248 hl->Draw("same"); c1->cd();
251 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
252 //hpt->Fit("gaus"); c1->cd();
253 hpt->Draw(); c1->cd();
256 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
259 hy->Draw("same"); c1->cd();
262 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
267 TCanvas *c2=new TCanvas("c2","",600,0,580,610);
271 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
272 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
273 hg->Divide(hfound,hgood,1,1.,"b");
274 hf->Divide(hfake,hgood,1,1.,"b");
276 hg->SetYTitle("Cascade reconstruction efficiency");
277 hg->SetXTitle("Pt (GeV/c)");
280 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
282 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
286 hf->SetFillStyle(3013);
289 hf->Draw("histsame");
290 TText *text = new TText(0.461176,0.248448,"Fake cascades");
291 text->SetTextSize(0.05);
293 text = new TText(0.453919,1.11408,"Good cascades");
294 text->SetTextSize(0.05);
299 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
300 //cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
303 Double_t max=cs->GetMaximum();
305 new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
308 new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
311 timer.Stop(); timer.Print();
318 Int_t good_cascades(GoodCascade *gc, Int_t max) {
320 /*** Get information about the cuts ***/
321 Double_t r2min=0.9*0.9;
322 Double_t r2max=2.9*2.9;
324 /*** Get labels of the "good" tracks ***/
325 Double_t dd; Int_t id, label[15000], ngt=0;
326 ifstream in("good_tracks_its");
328 cerr<<"Can't open the file good_tracks_its \n";
331 while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
334 cerr<<"Too many good ITS tracks !\n";
339 cerr<<"Read error (good_tracks_its) !\n";
343 /*** Get an access to the kinematics ***/
344 if (gAlice) {delete gAlice; gAlice=0;}
346 TFile *file=TFile::Open("galice.root");
347 if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
348 if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
349 cerr<<"gAlice has not been found on galice.root !\n";
353 Int_t np=gAlice->GetEvent(0);
356 TParticle *cp=gAlice->Particle(np);
358 /*** only these cascades are "good" ***/
359 Int_t code=cp->GetPdgCode();
360 if (code!=kXiMinus) if (code!=kXiPlusBar)
361 if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue;
363 /*** daughter tracks must be "good" ***/
364 Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
365 if (v0lab==blab) continue;
366 if (v0lab<0) continue;
367 if (blab<0) continue;
369 TParticle *p0=gAlice->Particle(v0lab);
370 TParticle *bp=gAlice->Particle(blab);
371 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
372 TParticle *p=p0; p0=bp; bp=p;
373 Int_t i=v0lab; v0lab=blab; blab=i;
374 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar))
378 /** is the bachelor "good" ? **/
380 for (i=0; i<ngt; i++) if (label[i]==blab) break;
381 if (i==ngt) continue;
383 /** is the V0 "good" ? **/
384 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
385 if (nlab==plab) continue;
386 if (nlab<0) continue;
387 if (plab<0) continue;
389 for (i=0; i<ngt; i++) if (label[i]==nlab) break;
390 if (i==ngt) continue;
391 for (i=0; i<ngt; i++) if (label[i]==plab) break;
392 if (i==ngt) continue;
394 /*** fiducial volume ***/
395 Double_t x=bp->Vx(), y=bp->Vy(), r2=x*x+y*y; //bachelor
396 if (r2<r2min) continue;
397 if (r2>r2max) continue;
398 TParticle *pp=gAlice->Particle(plab);
399 x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y; //V0
400 if (r2<r2min) continue;
401 if (r2>r2max) continue;
403 if (gAlice->Particle(plab)->GetPDG()->Charge() < 0.) {
404 i=plab; plab=nlab; nlab=i;
408 gc[nc].plab=plab; gc[nc].nlab=nlab; gc[nc].blab=blab;
409 gc[nc].px=cp->Px(); gc[nc].py=cp->Py(); gc[nc].pz=cp->Pz();
410 gc[nc].x=bp->Vx(); gc[nc].y=bp->Vy(); gc[nc].z=bp->Vz();