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 ****************************************************************************/
10 #if !defined(__CINT__) || defined(__MAKECINT__)
11 #include <Riostream.h>
15 #include "AliHeader.h"
16 #include "AliRunLoader.h"
17 #include "AliITSLoader.h"
22 #include "TObjArray.h"
27 #include "TParticle.h"
28 #include "TStopwatch.h"
33 #include "AliCascadeVertex.h"
37 Int_t nlab,plab; // V0's daughter labels
38 Int_t blab; // Bachelor label
43 Int_t good_cascades(GoodCascade *gt, Int_t max);
45 extern AliRun *gAlice;
47 Int_t AliCascadeComparison(Int_t code=3312) {
48 //code= 3312; //kXiMinus
49 //code=-3312; //kXiPlusBar
50 //code= 3334; //kOmegaMinus
51 //code=-3334; //kOmegaPlusBar
53 cerr<<"Doing comparison...\n";
57 const Double_t cascadeWindow=0.05, cascadeWidth=0.015;
58 Double_t ptncut=0.12, ptpcut=0.33, kine0cut=0.003;
59 Double_t ptbcut=0.11, kinecut=0.002;
60 Double_t cascadeMass=1.32131;
65 ptncut=0.33; ptpcut=0.12;
70 ptbcut=0.22; kinecut=0.006;
75 ptncut=0.33; ptpcut=0.12;
76 ptbcut=0.22; kinecut=0.006;
78 default: cerr<<"Invalid PDG code !\n"; return 1;
81 /*** create some histograms ***/
82 TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution
83 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
84 TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
85 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
86 TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
87 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
89 TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res.
90 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
91 TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.); //y res
92 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
93 TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.); //z res.
94 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
96 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
97 TH1F *hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);
98 TH1F *hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
99 TH1F *hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
100 TH1F *hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
101 hg->SetLineColor(4); hg->SetLineWidth(2);
102 TH1F *hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
103 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
105 Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
106 TH1F *cs =new TH1F("cs","Cascade Effective Mass",40, mmin, mmax);
107 cs->SetXTitle("(GeV)");
108 cs->SetLineColor(4); cs->SetLineWidth(4);
109 TH1F *csf =new TH1F("csf","Fake Cascade Effective Mass",40, mmin, mmax);
110 csf->SetXTitle("(GeV)"); csf->SetFillColor(6);
113 delete gAlice->GetRunLoader();
117 AliRunLoader *rl = AliRunLoader::Open("galice.root");
119 cerr<<"AliV0Comparison.C :Can't start sesion !\n";
122 AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
124 cerr<<"AliV0Comparison.C : Can not find the ITSLoader\n";
129 /*** Load reconstructed cascades ***/
130 TObjArray carray(1000);
131 itsl->LoadCascades();
132 TTree *xTree=itsl->TreeX();
133 TBranch *branch=xTree->GetBranch("cascades");
134 Int_t nentr=(Int_t)xTree->GetEntries();
135 for (Int_t i=0; i<nentr; i++) {
136 AliCascadeVertex *iovertex=new AliCascadeVertex;
137 branch->SetAddress(&iovertex);
139 carray.AddLast(iovertex);
142 /*** Check if the file with the "good" cascades exists ***/
145 ifstream in("good_cascades");
147 cerr<<"Reading good cascades...\n";
148 while (in>>gc[ngood].nlab>>gc[ngood].plab>>
149 gc[ngood].blab>>gc[ngood].code>>
150 gc[ngood].px>>gc[ngood].py>>gc[ngood].pz>>
151 gc[ngood].x >>gc[ngood].y >>gc[ngood].z) {
155 cerr<<"Too many good cascades !\n";
159 if (!in.eof()) cerr<<"Read error (good_cascades) !\n";
161 /*** generate a file with the "good" cascades ***/
162 cerr<<"Marking good cascades (this will take a while)...\n";
163 ngood=good_cascades(gc,100);
164 ofstream out("good_cascades");
166 for (Int_t ngd=0; ngd<ngood; ngd++)
167 out<<gc[ngd].nlab<<' '<<gc[ngd].plab<<' '<<
168 gc[ngd].blab<<' '<<gc[ngd].code<<' '<<
169 gc[ngd].px<<' '<<gc[ngd].py<<' '<<gc[ngd].pz<<' '<<
170 gc[ngd].x <<' '<<gc[ngd].y <<' '<<gc[ngd].z <<endl;
171 } else cerr<<"Can not open file (good_cascades) !\n";
175 Double_t pxg=0.,pyg=0.,ptg=0.;
176 Int_t nlab=-1, plab=-1, blab=-1;
178 for (i=0; i<nentr; i++) {
179 AliCascadeVertex *cascade=(AliCascadeVertex*)carray.UncheckedAt(i);
180 nlab=TMath::Abs(cascade->GetNlabel());
181 plab=TMath::Abs(cascade->GetPlabel());
182 blab=TMath::Abs(cascade->GetBlabel());
184 /** Kinematical cuts **/
185 Double_t pxn,pyn,pzn; cascade->GetNPxPyPz(pxn,pyn,pzn);
186 Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
187 if (ptn < ptncut) continue;
188 Double_t pxp,pyp,pzp; cascade->GetPPxPyPz(pxp,pyp,pzp);
189 Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
190 if (ptp < ptpcut) continue;
191 Double_t pxb,pyb,pzb; cascade->GetBPxPyPz(pxb,pyb,pzb);
192 Double_t ptb=TMath::Sqrt(pxb*pxb + pyb*pyb);
193 if (ptb < ptbcut) continue;
195 Double_t kine=cascade->ChangeMassHypothesis(kine0,code);
196 if (TMath::Abs(kine0)>kine0cut) continue;
197 //if (TMath::Abs(kine)>kinecut) continue;
199 Double_t mass=cascade->GetEffMass();
203 if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
206 for (j=0; j<ngood; j++) {
207 if (gc[j].code != cascade->GetPdgCode()) continue;
208 if (gc[j].nlab == nlab)
209 if (gc[j].plab == plab)
210 if (gc[j].blab == blab) break;
213 Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
214 Double_t pt=TMath::Sqrt(px*px+py*py);
218 cerr<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
223 pxg=gc[j].px; pyg=gc[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
224 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
225 Double_t lamg=TMath::ATan2(gc[j].pz,ptg), lam=TMath::ATan2(pz,pt);
226 hp->Fill((phi - phig)*1000.);
227 hl->Fill((lam - lamg)*1000.);
228 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
230 Double_t x,y,z; cascade->GetXYZ(x,y,z);
231 hx->Fill((x-gc[j].x)*10);
232 hy->Fill((y-gc[j].y)*10);
233 hz->Fill((z-gc[j].z)*10);
239 for (i=0; i<ngood; i++) {
240 if (gc[i].code != code) continue;
241 pxg=gc[i].px; pyg=gc[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
243 nlab=gc[i].nlab; plab=gc[i].plab; blab=gc[i].blab;
244 if (nlab < 0) continue;
245 cerr<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
250 Stat_t ng=hgood->GetEntries();
251 Stat_t nf=hfound->GetEntries();
253 cerr<<"Number of found cascades: "<<nentr<<" ("<<nf<<" in the peak)\n";
254 cerr<<"Number of good cascades: "<<ng<<endl;
257 cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
259 gStyle->SetOptStat(111110);
260 gStyle->SetOptFit(1);
262 TCanvas *c1=new TCanvas("c1","",0,0,580,610);
266 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
269 hl->Draw("same"); c1->cd();
272 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
273 //hpt->Fit("gaus"); c1->cd();
274 hpt->Draw(); c1->cd();
277 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
280 hy->Draw("same"); c1->cd();
283 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
288 TCanvas *c2=new TCanvas("c2","",600,0,580,610);
292 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
293 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
294 hg->Divide(hfound,hgood,1,1.,"b");
295 hf->Divide(hfake,hgood,1,1.,"b");
297 hg->SetYTitle("Cascade reconstruction efficiency");
298 hg->SetXTitle("Pt (GeV/c)");
301 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
303 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
307 hf->SetFillStyle(3013);
310 hf->Draw("histsame");
311 TText *text = new TText(0.461176,0.248448,"Fake cascades");
312 text->SetTextSize(0.05);
314 text = new TText(0.453919,1.11408,"Good cascades");
315 text->SetTextSize(0.05);
320 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
321 //cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
324 Double_t max=cs->GetMaximum();
326 new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
329 new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
332 timer.Stop(); timer.Print();
340 Int_t good_cascades(GoodCascade *gc, Int_t max) {
342 /*** Get information about the cuts ***/
343 Double_t r2min=0.9*0.9;
344 Double_t r2max=2.9*2.9;
346 /*** Get labels of the "good" tracks ***/
347 Double_t dd; Int_t id, label[15000], ngt=0;
348 ifstream in("good_tracks_its");
350 cerr<<"Can't open the file good_tracks_its \n";
353 while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
356 cerr<<"Too many good ITS tracks !\n";
361 cerr<<"Read error (good_tracks_its) !\n";
365 /*** Get an access to the kinematics ***/
367 AliRunLoader::GetRunLoader(AliConfig::fgkDefaultEventFolderName);
369 ::Fatal("AliCascadeComparison.C::good_cascades","Can not find Run Loader !");
372 AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
374 cerr<<"AliITSComparisonV2.C : Can not find TPCLoader\n";
380 rl->LoadKinematics();
381 Int_t np = rl->GetHeader()->GetNtrack();
385 TParticle *cp=gAlice->Particle(np);
387 /*** only these cascades are "good" ***/
388 Int_t code=cp->GetPdgCode();
389 if (code!=kXiMinus) if (code!=kXiPlusBar)
390 if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue;
392 /*** daughter tracks must be "good" ***/
393 Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
394 if (v0lab==blab) continue;
395 if (v0lab<0) continue;
396 if (blab<0) continue;
398 TParticle *p0=gAlice->Particle(v0lab);
399 TParticle *bp=gAlice->Particle(blab);
400 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
401 TParticle *p=p0; p0=bp; bp=p;
402 Int_t i=v0lab; v0lab=blab; blab=i;
403 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar))
407 /** is the bachelor "good" ? **/
409 for (i=0; i<ngt; i++) if (label[i]==blab) break;
410 if (i==ngt) continue;
412 /** is the V0 "good" ? **/
413 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
414 if (nlab==plab) continue;
415 if (nlab<0) continue;
416 if (plab<0) continue;
418 for (i=0; i<ngt; i++) if (label[i]==nlab) break;
419 if (i==ngt) continue;
420 for (i=0; i<ngt; i++) if (label[i]==plab) break;
421 if (i==ngt) continue;
423 /*** fiducial volume ***/
424 Double_t x=bp->Vx(), y=bp->Vy(), r2=x*x+y*y; //bachelor
425 if (r2<r2min) continue;
426 if (r2>r2max) continue;
427 TParticle *pp=gAlice->Particle(plab);
428 x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y; //V0
429 if (r2<r2min) continue;
430 if (r2>r2max) continue;
432 if (gAlice->Particle(plab)->GetPDG()->Charge() < 0.) {
433 i=plab; plab=nlab; nlab=i;
437 gc[nc].plab=plab; gc[nc].nlab=nlab; gc[nc].blab=blab;
438 gc[nc].px=cp->Px(); gc[nc].py=cp->Py(); gc[nc].pz=cp->Pz();
439 gc[nc].x=bp->Vx(); gc[nc].y=bp->Vy(); gc[nc].z=bp->Vz();