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__)
13 #include <Riostream.h>
16 #include <TParticle.h>
20 #include <TBenchmark.h>
26 #include "AliHeader.h"
27 #include "AliTrackReference.h"
28 #include "AliRunLoader.h"
32 const Int_t kXiMinus = 3312;
33 const Int_t kXiPlusBar = -3312;
34 const Int_t kOmegaMinus = 3334;
35 const Int_t kOmegaPlusBar = -3334;
38 Int_t GoodCascades(const Char_t *dir=".");
40 extern AliRun *gAlice;
41 extern TBenchmark *gBenchmark;
44 static Int_t allgood=0;
45 static Int_t allfound=0;
47 Int_t AliCascadeComparison(Int_t code=3312, const Char_t *dir=".") {
48 //code= 3312; //kXiMinus
49 //code=-3312; //kXiPlusBar
50 //code= 3334; //kOmegaMinus
51 //code=-3334; //kOmegaPlusBar
52 gBenchmark->Start("AliCascadeComparison");
54 ::Info("AliCascadeComparison.C","Doing comparison...");
56 const Double_t cascadeWindow=0.05, cascadeWidth=0.015;
57 Double_t ptncut=0.12, ptpcut=0.33, kine0cut=0.003;
58 Double_t ptbcut=0.11, kinecut=0.002;
59 Double_t cascadeMass=1.32131;
64 ptncut=0.33; ptpcut=0.12;
69 ptbcut=0.22; kinecut=0.006;
74 ptncut=0.33; ptpcut=0.12;
75 ptbcut=0.22; kinecut=0.006;
77 default: cerr<<"Invalid PDG code !\n"; return 1;
81 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
82 if (!hp) hp=new TH1F("hp","Angular Resolution",30,-30.,30.);
83 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
85 TH1F *hl=(TH1F*)gROOT->FindObject("hl");
86 if (!hl) hl=new TH1F("hl","Lambda Resolution",30,-30,30);
87 hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
89 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
90 if (!hpt) hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
91 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
93 TH1F *hx=(TH1F*)gROOT->FindObject("hx");
94 if (!hx) hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.);
95 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
97 TH1F *hy=(TH1F*)gROOT->FindObject("hy");
98 if (!hy) hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.);
99 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
101 TH1F *hz=(TH1F*)gROOT->FindObject("hz");
102 if (!hz) hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.);
103 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
106 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
107 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
108 if (!hgood) hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);
110 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
111 if (!hfound) hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
113 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
114 if (!hfake) hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
116 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
117 if (!hg) hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
118 hg->SetLineColor(4); hg->SetLineWidth(2);
120 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
121 if (!hf) hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
122 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
125 Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
126 TH1F *cs=(TH1F*)gROOT->FindObject("cs");
127 if (!cs) cs=new TH1F("cs","Cascade Effective Mass",40, mmin, mmax);
128 cs->SetXTitle("(GeV)"); cs->SetLineColor(4); cs->SetLineWidth(4);
130 TH1F *csf=(TH1F*)gROOT->FindObject("csf");
131 if (!csf) csf=new TH1F("csf","Fake Cascade Effective Mass",40, mmin, mmax);
132 csf->SetXTitle("(GeV)"); csf->SetFillColor(6);
136 sprintf(fname,"%s/GoodCascades.root",dir);
138 TFile *refFile=TFile::Open(fname,"old");
139 if (!refFile || !refFile->IsOpen()) {
140 ::Info("AliCascadeComparison.C","Marking good cascades (will take a while)...");
141 if (GoodCascades(dir)) {
142 ::Error("AliCascadesComparison.C","Can't generate the reference file !");
146 refFile=TFile::Open(fname,"old");
147 if (!refFile || !refFile->IsOpen()) {
148 ::Error("AliCascadeComparison.C","Can't open the reference file !");
152 TTree *csTree=(TTree*)refFile->Get("csTree");
154 ::Error("AliCascadeComparison.C","Can't get the reference tree !");
157 TBranch *pbranch=csTree->GetBranch("positive");
159 ::Error("AliCascadeComparison.C","Can't get the positive daughter branch !");
162 TClonesArray dummy("AliTrackReference",1000), *prefs=&dummy;
163 pbranch->SetAddress(&prefs);
165 TBranch *nbranch=csTree->GetBranch("negative");
167 ::Error("AliCascadeComparison.C","Can't get the negative daughter branch !");
170 TClonesArray dumm("AliTrackReference",1000), *nrefs=&dumm;
171 nbranch->SetAddress(&nrefs);
173 TBranch *bbranch=csTree->GetBranch("bachelor");
175 ::Error("AliCascadeComparison.C","Can't get the bachelor branch !");
178 TClonesArray dum("AliTrackReference",1000), *brefs=&dum;
179 bbranch->SetAddress(&brefs);
183 sprintf(fname,"%s/AliESDs.root",dir);
184 TFile *ef=TFile::Open(fname);
185 if ((!ef)||(!ef->IsOpen())) {
186 sprintf(fname,"%s/AliESDcascade.root",dir);
187 ef=TFile::Open(fname);
188 if ((!ef)||(!ef->IsOpen())) {
189 ::Error("AliCascadeComparison.C","Can't open AliESDcascade.root !");
193 AliESD* event = new AliESD;
194 TTree* esdTree = (TTree*) ef->Get("esdTree");
196 ::Error("AliCascadeComparison.C", "no ESD tree found");
199 esdTree->SetBranchAddress("ESD", &event);
202 //******* Loop over events *********
204 while (esdTree->GetEvent(e)) {
205 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
207 Int_t nentr=event->GetNumberOfCascades();
211 if (csTree->GetEvent(e++)==0) {
212 cerr<<"No reconstructable cascades !\n";
216 Int_t ngood=prefs->GetEntriesFast(),ng=0;
218 Double_t pxg=0.,pyg=0.,pzg=0.,ptg=0.;
219 Int_t nlab=-1, plab=-1, blab=-1;
221 for (i=0; i<nentr; i++) {
222 AliESDcascade *cascade=event->GetCascade(i);
224 Int_t nidx=TMath::Abs(cascade->GetNindex());
225 Int_t pidx=TMath::Abs(cascade->GetPindex());
226 Int_t bidx=TMath::Abs(cascade->GetBindex());
228 AliESDtrack *ntrack=event->GetTrack(nidx);
229 AliESDtrack *ptrack=event->GetTrack(pidx);
230 AliESDtrack *btrack=event->GetTrack(bidx);
232 nlab=TMath::Abs(ntrack->GetLabel());
233 plab=TMath::Abs(ptrack->GetLabel());
234 blab=TMath::Abs(btrack->GetLabel());
236 /** Kinematical cuts **/
237 Double_t pxn,pyn,pzn; cascade->GetNPxPyPz(pxn,pyn,pzn);
238 Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
239 if (ptn < ptncut) continue;
240 Double_t pxp,pyp,pzp; cascade->GetPPxPyPz(pxp,pyp,pzp);
241 Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
242 if (ptp < ptpcut) continue;
243 Double_t pxb,pyb,pzb; cascade->GetBPxPyPz(pxb,pyb,pzb);
244 Double_t ptb=TMath::Sqrt(pxb*pxb + pyb*pyb);
245 if (ptb < ptbcut) continue;
247 Double_t kine=cascade->ChangeMassHypothesis(kine0,code);
248 if (TMath::Abs(kine0)>kine0cut) continue;
249 //if (TMath::Abs(kine)>kinecut) continue;
251 Double_t mass=cascade->GetEffMass();
255 AliTrackReference *nref=0, *pref=0, *bref=0;
257 for (j=0; j<ngood; j++) {
258 bref=(AliTrackReference*)brefs->UncheckedAt(j);
259 nref=(AliTrackReference*)nrefs->UncheckedAt(j);
260 pref=(AliTrackReference*)prefs->UncheckedAt(j);
261 if (bref->Label() == blab)
262 if (nref->Label() == nlab)
263 if (pref->Label() == plab) break;
266 if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
268 Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
269 Double_t pt=TMath::Sqrt(px*px+py*py);
273 cout<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
278 pxg=bref->Px()+nref->Px()+pref->Px();
279 pyg=bref->Px()+nref->Py()+pref->Py();
280 pzg=nref->Pz()+pref->Pz();
281 ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
282 Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
283 Double_t lamg=TMath::ATan2(pzg,ptg), lam=TMath::ATan2(pz,pt);
284 hp->Fill((phi - phig)*1000.);
285 hl->Fill((lam - lamg)*1000.);
286 hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
288 Double_t x,y,z; cascade->GetXYZ(x,y,z);
289 hx->Fill((x-nref->X())*10);
290 hy->Fill((y-nref->Y())*10);
291 hz->Fill((z-nref->Z())*10);
297 for (i=0; i<ngood; i++) {
298 AliTrackReference *bref=(AliTrackReference*)brefs->UncheckedAt(i);
299 AliTrackReference *nref=(AliTrackReference*)nrefs->UncheckedAt(i);
300 AliTrackReference *pref=(AliTrackReference*)prefs->UncheckedAt(i);
301 Int_t pdg=(Int_t)nref->GetLength(); //this is the cascade's PDG !
302 if (code!=pdg) continue;
304 pxg=bref->Px()+nref->Px()+pref->Px();
305 pyg=bref->Px()+nref->Py()+pref->Py();
306 ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
308 nlab=nref->Label(); plab=pref->Label(); blab=bref->Label();
309 if (nlab < 0) continue;
310 cout<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
314 cout<<"Number of found cascades : "<<nentr<<endl;
315 cout<<"Number of \"good\" cascades : "<<ng<<endl;
321 } //**** End of the loop over events
329 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
330 if (ng!=0) cout<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
332 "Total number of found cascades: "<<allfound<<" ("<<nf<<" in the peak)\n";
333 cout<<"Total number of \"good\" cascades: "<<allgood<<endl;
336 gStyle->SetOptStat(111110);
337 gStyle->SetOptFit(1);
339 TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
341 c1=new TCanvas("c1","",0,0,580,610);
348 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
349 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus");
350 hl->Draw("same"); c1->cd();
353 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
354 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus");
358 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
359 if (hx->GetEntries()<minc) hx->Draw(); else hx->Fit("gaus");
360 hy->Draw("same"); c1->cd();
363 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
364 if (hz->GetEntries()<minc) hz->Draw(); else hz->Fit("gaus");
368 TCanvas *c2=(TCanvas*)gROOT->FindObject("c2");
370 c2=new TCanvas("c2","",600,0,580,610);
375 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
376 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
377 hg->Divide(hfound,hgood,1,1.,"b");
378 hf->Divide(hfake,hgood,1,1.,"b");
380 hg->SetYTitle("Cascade reconstruction efficiency");
381 hg->SetXTitle("Pt (GeV/c)");
384 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
386 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
390 hf->SetFillStyle(3013);
393 hf->Draw("histsame");
394 TText *text = new TText(0.461176,0.248448,"Fake cascades");
395 text->SetTextSize(0.05);
397 text = new TText(0.453919,1.11408,"Good cascades");
398 text->SetTextSize(0.05);
403 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
404 if (cs->GetEntries()<minc) cs->Draw();
405 else cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
407 Double_t max=cs->GetMaximum();
409 new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
412 new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
417 TFile fc("AliCascadeComparison.root","RECREATE");
422 gBenchmark->Stop("AliCascadeComparison");
423 gBenchmark->Show("AliCascadeComparison");
430 Int_t GoodCascades(const Char_t *dir) {
432 delete gAlice->GetRunLoader();
438 sprintf(fname,"%s/galice.root",dir);
440 AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
442 ::Error("GoodCascades","Can't start session !");
448 rl->LoadKinematics();
451 Int_t nev=rl->GetNumberOfEvents();
452 ::Info("GoodCascades","Number of events : %d\n",nev);
455 sprintf(fname,"%s/GoodTracksITS.root",dir);
456 TFile *itsFile=TFile::Open(fname);
457 if ((!itsFile)||(!itsFile->IsOpen())) {
458 ::Error("GoodCAscades","Can't open the GoodTracksITS.root !");
462 TClonesArray dm("AliTrackReference",1000), *itsRefs=&dm;
463 TTree *itsTree=(TTree*)itsFile->Get("itsTree");
465 ::Error("GoodCascades","Can't get the ITS reference tree !");
469 TBranch *itsBranch=itsTree->GetBranch("ITS");
471 ::Error("GoodCascades","Can't get the ITS reference branch !");
475 itsBranch->SetAddress(&itsRefs);
478 sprintf(fname,"%s/GoodCascades.root",dir);
479 TFile *csFile=TFile::Open(fname,"recreate");
480 TClonesArray dummy("AliTrackReference",1000), *nrefs=&dummy;
481 TClonesArray dumm("AliTrackReference",1000), *prefs=&dumm;
482 TClonesArray dum("AliTrackReference",1000), *brefs=&dum;
483 TTree csTree("csTree","Tree with info about the reconstructable cascades");
484 csTree.Branch("negative",&nrefs);
485 csTree.Branch("positive",&prefs);
486 csTree.Branch("bachelor",&brefs);
489 // *** Get information about the cuts ***
490 Double_t r2min=0.9*0.9;
491 Double_t r2max=2.9*2.9;
494 //******** Loop over generated events
495 for (Int_t e=0; e<nev; e++) {
496 rl->GetEvent(e); csFile->cd();
498 Int_t np = rl->GetHeader()->GetNtrack();
499 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
501 itsTree->GetEvent(e);
502 Int_t nk=itsRefs->GetEntriesFast();
504 AliStack *stack=rl->Stack();
506 AliTrackReference *nref=0, *pref=0, *bref=0;
511 TParticle *cp=stack->Particle(np);
513 // *** only these cascades are "good" ***
514 Int_t code=cp->GetPdgCode();
515 if (code!=kXiMinus) if (code!=kXiPlusBar)
516 if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue;
518 // *** daughter tracks must be "good" ***
519 Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
520 if (v0lab==blab) continue;
521 if (v0lab<0) continue;
522 if (blab<0) continue;
524 TParticle *p0=stack->Particle(v0lab);
525 TParticle *bp=stack->Particle(blab);
527 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
528 TParticle *p=p0; p0=bp; bp=p;
529 i=v0lab; v0lab=blab; blab=i;
530 if ((p0->GetPdgCode()!=kLambda0)&&(p0->GetPdgCode()!=kLambda0Bar))
534 // ** is the bachelor "good" ? **
535 for (i=0; i<nk; i++) {
536 bref=(AliTrackReference*)itsRefs->UncheckedAt(i);
537 if (bref->Label()==blab) break;
541 // ** is the V0 "good" ? **
542 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
543 if (nlab==plab) continue;
544 if (nlab<0) continue;
545 if (plab<0) continue;
546 if (stack->Particle(plab)->GetPDG()->Charge() < 0.) {
547 i=plab; plab=nlab; nlab=i;
549 for (i=0; i<nk; i++) {
550 nref=(AliTrackReference*)itsRefs->UncheckedAt(i);
551 if (nref->Label()==nlab) break;
554 for (i=0; i<nk; i++) {
555 pref=(AliTrackReference*)itsRefs->UncheckedAt(i);
556 if (pref->Label()==plab) break;
561 // *** fiducial volume ***
562 Double_t x=bp->Vx(), y=bp->Vy(), r2=x*x+y*y; //bachelor
563 if (r2<r2min) continue;
564 if (r2>r2max) continue;
565 TParticle *pp=stack->Particle(plab);
566 x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y; //V0
567 if (r2<r2min) continue;
568 if (r2>r2max) continue;
570 Int_t pdg=cp->GetPdgCode();
571 nref->SetLength(pdg); //This will the cascade's PDG !
573 new((*nrefs)[nc]) AliTrackReference(*nref);
574 new((*prefs)[nc]) AliTrackReference(*pref);
575 new((*brefs)[nc]) AliTrackReference(*bref);
582 nrefs->Clear(); prefs->Clear(); brefs->Clear();
584 } //**** end of the loop over generated events