]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/UPGRADE/macros/AliTPCUComparison.C
Fixing a memory leak
[u/mrichter/AliRoot.git] / ITS / UPGRADE / macros / AliTPCUComparison.C
CommitLineData
5bc07054 1/****************************************************************************
2 * Legacy comparison macro, adapted to the upgraded ITS. *
3 * *
4 * Creates list of "trackable" tracks, *
5 * calculates efficiency, resolutions etc. *
6 * There is a possibility to run this macro over several events. *
7 * *
8 * Origin: I.Belikov, IPHC, Iouri.Belikov@iphc.cnrs.fr *
9 * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de *
10 ****************************************************************************/
11
12#if !defined(__CINT__) || defined(__MAKECINT__)
13 #include <TMath.h>
14 #include <TError.h>
15 #include <Riostream.h>
16 #include <TH1F.h>
17 #include <TH2F.h>
18 #include <TTree.h>
19 #include <TParticle.h>
20 #include <TCanvas.h>
21 #include <TLine.h>
22 #include <TText.h>
23 #include <TBenchmark.h>
24 #include <TStyle.h>
25 #include <TFile.h>
26 #include <TROOT.h>
27 #include <TGeoGlobalMagField.h>
28
29 #include "AliStack.h"
30 #include "AliHeader.h"
31 #include "AliTrackReference.h"
32 #include "AliRunLoader.h"
33 #include "AliRun.h"
34 #include "AliMagF.h"
35 #include "AliESDEvent.h"
36 #include "AliESDtrack.h"
37
38 #include "Base/AliSimDigits.h"
39 #include "Base/AliTPCParamSR.h"
40 #include "Base/AliTPCLoader.h"
41 #include "Base/AliTPCcalibDB.h"
42 #include "Sim/AliTPC.h"
43 #include "Rec/AliTPCClustersRow.h"
44
45 #include "AliCDBManager.h"
46#endif
47
48Int_t GoodTracksTPC(const Char_t *dir=".");
49
50extern TBenchmark *gBenchmark;
51extern TROOT *gROOT;
52
53static Int_t allgood=0;
54static Int_t allselected=0;
55static Int_t allfound=0;
56
57Int_t AliTPCUComparison
6ec85964 58(Float_t ptcutl=0., Float_t ptcuth=2., const Char_t *dir=".") {
5bc07054 59 gBenchmark->Start("AliTPCUComparison");
60
61 ::Info("AliTPCUComparison.C","Doing comparison...");
62
63
64
65 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
66 if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.);
67 hp->SetFillColor(4);
68
69 TH1F *hl=(TH1F*)gROOT->FindObject("hl");
70 if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
71 hl->SetFillColor(4);
72
73 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
74 if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
75 hpt->SetFillColor(2);
76
77 TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
78 if (!hmpt)
79 hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
80 hmpt->SetFillColor(6);
81
82
6ec85964 83 Int_t nb=100;
5bc07054 84 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
6ec85964 85 if (!hgood) hgood=new TH1F("hgood","Good tracks",nb,ptcutl,ptcuth);
5bc07054 86
87 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
6ec85964 88 if (!hfound) hfound=new TH1F("hfound","Found tracks",nb,ptcutl,ptcuth);
5bc07054 89
90 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
6ec85964 91 if (!hfake) hfake=new TH1F("hfake","Fake tracks",nb,ptcutl,ptcuth);
5bc07054 92
93 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
6ec85964 94 if (!hg) hg=new TH1F("hg","Efficiency for good tracks",nb,ptcutl,ptcuth);
5bc07054 95 hg->SetLineColor(4); hg->SetLineWidth(2);
96
97 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
6ec85964 98 if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",nb,ptcutl,ptcuth);
5bc07054 99 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
100
6ec85964 101
5bc07054 102 TH1F *he=(TH1F*)gROOT->FindObject("he");
103 if (!he)
104 he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
105
106 TH2F *hep=(TH2F*)gROOT->FindObject("hep");
107 if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
108 hep->SetMarkerStyle(8);
109 hep->SetMarkerSize(0.4);
110
111
112 Char_t fname[100];
113 sprintf(fname,"%s/GoodTracksTPC.root",dir);
114
115 TFile *refFile=TFile::Open(fname,"old");
116 if (!refFile || !refFile->IsOpen()) {
117 ::Info("AliTPCUComparison.C","Marking good tracks (will take a while)...");
118 if (GoodTracksTPC(dir)) {
119 ::Error("AliTPCUComparison.C","Can't generate the reference file !");
120 return 1;
121 }
122 }
123 refFile=TFile::Open(fname,"old");
124 if (!refFile || !refFile->IsOpen()) {
125 ::Error("AliTPCUComparison.C","Can't open the reference file !");
126 return 2;
127 }
128
129 TTree *tpcTree=(TTree*)refFile->Get("tpcTree");
130 if (!tpcTree) {
131 ::Error("AliTPCUComparison.C","Can't get the reference tree !");
132 return 3;
133 }
134 TBranch *branch=tpcTree->GetBranch("TPC");
135 if (!branch) {
136 ::Error("AliTPCUComparison.C","Can't get the TPC branch !");
137 return 4;
138 }
139 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
140 branch->SetAddress(&refs);
141
142
143 sprintf(fname,"%s/AliESDs.root",dir);
144 TFile *ef=TFile::Open(fname);
145 if ((!ef)||(!ef->IsOpen())) {
146 sprintf(fname,"%s/AliESDtpc.root",dir);
147 ef=TFile::Open(fname);
148 if ((!ef)||(!ef->IsOpen())) {
149 ::Error("AliTPCUComparison.C","Can't open AliESDtpc.root !");
150 return 5;
151 }
152 }
153 AliESDEvent* event = new AliESDEvent();
154 TTree* esdTree = (TTree*) ef->Get("esdTree");
155 if (!esdTree) {
156 ::Error("AliTPCUComparison.C", "no ESD tree found");
157 return 6;
158 }
159 event->ReadFromTree(esdTree);
160
161
162 //******* Loop over events *********
163 Int_t e=0;
164 while (esdTree->GetEvent(e)) {
165 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
166
167 Int_t nentr=event->GetNumberOfTracks();
168 allfound+=nentr;
169
170 if (tpcTree->GetEvent(e++)==0) {
171 cerr<<"No reconstructable tracks !\n";
172 continue;
173 }
174
175 Int_t ngood=refs->GetEntriesFast();
176 allgood+=ngood;
177
178 const Int_t MAX=15000;
179 //MI change
180 Int_t track_notfound[MAX], itrack_notfound=0;
181 Int_t track_fake[MAX], itrack_fake=0;
182 Int_t track_multifound[MAX],track_multifound_n[MAX],itrack_multifound=0;
183
184 Int_t i;
185 for (Int_t k=0; k<ngood; k++) {
186 AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k);
187 Int_t lab=ref->Label(), tlab=-1;
188 Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());
189
190 if (ptg<1e-33) continue; // for those not crossing 0 pad row
191
192 if (ptg<ptcutl) continue;
193 if (ptg>ptcuth) continue;
194
195 allselected++;
196
197 hgood->Fill(ptg);
198
199 AliESDtrack *track=0;
200 for (i=0; i<nentr; i++) {
201 track=event->GetTrack(i);
202 tlab=track->GetTPCLabel();
203 if (lab==TMath::Abs(tlab)) break;
204 }
205 if (i==nentr) {
206 track_notfound[itrack_notfound++]=lab;
207 continue;
208 }
209
210 //MI change - addition
211 Int_t micount=0;
212 Int_t mi;
213 AliESDtrack * mitrack;
214 for (mi=0; mi<nentr; mi++) {
215 mitrack=event->GetTrack(mi);
216 if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++;
217 }
218 if (micount>1) {
219 track_multifound[itrack_multifound]=lab;
220 track_multifound_n[itrack_multifound]=micount;
221 itrack_multifound++;
222 }
223 if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
224 if (lab==tlab) hfound->Fill(ptg);
225 else {
226 track_fake[itrack_fake++]=lab;
227 hfake->Fill(ptg);
228 }
229
230 Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz);
231 Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]);
232 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
233 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
234 Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]);
235 Float_t lam=TMath::ATan2(pxpypz[2],pt);
236 Float_t pt_1=1/pt;
237
238 Int_t pdg=(Int_t)ref->GetLength(); //this is particle's PDG !
239
240 if (TMath::Abs(pdg)==11 && ptg>4.) {//high momentum electrons
241 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
242 } else {
243 Float_t phig=TMath::ATan2(ref->Py(),ref->Px());
244 hp->Fill((phi - phig)*1000.);
245
246 Float_t lamg=TMath::ATan2(ref->Pz(),ptg);
247 hl->Fill((lam - lamg)*1000.);
248
249 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
250 }
251
252 Float_t mom=pt/TMath::Cos(lam);
253 Float_t dedx=track->GetTPCsignal();
254 hep->Fill(mom,dedx,1.);
255 if (TMath::Abs(pdg)==211) //pions
256 if (mom>0.4 && mom<0.5) {
257 he->Fill(dedx,1.);
258 }
259 }
260
261 cout<<"\nList of Not found tracks :\n";
262 for ( i = 0; i< itrack_notfound; i++){
263 cout<<track_notfound[i]<<"\t";
264 if ((i%5)==4) cout<<"\n";
265 }
266 cout<<"\nList of fake tracks :\n";
267 for ( i = 0; i< itrack_fake; i++){
268 cout<<track_fake[i]<<"\t";
269 if ((i%5)==4) cout<<"\n";
270 }
271 cout<<"\nList of multiple found tracks :\n";
272 for ( i=0; i<itrack_multifound; i++) {
273 cout<<"id. "<<track_multifound[i]
274 <<" found - "<<track_multifound_n[i]<<"times\n";
275 }
276
277 cout<<"Number of found tracks ="<<nentr<<endl;
278 cout<<"Number of \"good\" tracks ="<<ngood<<endl;
279
280 refs->Clear();
281 }// ***** End of the loop over events
282
283 delete event;
284 delete esdTree;
285 ef->Close();
286
287 delete tpcTree;
288 refFile->Close();
289
290 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
291 if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
292 cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl;
293 cout<<"Total number of found tracks ="<<allfound<<endl;
294 cout<<"Total number of \"good\" tracks ="<<allgood<<endl;
295 cout<<endl;
296
297 gStyle->SetOptStat(111110);
298 gStyle->SetOptFit(1);
299
300 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
301
302 Int_t minc=33;
303
304 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
305 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
306 hp->SetFillColor(4); hp->SetXTitle("(mrad)");
307 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();
308
309 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
310 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
311 hl->SetXTitle("(mrad)");
312 if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();
313
314 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
315 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
316 hpt->SetXTitle("(%)");
317 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();
318
319 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
320 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
321 hmpt->SetXTitle("(%)");
322 if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd();
323
324 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
325 p5->SetFillColor(41); p5->SetFrameFillColor(10);
326 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
327 hg->Divide(hfound,hgood,1,1.,"b");
328 hf->Divide(hfake,hgood,1,1.,"b");
329 hg->SetMaximum(1.4);
330 hg->SetYTitle("Tracking efficiency");
331 hg->SetXTitle("Pt (GeV/c)");
332 hg->Draw();
333
6ec85964 334 TLine *line1 = new TLine(ptcutl,1.0,ptcuth,1.0); line1->SetLineStyle(4);
5bc07054 335 line1->Draw("same");
6ec85964 336 TLine *line2 = new TLine(ptcutl,0.9,ptcuth,0.9); line2->SetLineStyle(4);
5bc07054 337 line2->Draw("same");
338
339 hf->SetFillColor(1);
340 hf->SetFillStyle(3013);
341 hf->SetLineColor(2);
342 hf->SetLineWidth(2);
343 hf->Draw("histsame");
344 TText *text = new TText(0.461176,0.248448,"Fake tracks");
345 text->SetTextSize(0.05);
346 text->Draw();
347 text = new TText(0.453919,1.11408,"Good tracks");
348 text->SetTextSize(0.05);
349 text->Draw();
350
351
352
353 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
354
355 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
356 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
357 he->SetFillColor(2); he->SetFillStyle(3005);
358 he->SetXTitle("Arbitrary Units");
359 if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();
360
361 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
362 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
363 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
364 hep->Draw(); c1->cd();
365
366 TFile fc("AliTPCUComparison.root","RECREATE");
367 c1->Write();
368 c2->Write();
369 fc.Close();
370
371 gBenchmark->Stop("AliTPCUComparison");
372 gBenchmark->Show("AliTPCUComparison");
373
374 return 0;
375}
376
377
378Int_t GoodTracksTPC(const Char_t *dir) {
379 Char_t fname[100];
380 sprintf(fname,"%s/galice.root",dir);
381
382 AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
383 if (!rl) {
384 ::Error("GoodTracksTPC","Can't start session !");
385 return 1;
386 }
387
388 rl->LoadgAlice();
389 rl->LoadHeader();
390 rl->LoadKinematics();
391 rl->LoadTrackRefs();
392
393 AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
394 if (tpcl == 0x0) {
395 ::Error("GoodTracksTPC","Can not find TPCLoader !");
396 delete rl;
397 return 2;
398 }
399
400 AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
401 Int_t ver = TPC->IsVersion();
402 cout<<"TPC version "<<ver<<" has been found !\n";
403 if (ver==1) tpcl->LoadRecPoints();
404 else if (ver==2) tpcl->LoadDigits();
405 else {
406 ::Error("GoodTracksTPC","Wrong TPC version !");
407 delete rl;
408 return 3;
409 }
410
411 TGeoGlobalMagField::Instance()->SetField(
412 new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG));
413
414 AliCDBManager *man=AliCDBManager::Instance();
415 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
416 man->SetRun(0);
417 AliTPCParamSR *digp=
418 (AliTPCParamSR*)(AliTPCcalibDB::Instance()->GetParameters());
419 if (!digp) {
420 ::Error("AliTPCUComparison.C","TPC parameters have not been found !");
421 delete rl;
422 return 4;
423 }
424
425 Int_t nrow_up=digp->GetNRowUp();
426 Int_t nrows=digp->GetNRowLow()+nrow_up;
427 Int_t zero=digp->GetZeroSup();
428 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
429 Int_t good_number=Int_t(0.4*nrows);
430
431 Int_t nev=rl->GetNumberOfEvents();
432 ::Info("GoodTracksTPC","Number of events : %d\n",nev);
433
434
435 sprintf(fname,"%s/GoodTracksTPC.root",dir);
436 TFile *file=TFile::Open(fname,"recreate");
437
438 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
439 TTree tpcTree("tpcTree","Tree with info about the reconstructable TPC tracks");
440 tpcTree.Branch("TPC",&refs);
441
442 //******** Loop over generated events
443 for (Int_t e=0; e<nev; e++) {
444 Int_t nt=0;
445 refs->Clear();
446
447 Int_t i;
448
449 rl->GetEvent(e); file->cd();
450
451 Int_t np = rl->GetHeader()->GetNtrack();
452 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
453
454 //******** Fill the "good" masks
455 Int_t *good=new Int_t[np]; for (i=0; i<np; i++) good[i]=0;
456 switch (ver) {
457 case 1:
458 /*
459 {
460 AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
461 ca.Setup(digp);
462 ca.SetClusterType("AliTPCcluster");
463 ca.ConnectTree(tpcl->TreeR());
464 Int_t nrows=Int_t(ca.GetTree()->GetEntries());
465 for (Int_t n=0; n<nrows; n++) {
466 AliSegmentID *s=ca.LoadEntry(n);
467 Int_t sec,row;
468 digp->AdjustSectorRow(s->GetID(),sec,row);
469 AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
470 Int_t ncl=clrow.GetArray()->GetEntriesFast();
471 while (ncl--) {
472 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
473 Int_t lab=c->GetLabel(0);
474 if (lab<0) continue; //noise cluster
475 lab=TMath::Abs(lab);
476
477 if (sec>=digp->GetNInnerSector())
478 if (row==nrow_up-1) good[lab]|=0x4000;
479 if (sec>=digp->GetNInnerSector())
480 if (row==nrow_up-1-gap) good[lab]|=0x1000;
481
482 if (sec>=digp->GetNInnerSector())
483 if (row==nrow_up-1-shift) good[lab]|=0x2000;
484 if (sec>=digp->GetNInnerSector())
485 if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
486
487 good[lab]++;
488 }
489 ca.ClearRow(sec,row);
490 }
491 delete pca;
492 }
493 break;
494 */
495 case 2:
496 {
497 TTree *TD=tpcl->TreeD();
498
499 AliSimDigits da, *digits=&da;
500 TD->GetBranch("Segment")->SetAddress(&digits);
501
502 Int_t *count = new Int_t[np];
503 Int_t i;
504 for (i=0; i<np; i++) count[i]=0;
505
506 Int_t sectors_by_rows=(Int_t)TD->GetEntries();
507 for (i=0; i<sectors_by_rows; i++) {
508 if (!TD->GetEvent(i)) continue;
509 Int_t sec,row;
510 digp->AdjustSectorRow(digits->GetID(),sec,row);
511 //cerr<<sec<<' '<<row<<'\r';
512 digits->First();
513 do { //Many thanks to J.Chudoba who noticed this
514 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
515 Short_t dig = digits->GetDigit(it,ip);
516 Int_t idx0=digits->GetTrackID(it,ip,0);
517 Int_t idx1=digits->GetTrackID(it,ip,1);
518 Int_t idx2=digits->GetTrackID(it,ip,2);
519 if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
520 if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
521 if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
522 } while (digits->Next());
523 for (Int_t j=0; j<np; j++) {
524 if (count[j]>1) {
525 if (sec>=digp->GetNInnerSector())
526 if (row==nrow_up-1 ) good[j]|=0x4000;
527 if (sec>=digp->GetNInnerSector())
528 if (row==nrow_up-1-gap) good[j]|=0x1000;
529
530 if (sec>=digp->GetNInnerSector())
531 if (row==nrow_up-1-shift) good[j]|=0x2000;
532 if (sec>=digp->GetNInnerSector())
533 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
534 good[j]++;
535 }
536 count[j]=0;
537 }
538 }
539 delete[] count;
540 }
541 break;
542 }
543
544
545 //****** select tracks which are "good" enough
546 AliStack* stack = rl->Stack();
547 for (i=0; i<np; i++) {
548 if ((good[i]&0x5000) != 0x5000)
549 if ((good[i]&0x2800) != 0x2800) continue;
550 if ((good[i]&0x7FF ) < good_number) continue;
5bc07054 551 TParticle *p = (TParticle*)stack->Particle(i);
552 if (p == 0x0) {
553 cerr<<"Can not get particle "<<i<<endl;
554 continue;
555 }
6ec85964 556 if (p->Pt() <= 0.) continue;
5bc07054 557 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
558
559 Double_t vx=p->Vx(),vy=p->Vy(),vz=p->Vz();
560 if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue;
561 if (TMath::Abs(vz) > 50.) continue;
562
563 AliTrackReference *ref=new((*refs)[nt]) AliTrackReference();
564
565 ref->SetLabel(i);
566 Int_t pdg=p->GetPdgCode();
567 ref->SetLength(pdg); //This will the particle's PDG !
568 ref->SetMomentum(0.,0.,0.);
569 ref->SetPosition(0.,0.,0.);
570
571 nt++;
572 }
573
574 //**** check if there is also information at the entrance of the TPC
575 TTree *TR=rl->TreeTR();
576 TBranch *branch=TR->GetBranch("TrackReferences");
577 if (branch==0) {
578 ::Error("GoodTracksTPC","No track references !");
579 delete rl;
580 return 5;
581 }
582 TClonesArray tpcdummy("AliTrackReference",1000), *tpcRefs=&tpcdummy;
583 branch->SetAddress(&tpcRefs);
584 Int_t nr=(Int_t)TR->GetEntries();
585 for (Int_t r=0; r<nr; r++) {
586 //cerr<<r<<' '<<nr<<'\r';
587 TR->GetEvent(r);
588
589 Int_t nref = tpcRefs->GetEntriesFast();
590 if (!nref) continue;
591 AliTrackReference *tpcRef= 0x0;
592 for (Int_t iref=0; iref<nref; ++iref) {
593 tpcRef = (AliTrackReference*)tpcRefs->UncheckedAt(iref);
594 if (tpcRef->DetectorId() == AliTrackReference::kTPC) break;
595 tpcRef = 0x0;
596 }
597
598 if (!tpcRef) continue;
599
600 Int_t j;
601 AliTrackReference *ref=0;
602 for (j=0; j<nt; j++) {
603 ref=(AliTrackReference *)refs->UncheckedAt(j);
604 if (ref->Label()==tpcRef->Label()) break;
605 ref=0;
606 }
607 if (ref==0) continue;
608
609 ref->SetMomentum(tpcRef->Px(),tpcRef->Py(),tpcRef->Pz());
610 ref->SetPosition(tpcRef->LocalX(),tpcRef->LocalY(),tpcRef->Z());
611
612 tpcRefs->Clear();
613 }
614
615 tpcTree.Fill();
616
617 delete[] good;
618 } //****** end of the loop over generated events
619
620
621 tpcTree.Write();
622 file->Close();
623
624 delete rl;
625 return 0;
626}