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