Little changes needed for running with new AliESDEvent
[u/mrichter/AliRoot.git] / STEER / AliCascadeComparison.C
CommitLineData
ca28c5f5 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
566abf75 10#if !defined(__CINT__) || defined(__MAKECINT__)
4fa1e629 11 #include <TMath.h>
12 #include <TError.h>
ddae8318 13 #include <Riostream.h>
4fa1e629 14 #include <TH1F.h>
15 #include <TTree.h>
16 #include <TParticle.h>
17 #include <TCanvas.h>
18 #include <TLine.h>
19 #include <TText.h>
20 #include <TBenchmark.h>
21 #include <TStyle.h>
8b462fd8 22 #include <TFile.h>
4fa1e629 23 #include <TROOT.h>
24
25 #include "AliStack.h"
566abf75 26 #include "AliHeader.h"
4fa1e629 27 #include "AliTrackReference.h"
566abf75 28 #include "AliRunLoader.h"
ca28c5f5 29 #include "AliRun.h"
4514157e 30 #include "AliESDEvent.h"
4fa1e629 31#else
32const Int_t kXiMinus = 3312;
33const Int_t kXiPlusBar = -3312;
34const Int_t kOmegaMinus = 3334;
35const Int_t kOmegaPlusBar = -3334;
ca28c5f5 36#endif
37
4fa1e629 38Int_t GoodCascades(const Char_t *dir=".");
ca28c5f5 39
566abf75 40extern AliRun *gAlice;
4fa1e629 41extern TBenchmark *gBenchmark;
42extern TROOT *gROOT;
43
44static Int_t allgood=0;
45static Int_t allfound=0;
566abf75 46
4fa1e629 47Int_t AliCascadeComparison(Int_t code=3312, const Char_t *dir=".") {
ca28c5f5 48 //code= 3312; //kXiMinus
49 //code=-3312; //kXiPlusBar
50 //code= 3334; //kOmegaMinus
51 //code=-3334; //kOmegaPlusBar
4fa1e629 52 gBenchmark->Start("AliCascadeComparison");
5102bab6 53
4fa1e629 54 ::Info("AliCascadeComparison.C","Doing comparison...");
5102bab6 55
ca28c5f5 56 const Double_t cascadeWindow=0.05, cascadeWidth=0.015;
04b2a5f1 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;
ca28c5f5 60 switch (code) {
61 case kXiMinus:
04b2a5f1 62 break;
63 case kXiPlusBar:
64 ptncut=0.33; ptpcut=0.12;
65 break;
ca28c5f5 66 case kOmegaMinus:
04b2a5f1 67 cascadeMass=1.67245;
68 kine0cut=0.001;
69 ptbcut=0.22; kinecut=0.006;
70 break;
71 case kOmegaPlusBar:
72 cascadeMass=1.67245;
73 kine0cut=0.001;
74 ptncut=0.33; ptpcut=0.12;
75 ptbcut=0.22; kinecut=0.006;
76 break;
ca28c5f5 77 default: cerr<<"Invalid PDG code !\n"; return 1;
78 }
79
4fa1e629 80
81 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
82 if (!hp) hp=new TH1F("hp","Angular Resolution",30,-30.,30.);
566abf75 83 hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
4fa1e629 84
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);
88
89 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
90 if (!hpt) hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.);
566abf75 91 hpt->SetXTitle("(%)"); hpt->SetFillColor(2);
92
4fa1e629 93 TH1F *hx=(TH1F*)gROOT->FindObject("hx");
94 if (!hx) hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.);
566abf75 95 hx->SetXTitle("(mm)"); hx->SetFillColor(6);
4fa1e629 96
97 TH1F *hy=(TH1F*)gROOT->FindObject("hy");
98 if (!hy) hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.);
566abf75 99 hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
4fa1e629 100
101 TH1F *hz=(TH1F*)gROOT->FindObject("hz");
102 if (!hz) hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.);
566abf75 103 hz->SetXTitle("(mm)"); hz->SetFillColor(6);
104
4fa1e629 105
566abf75 106 Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
4fa1e629 107 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
108 if (!hgood) hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);
109
110 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
111 if (!hfound) hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
112
113 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
114 if (!hfake) hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
115
116 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
117 if (!hg) hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
566abf75 118 hg->SetLineColor(4); hg->SetLineWidth(2);
4fa1e629 119
120 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
121 if (!hf) hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
566abf75 122 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
123
4fa1e629 124
566abf75 125 Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
4fa1e629 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);
129
130 TH1F *csf=(TH1F*)gROOT->FindObject("csf");
131 if (!csf) csf=new TH1F("csf","Fake Cascade Effective Mass",40, mmin, mmax);
566abf75 132 csf->SetXTitle("(GeV)"); csf->SetFillColor(6);
133
ca28c5f5 134
4fa1e629 135 Char_t fname[100];
136 sprintf(fname,"%s/GoodCascades.root",dir);
137
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 !");
143 return 1;
144 }
ca28c5f5 145 }
4fa1e629 146 refFile=TFile::Open(fname,"old");
147 if (!refFile || !refFile->IsOpen()) {
148 ::Error("AliCascadeComparison.C","Can't open the reference file !");
149 return 1;
150 }
151
152 TTree *csTree=(TTree*)refFile->Get("csTree");
153 if (!csTree) {
154 ::Error("AliCascadeComparison.C","Can't get the reference tree !");
155 return 2;
ca28c5f5 156 }
4fa1e629 157 TBranch *pbranch=csTree->GetBranch("positive");
158 if (!pbranch) {
159 ::Error("AliCascadeComparison.C","Can't get the positive daughter branch !");
160 return 3;
161 }
162 TClonesArray dummy("AliTrackReference",1000), *prefs=&dummy;
163 pbranch->SetAddress(&prefs);
164
165 TBranch *nbranch=csTree->GetBranch("negative");
166 if (!nbranch) {
167 ::Error("AliCascadeComparison.C","Can't get the negative daughter branch !");
168 return 4;
169 }
170 TClonesArray dumm("AliTrackReference",1000), *nrefs=&dumm;
171 nbranch->SetAddress(&nrefs);
172
173 TBranch *bbranch=csTree->GetBranch("bachelor");
174 if (!nbranch) {
175 ::Error("AliCascadeComparison.C","Can't get the bachelor branch !");
176 return 4;
177 }
178 TClonesArray dum("AliTrackReference",1000), *brefs=&dum;
179 bbranch->SetAddress(&brefs);
180
181
182
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 !");
190 return 5;
191 }
192 }
4514157e 193 AliESDEvent* event = new AliESDEvent();
8b462fd8 194 TTree* esdTree = (TTree*) ef->Get("esdTree");
195 if (!esdTree) {
196 ::Error("AliCascadeComparison.C", "no ESD tree found");
197 return 6;
198 }
4514157e 199 event->ReadFromTree(esdTree);
4fa1e629 200
ca28c5f5 201
4fa1e629 202 //******* Loop over events *********
203 Int_t e=0;
8b462fd8 204 while (esdTree->GetEvent(e)) {
4fa1e629 205 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
206
4fa1e629 207 Int_t nentr=event->GetNumberOfCascades();
208 allfound+=nentr;
209
210
211 if (csTree->GetEvent(e++)==0) {
212 cerr<<"No reconstructable cascades !\n";
213 continue;
214 }
215
216 Int_t ngood=prefs->GetEntriesFast(),ng=0;
217
218 Double_t pxg=0.,pyg=0.,pzg=0.,ptg=0.;
219 Int_t nlab=-1, plab=-1, blab=-1;
220 Int_t i;
221 for (i=0; i<nentr; i++) {
222 AliESDcascade *cascade=event->GetCascade(i);
223
224 Int_t nidx=TMath::Abs(cascade->GetNindex());
225 Int_t pidx=TMath::Abs(cascade->GetPindex());
226 Int_t bidx=TMath::Abs(cascade->GetBindex());
227
228 AliESDtrack *ntrack=event->GetTrack(nidx);
229 AliESDtrack *ptrack=event->GetTrack(pidx);
230 AliESDtrack *btrack=event->GetTrack(bidx);
231
232 nlab=TMath::Abs(ntrack->GetLabel());
233 plab=TMath::Abs(ptrack->GetLabel());
234 blab=TMath::Abs(btrack->GetLabel());
235
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;
246 Double_t kine0;
247 Double_t kine=cascade->ChangeMassHypothesis(kine0,code);
248 if (TMath::Abs(kine0)>kine0cut) continue;
249 //if (TMath::Abs(kine)>kinecut) continue;
250
251 Double_t mass=cascade->GetEffMass();
252 cs->Fill(mass);
253 csf->Fill(mass);
254
255 AliTrackReference *nref=0, *pref=0, *bref=0;
256 Int_t j;
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;
264 }
265
266 if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
267
268 Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
269 Double_t pt=TMath::Sqrt(px*px+py*py);
270
271 if (j==ngood) {
272 hfake->Fill(pt);
273 cout<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
274 continue;
275 }
276 csf->Fill(mass,-1);
277
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.);
287
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);
292
293 hfound->Fill(ptg);
294 nref->SetLabel(-1);
295
296 }
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);
18b4df7a 301 Int_t pdg=(Int_t)nref->GetLength(); //this is the cascade's PDG !
4fa1e629 302 if (code!=pdg) continue;
303 ng++;
304 pxg=bref->Px()+nref->Px()+pref->Px();
305 pyg=bref->Px()+nref->Py()+pref->Py();
306 ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
307 hgood->Fill(ptg);
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";
311 }
312 allgood+=ng;
313
314 cout<<"Number of found cascades : "<<nentr<<endl;
315 cout<<"Number of \"good\" cascades : "<<ng<<endl;
316
317 brefs->Clear();
318 prefs->Clear();
319 nrefs->Clear();
4fa1e629 320
321 } //**** End of the loop over events
322
8b462fd8 323 delete event;
4514157e 324 delete esdTree;
4fa1e629 325 ef->Close();
ca28c5f5 326
4fa1e629 327 delete csTree;
328 refFile->Close();
329
330 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
331 if (ng!=0) cout<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
332 cout<<
333 "Total number of found cascades: "<<allfound<<" ("<<nf<<" in the peak)\n";
334 cout<<"Total number of \"good\" cascades: "<<allgood<<endl;
ca28c5f5 335
ca28c5f5 336
337 gStyle->SetOptStat(111110);
338 gStyle->SetOptFit(1);
339
4fa1e629 340 TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
341 if (!c1) {
342 c1=new TCanvas("c1","",0,0,580,610);
343 c1->Divide(2,2);
344 }
345
346 Int_t minc=33;
ca28c5f5 347
348 c1->cd(1);
349 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
4fa1e629 350 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus");
ca28c5f5 351 hl->Draw("same"); c1->cd();
352
353 c1->cd(2);
354 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
4fa1e629 355 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus");
356 c1->cd();
ca28c5f5 357
358 c1->cd(3);
359 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
4fa1e629 360 if (hx->GetEntries()<minc) hx->Draw(); else hx->Fit("gaus");
ca28c5f5 361 hy->Draw("same"); c1->cd();
362
363 c1->cd(4);
364 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
4fa1e629 365 if (hz->GetEntries()<minc) hz->Draw(); else hz->Fit("gaus");
ca28c5f5 366
4fa1e629 367 c1->Update();
ca28c5f5 368
4fa1e629 369 TCanvas *c2=(TCanvas*)gROOT->FindObject("c2");
370 if (!c2) {
371 c2=new TCanvas("c2","",600,0,580,610);
372 c2->Divide(1,2);
373 }
ca28c5f5 374
375 c2->cd(1);
376 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
377 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
378 hg->Divide(hfound,hgood,1,1.,"b");
379 hf->Divide(hfake,hgood,1,1.,"b");
380 hg->SetMaximum(1.4);
381 hg->SetYTitle("Cascade reconstruction efficiency");
382 hg->SetXTitle("Pt (GeV/c)");
383 hg->Draw();
384
385 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
386 line1->Draw("same");
387 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
388 line2->Draw("same");
389
390 hf->SetFillColor(1);
391 hf->SetFillStyle(3013);
392 hf->SetLineColor(2);
393 hf->SetLineWidth(2);
394 hf->Draw("histsame");
395 TText *text = new TText(0.461176,0.248448,"Fake cascades");
396 text->SetTextSize(0.05);
397 text->Draw();
398 text = new TText(0.453919,1.11408,"Good cascades");
399 text->SetTextSize(0.05);
400 text->Draw();
401
402
403 c2->cd(2);
404 gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
4fa1e629 405 if (cs->GetEntries()<minc) cs->Draw();
406 else cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
04b2a5f1 407 csf->Draw("same");
ca28c5f5 408 Double_t max=cs->GetMaximum();
409 TLine *line3 =
410 new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
411 line3->Draw("same");
412 TLine *line4 =
413 new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
414 line4->Draw("same");
415
4fa1e629 416 c2->Update();
ca28c5f5 417
4fa1e629 418 TFile fc("AliCascadeComparison.root","RECREATE");
419 c1->Write();
420 c2->Write();
421 fc.Close();
ca28c5f5 422
4fa1e629 423 gBenchmark->Stop("AliCascadeComparison");
424 gBenchmark->Show("AliCascadeComparison");
ca28c5f5 425
ca28c5f5 426
4fa1e629 427 return 0;
428}
ca28c5f5 429
4fa1e629 430
431Int_t GoodCascades(const Char_t *dir) {
5102bab6 432 if (gAlice) {
433 delete gAlice->GetRunLoader();
434 delete gAlice;
435 gAlice=0;
436 }
4fa1e629 437
438 Char_t fname[100];
439 sprintf(fname,"%s/galice.root",dir);
440
441 AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
5102bab6 442 if (!rl) {
4fa1e629 443 ::Error("GoodCascades","Can't start session !");
566abf75 444 return 1;
ca28c5f5 445 }
5102bab6 446
566abf75 447 rl->LoadgAlice();
448 rl->LoadHeader();
449 rl->LoadKinematics();
4fa1e629 450
451
452 Int_t nev=rl->GetNumberOfEvents();
453 ::Info("GoodCascades","Number of events : %d\n",nev);
454
455
456 sprintf(fname,"%s/GoodTracksITS.root",dir);
457 TFile *itsFile=TFile::Open(fname);
458 if ((!itsFile)||(!itsFile->IsOpen())) {
459 ::Error("GoodCAscades","Can't open the GoodTracksITS.root !");
460 delete rl;
461 return 5;
462 }
463 TClonesArray dm("AliTrackReference",1000), *itsRefs=&dm;
464 TTree *itsTree=(TTree*)itsFile->Get("itsTree");
465 if (!itsTree) {
466 ::Error("GoodCascades","Can't get the ITS reference tree !");
467 delete rl;
468 return 6;
469 }
470 TBranch *itsBranch=itsTree->GetBranch("ITS");
471 if (!itsBranch) {
472 ::Error("GoodCascades","Can't get the ITS reference branch !");
473 delete rl;
474 return 7;
475 }
476 itsBranch->SetAddress(&itsRefs);
477
478
479 sprintf(fname,"%s/GoodCascades.root",dir);
480 TFile *csFile=TFile::Open(fname,"recreate");
481 TClonesArray dummy("AliTrackReference",1000), *nrefs=&dummy;
482 TClonesArray dumm("AliTrackReference",1000), *prefs=&dumm;
483 TClonesArray dum("AliTrackReference",1000), *brefs=&dum;
484 TTree csTree("csTree","Tree with info about the reconstructable cascades");
485 csTree.Branch("negative",&nrefs);
486 csTree.Branch("positive",&prefs);
487 csTree.Branch("bachelor",&brefs);
488
489
490 // *** Get information about the cuts ***
491 Double_t r2min=0.9*0.9;
492 Double_t r2max=2.9*2.9;
493
494
495 //******** Loop over generated events
496 for (Int_t e=0; e<nev; e++) {
497 rl->GetEvent(e); csFile->cd();
498
499 Int_t np = rl->GetHeader()->GetNtrack();
500 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
501
502 itsTree->GetEvent(e);
503 Int_t nk=itsRefs->GetEntriesFast();
504
505 AliStack *stack=rl->Stack();
506
507 AliTrackReference *nref=0, *pref=0, *bref=0;
508
509 Int_t nc=0;
510 while (np--) {
511 //cerr<<np<<'\r';
512 TParticle *cp=stack->Particle(np);
513
514 // *** only these cascades are "good" ***
515 Int_t code=cp->GetPdgCode();
516 if (code!=kXiMinus) if (code!=kXiPlusBar)
517 if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue;
518
519 // *** daughter tracks must be "good" ***
520 Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
521 if (v0lab==blab) continue;
522 if (v0lab<0) continue;
523 if (blab<0) continue;
524
525 TParticle *p0=stack->Particle(v0lab);
526 TParticle *bp=stack->Particle(blab);
527 Int_t i;
528 if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
529 TParticle *p=p0; p0=bp; bp=p;
530 i=v0lab; v0lab=blab; blab=i;
531 if ((p0->GetPdgCode()!=kLambda0)&&(p0->GetPdgCode()!=kLambda0Bar))
ca28c5f5 532 continue;
4fa1e629 533 }
ca28c5f5 534
4fa1e629 535 // ** is the bachelor "good" ? **
536 for (i=0; i<nk; i++) {
537 bref=(AliTrackReference*)itsRefs->UncheckedAt(i);
538 if (bref->Label()==blab) break;
539 }
540 if (i==nk) continue;
541
542 // ** is the V0 "good" ? **
543 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
544 if (nlab==plab) continue;
545 if (nlab<0) continue;
546 if (plab<0) continue;
547 if (stack->Particle(plab)->GetPDG()->Charge() < 0.) {
548 i=plab; plab=nlab; nlab=i;
549 }
550 for (i=0; i<nk; i++) {
551 nref=(AliTrackReference*)itsRefs->UncheckedAt(i);
552 if (nref->Label()==nlab) break;
553 }
554 if (i==nk) continue;
555 for (i=0; i<nk; i++) {
556 pref=(AliTrackReference*)itsRefs->UncheckedAt(i);
557 if (pref->Label()==plab) break;
558 }
559 if (i==nk) continue;
560
561
562 // *** fiducial volume ***
563 Double_t x=bp->Vx(), y=bp->Vy(), r2=x*x+y*y; //bachelor
564 if (r2<r2min) continue;
565 if (r2>r2max) continue;
566 TParticle *pp=stack->Particle(plab);
567 x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y; //V0
568 if (r2<r2min) continue;
569 if (r2>r2max) continue;
570
571 Int_t pdg=cp->GetPdgCode();
572 nref->SetLength(pdg); //This will the cascade's PDG !
573
574 new((*nrefs)[nc]) AliTrackReference(*nref);
575 new((*prefs)[nc]) AliTrackReference(*pref);
576 new((*brefs)[nc]) AliTrackReference(*bref);
577
578 nc++;
ca28c5f5 579 }
4fa1e629 580 itsRefs->Clear();
ca28c5f5 581
4fa1e629 582 csTree.Fill();
583 nrefs->Clear(); prefs->Clear(); brefs->Clear();
ca28c5f5 584
4fa1e629 585 } //**** end of the loop over generated events
ca28c5f5 586
4fa1e629 587 csTree.Write();
588 csFile->Close();
ca28c5f5 589
4fa1e629 590 delete itsTree;
591 itsFile->Close();
592
593 delete rl;
594 return 0;
ca28c5f5 595}