]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCComparison.C
Comparison working with ESD (Yu.Belikov)
[u/mrichter/AliRoot.git] / TPC / AliTPCComparison.C
CommitLineData
7f6ddf58 1/****************************************************************************
2 * Very important, delicate and rather obscure macro. *
3 * *
4 * Creates list of "trackable" tracks, *
7f6ddf58 5 * calculates efficiency, resolutions etc. *
5102bab6 6 * The ESD tracks must be in an appropriate state (kTPCin or kTPCrefit) *
7f6ddf58 7 * *
024a7fe9 8 * Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch *
7f6ddf58 9 * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de *
10 ****************************************************************************/
11
c630aafd 12#if !defined(__CINT__) || defined(__MAKECINT__)
5102bab6 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 <TPad.h>
21 #include <TCanvas.h>
22 #include <TLine.h>
23 #include <TText.h>
24 #include <TBenchmark.h>
25 #include <TStyle.h>
26 #include <TKey.h>
27
28 #include "AliStack.h"
29 #include "AliHeader.h"
30 #include "AliTrackReference.h"
88cb7938 31 #include "AliRunLoader.h"
32 #include "AliMagF.h"
c630aafd 33 #include "AliRun.h"
5102bab6 34 #include "AliESD.h"
35
36 #include "AliSimDigits.h"
37 #include "AliTPC.h"
38 #include "AliTPCParamSR.h"
39 #include "AliTPCClustersArray.h"
40 #include "AliTPCClustersRow.h"
41 #include "AliTPCtracker.h"
42 #include "AliTPCLoader.h"
be9c5115 43#endif
44
f38c8ae5 45struct GoodTrackTPC {
cc80f89e 46 Int_t lab;
47 Int_t code;
48 Float_t px,py,pz;
49 Float_t x,y,z;
50};
cc80f89e 51
5102bab6 52Int_t good_tracks_tpc(GoodTrackTPC *gt, const Int_t max,
53const char* evfoldname=AliConfig::fgkDefaultEventFolderName);
88cb7938 54
c630aafd 55extern AliRun *gAlice;
5102bab6 56extern TBenchmark *gBenchmark;
fbb58996 57
5102bab6 58Int_t AliTPCComparison() {
59 gBenchmark->Start("AliTPCComparison");
60
61 ::Info("AliTPCComparison.C","Doing comparison...");
c630aafd 62
63 if (gAlice) {
88cb7938 64 delete gAlice->GetRunLoader();
65 delete gAlice;//if everything was OK here it is already NULL
66 gAlice = 0x0;
67 }
68
88cb7938 69 AliRunLoader *rl = AliRunLoader::Open("galice.root","COMPARISON");
c630aafd 70 if (!rl) {
5102bab6 71 ::Fatal("AliTPCComparison.C","Can't start session !");
c630aafd 72 }
7f6ddf58 73
88cb7938 74 rl->LoadgAlice();
3cb3ee52 75 if (rl->GetAliRun())
c630aafd 76 AliKalmanTrack::SetConvConst(
77 1000/0.299792458/rl->GetAliRun()->Field()->SolenoidField()
78 );
79 else {
5102bab6 80 ::Fatal("AliTPCComparison.C","Can't get AliRun !");
c630aafd 81 }
3cb3ee52 82 //rl->UnloadgAlice();
88cb7938 83
afc42102 84
7f6ddf58 85 /* Generate a list of "good" tracks */
5102bab6 86 const Int_t MAX=20000;
f38c8ae5 87 GoodTrackTPC gt[MAX];
cc80f89e 88 Int_t ngood=0;
be9c5115 89 ifstream in("good_tracks_tpc");
cc80f89e 90 if (in) {
5102bab6 91 ::Info("AliTPCComparison.C","Reading good tracks...");
7f6ddf58 92 while (in>>gt[ngood].lab>>gt[ngood].code>>
be9c5115 93 gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
94 gt[ngood].x >>gt[ngood].y >>gt[ngood].z) {
7ac477db 95 ngood++;
7ac477db 96 if (ngood==MAX) {
5102bab6 97 ::Warning("AliTPCComparison.C","Too many good tracks !");
7ac477db 98 break;
99 }
cc80f89e 100 }
5102bab6 101 if (!in.eof())
102 ::Fatal("AliTPCComparison.C","Read error (good_tracks_tpc) !");
cc80f89e 103 } else {
5102bab6 104 ::Info
105 ("AliTPCComparison","Marking good tracks (this will take a while)...");
88cb7938 106 ngood=good_tracks_tpc(gt,MAX,"COMPARISON");
be9c5115 107 ofstream out("good_tracks_tpc");
cc80f89e 108 if (out) {
109 for (Int_t ngd=0; ngd<ngood; ngd++)
7f6ddf58 110 out<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
be9c5115 111 gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
112 gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
5102bab6 113 } else
114 ::Fatal("AliTPCComparison.C","Can not open file (good_tracks_tpc) !");
cc80f89e 115 out.close();
116 }
7f6ddf58 117
5102bab6 118 Int_t i=0; TObjArray tarray(MAX);
3cb3ee52 119 { /*Load tracks*/
5102bab6 120 TFile *ef=TFile::Open("AliESDtpc.root");
121 if ((!ef)||(!ef->IsOpen())) {
122 ::Fatal("AliTPCComparison.C","Can't open AliESDtpc.root !");
123 }
124 TKey *key=0;
125 TIter next(ef->GetListOfKeys());
126 if ((key=(TKey*)next())!=0) {
127 AliESD *event=(AliESD*)key->ReadObj();
128 Int_t ntrk=event->GetNumberOfTracks();
129 for (Int_t i=0; i<ntrk; i++) {
130 AliESDtrack *t=event->GetTrack(i);
131 if ((t->GetStatus()&AliESDtrack::kTPCin)==0) continue;
132 AliTPCtrack *iotrack=new AliTPCtrack(*t);
133 tarray.AddLast(iotrack);
134 }
135 delete event;
136 }
137 ef->Close();
3cb3ee52 138 }
5102bab6 139 Int_t nentr=tarray.GetEntriesFast();
7f6ddf58 140
cc80f89e 141 TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4);
142 TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4);
8c555625 143 TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
144 hpt->SetFillColor(2);
cc80f89e 145 TH1F *hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
146 hmpt->SetFillColor(6);
8c555625 147
5102bab6 148 TH1F *hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
149 TH1F *hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);
150 TH1F *hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);
151 TH1F *hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
8c555625 152 hg->SetLineColor(4); hg->SetLineWidth(2);
5102bab6 153 TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
8c555625 154 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
155
cc80f89e 156 TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
7f6ddf58 157 TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
be9c5115 158 hep->SetMarkerStyle(8);
159 hep->SetMarkerSize(0.4);
8c555625 160
7f6ddf58 161 //MI change
162 Int_t mingood=ngood;
163 Int_t track_notfound[MAX], itrack_notfound=0;
164 Int_t track_fake[MAX], itrack_fake=0;
165 Int_t track_multifound[MAX], track_multifound_n[MAX], itrack_multifound=0;
cf98c13f 166
cc80f89e 167 while (ngood--) {
cf98c13f 168 Int_t lab=gt[ngood].lab,tlab=-1;
cc80f89e 169 Float_t ptg=
170 TMath::Sqrt(gt[ngood].px*gt[ngood].px + gt[ngood].py*gt[ngood].py);
8c555625 171
7f6ddf58 172 if (ptg<1e-33) continue; // for those not crossing 0 pad row
b9de75e1 173
8c555625 174 hgood->Fill(ptg);
cc80f89e 175
cf98c13f 176 AliTPCtrack *track=0;
7f6ddf58 177 for (i=0; i<nentr; i++) {
cf98c13f 178 track=(AliTPCtrack*)tarray.UncheckedAt(i);
73042f01 179 tlab=track->GetLabel();
cc80f89e 180 if (lab==TMath::Abs(tlab)) break;
8c555625 181 }
7f6ddf58 182 if (i==nentr) {
cf98c13f 183 track_notfound[itrack_notfound++]=lab;
cc80f89e 184 continue;
185 }
cf98c13f 186
187 //MI change - addition
188 Int_t micount=0;
189 Int_t mi;
190 AliTPCtrack * mitrack;
7f6ddf58 191 for (mi=0; mi<nentr; mi++) {
cf98c13f 192 mitrack=(AliTPCtrack*)tarray.UncheckedAt(mi);
193 if (lab==TMath::Abs(mitrack->GetLabel())) micount++;
194 }
195 if (micount>1) {
cf98c13f 196 track_multifound[itrack_multifound]=lab;
197 track_multifound_n[itrack_multifound]=micount;
198 itrack_multifound++;
199 }
200
7f6ddf58 201 Double_t xk=gt[ngood].x;
cc80f89e 202 track->PropagateTo(xk);
3c0f9266 203
cc80f89e 204 if (lab==tlab) hfound->Fill(ptg);
cf98c13f 205 else {
206 track_fake[itrack_fake++]=lab;
207 hfake->Fill(ptg);
cf98c13f 208 }
cc80f89e 209
be9c5115 210 Double_t par[5]; track->GetExternalParameters(xk,par);
211 Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
212 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
213 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
214 Float_t lam=TMath::ATan(par[3]);
215 Float_t pt_1=TMath::Abs(par[4]);
cc80f89e 216
217 if (TMath::Abs(gt[ngood].code)==11 && ptg>4.) {
be9c5115 218 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
cc80f89e 219 } else {
220 Float_t phig=TMath::ATan2(gt[ngood].py,gt[ngood].px);
cc80f89e 221 hp->Fill((phi - phig)*1000.);
3c0f9266 222
cc80f89e 223 Float_t lamg=TMath::ATan2(gt[ngood].pz,ptg);
cc80f89e 224 hl->Fill((lam - lamg)*1000.);
225
be9c5115 226 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
cc80f89e 227 }
73042f01 228
be9c5115 229 Float_t mom=1./(pt_1*TMath::Cos(lam));
73042f01 230 Float_t dedx=track->GetdEdx();
cc80f89e 231 hep->Fill(mom,dedx,1.);
232 if (TMath::Abs(gt[ngood].code)==211)
233 if (mom>0.4 && mom<0.5) {
234 he->Fill(dedx,1.);
235 }
73042f01 236
cc80f89e 237 }
238
cf98c13f 239 cout<<"\nList of Not found tracks :\n";
cf98c13f 240 for ( i = 0; i< itrack_notfound; i++){
241 cout<<track_notfound[i]<<"\t";
242 if ((i%5)==4) cout<<"\n";
243 }
244 cout<<"\nList of fake tracks :\n";
245 for ( i = 0; i< itrack_fake; i++){
246 cout<<track_fake[i]<<"\t";
247 if ((i%5)==4) cout<<"\n";
248 }
249 cout<<"\nList of multiple found tracks :\n";
250 for ( i=0; i<itrack_multifound; i++) {
7f6ddf58 251 cout<<"id. "<<track_multifound[i]
252 <<" found - "<<track_multifound_n[i]<<"times\n";
cf98c13f 253 }
8c555625 254
7f6ddf58 255 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
256 if (ng!=0) cerr<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
257 cout<<"Total number of found tracks ="<<nentr<<endl;
258 cout<<"Total number of \"good\" tracks ="
259 <<mingood<<" (selected for comparison: "<<ng<<')'<<endl<<endl;
260
261
8c555625 262 gStyle->SetOptStat(111110);
263 gStyle->SetOptFit(1);
264
265 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
266
267 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
268 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
269 hp->SetFillColor(4); hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd();
270
271 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
272 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
273 hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd();
274
275 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
276 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
277 hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd();
278
279 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
280 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
cc80f89e 281 hmpt->SetXTitle("(%)"); hmpt->Fit("gaus"); c1->cd();
8c555625 282
283 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
3c0f9266 284 p5->SetFillColor(41); p5->SetFrameFillColor(10);
8c555625 285 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
286 hg->Divide(hfound,hgood,1,1.,"b");
287 hf->Divide(hfake,hgood,1,1.,"b");
288 hg->SetMaximum(1.4);
289 hg->SetYTitle("Tracking efficiency");
290 hg->SetXTitle("Pt (GeV/c)");
291 hg->Draw();
292
7f6ddf58 293 TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
8c555625 294 line1->Draw("same");
7f6ddf58 295 TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
8c555625 296 line2->Draw("same");
297
298 hf->SetFillColor(1);
299 hf->SetFillStyle(3013);
300 hf->SetLineColor(2);
301 hf->SetLineWidth(2);
302 hf->Draw("histsame");
303 TText *text = new TText(0.461176,0.248448,"Fake tracks");
304 text->SetTextSize(0.05);
305 text->Draw();
306 text = new TText(0.453919,1.11408,"Good tracks");
307 text->SetTextSize(0.05);
308 text->Draw();
3c0f9266 309
310
311
312 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
313
314 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
315 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
316 he->SetFillColor(2); he->SetFillStyle(3005);
317 he->SetXTitle("Arbitrary Units");
318 he->Fit("gaus"); c2->cd();
319
320 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
321 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
322 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
323 hep->Draw(); c1->cd();
324
cf98c13f 325 gBenchmark->Stop("AliTPCComparison");
326 gBenchmark->Show("AliTPCComparison");
327
88cb7938 328 delete rl;
73042f01 329 return 0;
cc80f89e 330}
3c0f9266 331
cc80f89e 332
c630aafd 333Int_t
334good_tracks_tpc(GoodTrackTPC *gt, const Int_t max, const char* evfoldname) {
7f6ddf58 335 Int_t nt=0;
cc80f89e 336
c630aafd 337 AliRunLoader *rl = AliRunLoader::GetRunLoader(evfoldname);
338 if (rl == 0x0) {
88cb7938 339 ::Fatal("AliTPCComparison.C::good_tracks_tpc",
340 "Can not find Run Loader in Folder Named %s",
341 evfoldname);
c630aafd 342 }
343 AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
344 if (tpcl == 0x0) {
5102bab6 345 ::Fatal("AliTPCHits2Digits.C","Can not find TPCLoader !");
c630aafd 346 }
88cb7938 347
348 rl->LoadgAlice();
349
350 AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
73042f01 351 Int_t ver = TPC->IsVersion();
352 cerr<<"TPC version "<<ver<<" has been found !\n";
353
88cb7938 354 rl->CdGAFile();
355 AliTPCParamSR *digp=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
5102bab6 356 if (!digp) {
357 ::Fatal("AliTPCHits2Digits.C","TPC parameters have not been found !");
358 }
73042f01 359 TPC->SetParam(digp);
360
88cb7938 361 rl->LoadHeader();
362
363 Int_t np = rl->GetHeader()->GetNtrack();
364
024a7fe9 365
73042f01 366 Int_t nrow_up=digp->GetNRowUp();
367 Int_t nrows=digp->GetNRowLow()+nrow_up;
368 Int_t zero=digp->GetZeroSup();
7f6ddf58 369 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
cc80f89e 370 Int_t good_number=Int_t(0.4*nrows);
371
7f6ddf58 372 Int_t *good=new Int_t[np];
373 Int_t i;
374 for (i=0; i<np; i++) good[i]=0;
375
376
377 switch (ver) {
378 case 1:
379 {
88cb7938 380 tpcl->LoadRecPoints();
7f6ddf58 381 AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
382 ca.Setup(digp);
383 ca.SetClusterType("AliTPCcluster");
88cb7938 384 ca.ConnectTree(tpcl->TreeR());
7f6ddf58 385 Int_t nrows=Int_t(ca.GetTree()->GetEntries());
386 for (Int_t n=0; n<nrows; n++) {
387 AliSegmentID *s=ca.LoadEntry(n);
388 Int_t sec,row;
389 digp->AdjustSectorRow(s->GetID(),sec,row);
390 AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
391 Int_t ncl=clrow.GetArray()->GetEntriesFast();
392 while (ncl--) {
393 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
394 Int_t lab=c->GetLabel(0);
395 if (lab<0) continue; //noise cluster
396 lab=TMath::Abs(lab);
397
398 if (sec>=digp->GetNInnerSector())
399 if (row==nrow_up-1) good[lab]|=0x4000;
400 if (sec>=digp->GetNInnerSector())
401 if (row==nrow_up-1-gap) good[lab]|=0x1000;
402
403 if (sec>=digp->GetNInnerSector())
404 if (row==nrow_up-1-shift) good[lab]|=0x2000;
405 if (sec>=digp->GetNInnerSector())
406 if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
407
408 good[lab]++;
409 }
410 ca.ClearRow(sec,row);
411 }
412 delete pca;
88cb7938 413 tpcl->UnloadRecPoints();
7f6ddf58 414 }
afc42102 415 break;
7f6ddf58 416 case 2:
417 {
88cb7938 418 tpcl->LoadDigits();
419 TTree *TD=tpcl->TreeD();
420
7f6ddf58 421 AliSimDigits da, *digits=&da;
422 TD->GetBranch("Segment")->SetAddress(&digits);
423
424 Int_t *count = new Int_t[np];
425 Int_t i;
426 for (i=0; i<np; i++) count[i]=0;
427
428 Int_t sectors_by_rows=(Int_t)TD->GetEntries();
429 for (i=0; i<sectors_by_rows; i++) {
430 if (!TD->GetEvent(i)) continue;
431 Int_t sec,row;
432 digp->AdjustSectorRow(digits->GetID(),sec,row);
433 cerr<<sec<<' '<<row<<" \r";
434 digits->First();
435 do { //Many thanks to J.Chudoba who noticed this
436 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
437 Short_t dig = digits->GetDigit(it,ip);
438 Int_t idx0=digits->GetTrackID(it,ip,0);
439 Int_t idx1=digits->GetTrackID(it,ip,1);
440 Int_t idx2=digits->GetTrackID(it,ip,2);
8ded1b7a 441 if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
442 if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
443 if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
afc42102 444 } while (digits->Next());
7f6ddf58 445 for (Int_t j=0; j<np; j++) {
446 if (count[j]>1) {
447 if (sec>=digp->GetNInnerSector())
448 if (row==nrow_up-1 ) good[j]|=0x4000;
449 if (sec>=digp->GetNInnerSector())
450 if (row==nrow_up-1-gap) good[j]|=0x1000;
451
452 if (sec>=digp->GetNInnerSector())
453 if (row==nrow_up-1-shift) good[j]|=0x2000;
454 if (sec>=digp->GetNInnerSector())
455 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
456 good[j]++;
457 }
458 count[j]=0;
459 }
460 }
461 delete[] count;
88cb7938 462 tpcl->UnloadDigits();
afc42102 463 }
7f6ddf58 464 break;
465 default:
5102bab6 466 ::Fatal("AliTPCComparison.C","Invalid TPC version !");
cc80f89e 467 }
7f6ddf58 468
88cb7938 469 rl->LoadKinematics();
470 AliStack* stack = rl->Stack();
7f6ddf58 471 /** select tracks which are "good" enough **/
472 for (i=0; i<np; i++) {
473 if ((good[i]&0x5000) != 0x5000)
474 if ((good[i]&0x2800) != 0x2800) continue;
475 if ((good[i]&0x7FF ) < good_number) continue;
476
88cb7938 477 TParticle *p = (TParticle*)stack->Particle(i);
478 if (p == 0x0)
479 {
480 cerr<<"Can not get particle "<<i<<endl;
481 continue;
482 }
7f6ddf58 483 if (p->Pt()<0.100) continue;
484 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
485
486 Int_t j=p->GetFirstMother();
487 if (j>=0) {
88cb7938 488 TParticle *pp = (TParticle*)stack->Particle(j);
489 if (pp == 0x0)
490 {
491 cerr<<"Can not get particle "<<j<<endl;
492 continue;
493 }
024a7fe9 494 if (pp->GetFirstMother()>=0) continue;//only one decay is allowed
495 /* for cascade hyperons only
496 Int_t jj=pp->GetFirstMother();
497 if (jj>=0) {
88cb7938 498 TParticle *ppp = (TParticle*)stack->Particle(jj);
024a7fe9 499 if (ppp->GetFirstMother()>=0) continue;//two decays are allowed
500 }
501 */
7f6ddf58 502 }
503
504 gt[nt].lab=i;
505 gt[nt].code=p->GetPdgCode();
506 gt[nt].px=0.; gt[nt].py=0.; gt[nt].pz=0.;
527abbe9 507 gt[nt].x=0.; gt[nt].y=0.; gt[nt].z=0.;
7f6ddf58 508 nt++;
509 if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
510 cerr<<np-i<<" \r";
511 }
88cb7938 512 rl->UnloadKinematics();
7f6ddf58 513
514 /** check if there is also information at the entrance of the TPC **/
88cb7938 515
5102bab6 516 rl->LoadTrackRefs();
517 TTree *TR=rl->TreeTR();
518 TBranch *branch=TR->GetBranch("TPC");
519 if (branch==0) {
520 ::Fatal("AliTPCComparison.C::good_tracks_tpc","No track references !");
521 }
522 TClonesArray *references=new TClonesArray("AliTrackReference",1000);
523 branch->SetAddress(&references);
524 np=(Int_t)TR->GetEntries();
7f6ddf58 525 for (i=0; i<np; i++) {
5102bab6 526 references->Clear();
527 TR->GetEvent(i);
528 Int_t nref=references->GetEntriesFast();
529 if (nref==0) continue;
530 AliTrackReference *ref=(AliTrackReference*)references->UncheckedAt(0);
531
532 Int_t j, lab=ref->Label();
533 for (j=0; j<nt; j++) {if (gt[j].lab==lab) break;}
534 if (j==nt) continue;
535
536 // (px,py,pz) - in global coordinate system, (x,y,z) - in local !
537 gt[j].px=ref->Px(); gt[j].py=ref->Py(); gt[j].pz=ref->Pz();
538 gt[j].x = ref->LocalX();
539 gt[j].y = ref->LocalY();
540 gt[j].z = ref->Z();
7f6ddf58 541 }
542
5102bab6 543 delete references;
7f6ddf58 544 delete[] good;
88cb7938 545
546 tpcl->UnloadHits();
5102bab6 547 rl->UnloadTrackRefs();
88cb7938 548 rl->UnloadgAlice();
7f6ddf58 549
cf98c13f 550 gBenchmark->Stop("AliTPCComparison");
551 gBenchmark->Show("AliTPCComparison");
cc80f89e 552 return nt;
8c555625 553}