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