]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliCascadeComparison.C
Bug fix (Yu.Belikov)
[u/mrichter/AliRoot.git] / ITS / 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 #if !defined(__CINT__) || defined(__MAKECINT__)
11   #include <TMath.h>
12   #include <TError.h>
13   #include <Riostream.h>
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>
22   #include <TKey.h>
23   #include <TROOT.h>
24
25   #include "AliStack.h"
26   #include "AliHeader.h"
27   #include "AliTrackReference.h"
28   #include "AliRunLoader.h"
29   #include "AliRun.h"
30   #include "AliESD.h"
31 #else
32 const Int_t kXiMinus    = 3312;
33 const Int_t kXiPlusBar  = -3312;
34 const Int_t kOmegaMinus = 3334;
35 const Int_t kOmegaPlusBar = -3334;
36 #endif
37
38 Int_t GoodCascades(const Char_t *dir=".");
39
40 extern AliRun *gAlice;
41 extern TBenchmark *gBenchmark;
42 extern TROOT *gROOT;
43
44 static Int_t allgood=0;
45 static Int_t allfound=0;
46
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");
53
54    ::Info("AliCascadeComparison.C","Doing comparison...");
55
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;
60    switch (code) {
61    case kXiMinus:
62         break;
63    case kXiPlusBar:
64         ptncut=0.33; ptpcut=0.12;
65         break;
66    case kOmegaMinus: 
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;
77    default: cerr<<"Invalid PDG code !\n"; return 1;
78    }
79
80
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);
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.); 
91    hpt->SetXTitle("(%)"); hpt->SetFillColor(2); 
92
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);
96
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);
100
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);
104
105
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);
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);
118    hg->SetLineColor(4); hg->SetLineWidth(2);
119
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);
123
124
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);
129
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);
133
134
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      }
145    }
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;
156    }
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    }
193    TKey *key=0;
194    TIter next(ef->GetListOfKeys());
195
196
197
198    //******* Loop over events *********
199    Int_t e=0;
200    while ((key=(TKey*)next())!=0) {
201       cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
202  
203       AliESD *event=(AliESD*)key->ReadObj();
204
205       Int_t nentr=event->GetNumberOfCascades();
206       allfound+=nentr;
207
208
209       if (csTree->GetEvent(e++)==0) {
210          cerr<<"No reconstructable cascades !\n";
211          continue;
212       }
213
214       Int_t ngood=prefs->GetEntriesFast(),ng=0; 
215
216       Double_t pxg=0.,pyg=0.,pzg=0.,ptg=0.;
217       Int_t nlab=-1, plab=-1, blab=-1;
218       Int_t i;
219       for (i=0; i<nentr; i++) {
220           AliESDcascade *cascade=event->GetCascade(i);
221
222           Int_t nidx=TMath::Abs(cascade->GetNindex());
223           Int_t pidx=TMath::Abs(cascade->GetPindex());
224           Int_t bidx=TMath::Abs(cascade->GetBindex());
225
226           AliESDtrack *ntrack=event->GetTrack(nidx);
227           AliESDtrack *ptrack=event->GetTrack(pidx);
228           AliESDtrack *btrack=event->GetTrack(bidx);
229
230           nlab=TMath::Abs(ntrack->GetLabel()); 
231           plab=TMath::Abs(ptrack->GetLabel());
232           blab=TMath::Abs(btrack->GetLabel());
233
234           /** Kinematical cuts **/
235           Double_t pxn,pyn,pzn; cascade->GetNPxPyPz(pxn,pyn,pzn); 
236           Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
237           if (ptn < ptncut) continue;
238           Double_t pxp,pyp,pzp; cascade->GetPPxPyPz(pxp,pyp,pzp); 
239           Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
240           if (ptp < ptpcut) continue;
241           Double_t pxb,pyb,pzb; cascade->GetBPxPyPz(pxb,pyb,pzb); 
242           Double_t ptb=TMath::Sqrt(pxb*pxb + pyb*pyb);
243           if (ptb < ptbcut) continue;
244           Double_t kine0;
245           Double_t kine=cascade->ChangeMassHypothesis(kine0,code);
246           if (TMath::Abs(kine0)>kine0cut) continue;
247           //if (TMath::Abs(kine)>kinecut) continue;
248
249           Double_t mass=cascade->GetEffMass();
250           cs->Fill(mass);
251           csf->Fill(mass);
252
253           AliTrackReference *nref=0, *pref=0, *bref=0;
254           Int_t j;
255           for (j=0; j<ngood; j++) {
256               bref=(AliTrackReference*)brefs->UncheckedAt(j);
257               nref=(AliTrackReference*)nrefs->UncheckedAt(j);
258               pref=(AliTrackReference*)prefs->UncheckedAt(j);
259               if (bref->Label() == blab)
260               if (nref->Label() == nlab)
261               if (pref->Label() == plab) break;
262           }
263
264           if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
265
266           Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
267           Double_t pt=TMath::Sqrt(px*px+py*py);
268
269           if (j==ngood) {
270              hfake->Fill(pt);
271              cout<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
272              continue;
273           }
274           csf->Fill(mass,-1);
275
276           pxg=bref->Px()+nref->Px()+pref->Px(); 
277           pyg=bref->Px()+nref->Py()+pref->Py(); 
278           pzg=nref->Pz()+pref->Pz(); 
279           ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
280           Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
281           Double_t lamg=TMath::ATan2(pzg,ptg), lam=TMath::ATan2(pz,pt);
282           hp->Fill((phi - phig)*1000.);
283           hl->Fill((lam - lamg)*1000.);
284           hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
285
286           Double_t x,y,z; cascade->GetXYZ(x,y,z);
287           hx->Fill((x-nref->X())*10);
288           hy->Fill((y-nref->Y())*10);
289           hz->Fill((z-nref->Z())*10);
290
291           hfound->Fill(ptg);
292           nref->SetLabel(-1);
293
294       }
295       for (i=0; i<ngood; i++) {
296          AliTrackReference *bref=(AliTrackReference*)brefs->UncheckedAt(i);
297          AliTrackReference *nref=(AliTrackReference*)nrefs->UncheckedAt(i);
298          AliTrackReference *pref=(AliTrackReference*)prefs->UncheckedAt(i);
299          Int_t pdg=(Int_t)nref->GetLength();  //this is the cascade's PDG !
300          if (code!=pdg) continue;
301          ng++;
302          pxg=bref->Px()+nref->Px()+pref->Px(); 
303          pyg=bref->Px()+nref->Py()+pref->Py(); 
304          ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
305          hgood->Fill(ptg);
306          nlab=nref->Label(); plab=pref->Label(); blab=bref->Label();
307          if (nlab < 0) continue;
308          cout<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
309       }
310       allgood+=ng;
311
312       cout<<"Number of found cascades : "<<nentr<<endl;
313       cout<<"Number of \"good\" cascades : "<<ng<<endl;
314
315       brefs->Clear();
316       prefs->Clear();
317       nrefs->Clear();
318       delete event;
319
320    } //**** End of the loop over events
321
322    ef->Close();
323
324    delete csTree;
325    refFile->Close();
326
327    Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
328    if (ng!=0) cout<<"Integral efficiency is about "<<nf/ng*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 gAlice->GetRunLoader();
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 }