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