]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCComparison2.C
New macroc for handling signal and background events simultanously
[u/mrichter/AliRoot.git] / TPC / AliTPCComparison2.C
1 /****************************************************************************
2  *           Very important, delicate and rather obscure macro.             *
3  *                                                                          *
4  *               Creates list of "trackable" tracks,                        *
5  *             sorts tracks for matching with the ITS,                      *
6  *             calculates efficiency, resolutions etc.                      *
7  *                                                                          *
8  *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
9  * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de        *
10  ****************************************************************************/
11
12 #if !defined(__CINT__) || defined(__MAKECINT__)
13 #include<fstream.h>
14 #include<TROOT.h>
15 #include<TCanvas.h>
16 #include<TClassTable.h>
17 #include<TClonesArray.h>
18 #include<TFile.h>
19 #include<TFolder.h>
20 #include<TH1.h>
21 #include<TObject.h>
22 #include<TObjArray.h>
23 #include<TParticle.h>
24 #include<TTree.h>
25 #include <AliMagF.h>
26 #include <AliRun.h>
27 #include <AliTPCtracker.h>
28 #include <alles.h>
29 #endif
30
31 /*
32   multievent comparison of reconstructed tracks with "exact tracks" 
33   from silulation
34
35 */
36
37
38 struct GoodTrackTPC {
39   Int_t fEventN; //event number
40   Int_t lab;
41   Int_t code;
42   Float_t px,py,pz;
43   Float_t x,y,z;
44 };
45
46 enum tagprimary {kPrimaryCharged = 0x4000};
47
48 Int_t good_tracks(GoodTrackTPC *gt, Int_t max, Int_t firstev, Int_t eventn);
49
50 Int_t FindPrimaries(Int_t nparticles);
51
52 Int_t AliTPCComparison2(Int_t firstev=0, Int_t eventn=1) {
53
54   const Int_t MAX = 20000;
55   const Bool_t kOLD = kFALSE; // if tracking was done with a previous version
56   /***********************************************************************/
57
58   TFile *inkin = TFile::Open("rfio:galice.root");
59   //  if(gAlice)delete gAlice;   COMMENTED BECAUSE OF A BUG (IN COMPILED MODE)
60   gAlice = (AliRun*)inkin->Get("gAlice");
61   cout<<"AliRun object found on file "<<gAlice<<endl;
62   Float_t fifac=gAlice->Field()->Factor();
63   AliKalmanTrack::SetConvConst(100/0.299792458/0.2/fifac);
64   cout<<fifac*0.2<<" T Magnetic field is used\n";
65   inkin->Close();
66   /*
67     delete gAlice;  COMMENTED BECAUSE OF A BUG IN COMPILED MODE
68     gAlice = 0;
69   */
70   /***********************************************************************/
71   cerr<<"Doing comparison...\n";
72   Int_t i;
73   gBenchmark->Start("AliTPCComparison2");
74
75   TFile *cf=0; 
76   AliTPCParam *digp=0;
77   if(kOLD){
78     cf=TFile::Open("AliTPCclusters.root");
79     if (!cf->IsOpen()) {cerr<<"Can't open AliTPCclusters.root !\n"; return 1;}
80     digp= (AliTPCParam*)cf->Get("75x40_100x60");
81     if (!digp) { cerr<<"TPC parameters have not been found !\n"; return 2; }
82   }
83   ///////////
84   AliTPCtracker *tracker =0; 
85   TObjArray tarray(MAX);
86   AliTPCtrack *iotrack=0;
87   Int_t nentr= 0;
88   Int_t eventptr[1000];
89   TFile *tf=TFile::Open("AliTPCtracks.root");
90   TTree *tracktree=0;
91
92   eventptr[firstev]=0;
93   eventptr[firstev+1]=0;
94   for (Int_t event=firstev;event<eventn; event++){
95     cout<<"================================================"<<endl;
96     cout<<"Processing event "<<event<<endl;
97     cout<<"================================================"<<endl;
98     if(kOLD){
99       cf->cd();
100       tracker = new AliTPCtracker(digp,event);
101       tracker->LoadInnerSectors();
102       tracker->LoadOuterSectors();
103     }
104     
105     char   tname[100];
106     if (eventn==-1) {
107       sprintf(tname,"TreeT_TPC");
108     }
109     else {
110       sprintf(tname,"TreeT_TPC_%d",event);
111     }
112   
113     // Load tracks
114     if (!tf->IsOpen()) {cerr<<"Can't open AliTPCtracks.root !\n"; return 3;}
115
116     tracktree=(TTree*)tf->Get(tname);
117     if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n"; return 4;}
118     TBranch *tbranch=tracktree->GetBranch("tracks");
119     Int_t nentr0=(Int_t)tracktree->GetEntries();
120     nentr+=nentr0;
121     for (i=0; i<nentr0; i++) {
122       iotrack=new AliTPCtrack;
123       tbranch->SetAddress(&iotrack);
124       tracktree->GetEvent(i);
125       if(kOLD)tracker->CookLabel(iotrack,0.1);
126       tarray.AddLast(iotrack);
127     }   
128     eventptr[event+1] = nentr;  //store end of the event
129     delete tracker;
130     delete tracktree;
131   }
132   tf->Close();
133   if(kOLD)cf->Close();
134
135   GoodTrackTPC gt[MAX];
136
137   Int_t ngood=0;
138   ifstream in("good_tracks_tpc");
139   if (in) {
140     cerr<<"Reading good tracks...\n";
141     while (in>>gt[ngood].fEventN>>gt[ngood].lab>>gt[ngood].code>>
142            gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
143            gt[ngood].x >>gt[ngood].y >>gt[ngood].z) {
144       ngood++;
145       cerr<<ngood<<"\r";
146       if (ngood==MAX) {
147         cerr<<"Too many good tracks !\n";
148         break;
149       }
150     }
151     cout<<endl;
152     if (!in.eof()) cerr<<"Read error (good_tracks_tpc) !\n";
153   } else {
154     cerr<<"Marking good tracks (this will take a while)...\n";
155     ngood=good_tracks(gt,45000,firstev,eventn); 
156     ofstream out("good_tracks_tpc");
157     ofstream out2("good_tracks_tpc_par");
158
159     if (out) {
160       cout<<"File good_tracks_tpc opened\n";
161       for (Int_t ngd=0; ngd<ngood; ngd++) {
162         Float_t pt =  TMath::Sqrt(gt[ngd].px*gt[ngd].px+gt[ngd].py*gt[ngd].py);
163         Float_t angle = gt[ngd].pz/pt;
164         out<<gt[ngd].fEventN<<' '<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
165           gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
166           gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
167         out2<<gt[ngd].fEventN<<"\t"<<gt[ngd].lab<<"\t"<<gt[ngd].code<<"\t"<<
168           pt<<"\t"<<angle<<"\t"<<endl;
169       }
170     } else cerr<<"Can not open file (good_tracks_tpc) !\n";
171     out.close();
172     out2.close();
173   }
174   cerr<<"Number of good tracks : "<<ngood<<endl;
175   if(ngood==0)return 5;
176   TH1F *hp=new TH1F("hp","PHI resolution",50,-200.,200.); hp->SetFillColor(4);
177   TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-200,200);hl->SetFillColor(4);
178   TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
179   hpt->SetFillColor(2); 
180   TH1F *hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60); 
181   hmpt->SetFillColor(6);
182
183   AliTPCtrack *trk=(AliTPCtrack*)tarray.UncheckedAt(0);
184
185   Double_t pmin=0.1*(100/0.299792458/0.2/trk->GetConvConst());
186   Double_t pmax=6.0+pmin;
187
188   TH1F *hgood=new TH1F("hgood","Good tracks",30,pmin,pmax);    
189   TH1F *hfound=new TH1F("hfound","Found tracks",30,pmin,pmax);
190   TH1F *hfake=new TH1F("hfake","Fake tracks",30,pmin,pmax);
191   TH1F *hg=new TH1F("hg","",30,pmin,pmax); //efficiency for good tracks
192   hg->SetLineColor(4); hg->SetLineWidth(2);
193   TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,pmin,pmax);
194   hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
195
196   TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
197   TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,250.);
198   hep->SetMarkerStyle(8);
199   hep->SetMarkerSize(0.4);
200
201   Int_t mingood=ngood;
202   Int_t track_notfound[MAX], itrack_notfound=0;
203   Int_t track_fake[MAX], itrack_fake=0;
204   Int_t track_multifound[MAX], track_multifound_n[MAX], itrack_multifound=0;
205   while (ngood--) {
206     Int_t lab=gt[ngood].lab,tlab=-1;
207     Float_t ptg=
208       TMath::Sqrt(gt[ngood].px*gt[ngood].px + gt[ngood].py*gt[ngood].py);
209
210     if (ptg<pmin) continue;  //  Should be 1e-33 instead pmin?
211
212     hgood->Fill(ptg);
213     Int_t ievent = gt[ngood].fEventN;
214
215     AliTPCtrack *track=0;
216     for (i=eventptr[ievent]; i<eventptr[ievent+1]; i++) {
217
218       track=(AliTPCtrack*)tarray.UncheckedAt(i);
219       tlab=track->GetLabel();
220       if (lab==TMath::Abs(tlab)) break;
221     }
222     if (i==eventptr[ievent+1]) {
223       track_notfound[itrack_notfound++]=lab;
224       continue;
225     }
226       
227     Int_t micount=0;
228     Int_t mi;
229     AliTPCtrack * mitrack;
230     for (mi=eventptr[ievent]; mi<eventptr[ievent+1]; mi++) {
231
232       mitrack=(AliTPCtrack*)tarray.UncheckedAt(mi);          
233       if (lab==TMath::Abs(mitrack->GetLabel())) micount++;
234     }
235     if (micount>1) {
236       track_multifound[itrack_multifound]=lab;
237       track_multifound_n[itrack_multifound]=micount;
238       itrack_multifound++;
239     }
240       
241     //
242     Double_t xk=gt[ngood].x;
243     printf("Track =%p\n",track);
244     track->PropagateTo(xk);
245
246     if (lab==tlab) hfound->Fill(ptg);
247     else { 
248       track_fake[itrack_fake++]=lab;
249       hfake->Fill(ptg); 
250     }
251     Double_t par[5]; track->GetExternalParameters(xk,par);
252     Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
253     if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
254     if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
255     Float_t lam=TMath::ATan(par[3]); 
256     Float_t pt_1=TMath::Abs(par[4]);
257
258     if (TMath::Abs(gt[ngood].code)==11 && ptg>4.) {
259       hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
260     } else {
261       Float_t phig=TMath::ATan2(gt[ngood].py,gt[ngood].px);
262       hp->Fill((phi - phig)*1000.);
263
264       Float_t lamg=TMath::ATan2(gt[ngood].pz,ptg);
265       hl->Fill((lam - lamg)*1000.);
266
267       hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
268     }
269
270     Float_t mom=1./(pt_1*TMath::Cos(lam));
271     Float_t dedx=track->GetdEdx();
272     hep->Fill(mom,dedx,1.);
273     if (TMath::Abs(gt[ngood].code)==211)
274       if (mom>0.4 && mom<0.5) {
275         he->Fill(dedx,1.);
276       }
277
278   }
279
280   cout<<"\nList of Not found tracks :\n";
281   for ( i = 0; i< itrack_notfound; i++){
282     cout<<track_notfound[i]<<"\t";
283     if ((i%5)==4) cout<<"\n";
284   }
285   cout<<"\nList of fake  tracks :\n";
286   for ( i = 0; i< itrack_fake; i++){
287     cout<<track_fake[i]<<"\t";
288     if ((i%5)==4) cout<<"\n";
289   }
290   cout<<"\nList of multiple found tracks :\n";
291   for ( i=0; i<itrack_multifound; i++) {
292     cout<<"id.   "<<track_multifound[i]
293         <<"     found - "<<track_multifound_n[i]<<"times\n";
294   }
295   
296   Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
297   if (ng!=0) cerr<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
298   cout<<"Total number of found tracks ="<<nentr<<endl;
299   cout<<"Total number of \"good\" tracks ="
300       <<mingood<<"   (selected for comparison: "<<ng<<')'<<endl<<endl;
301
302    
303   gStyle->SetOptStat(111110);
304   gStyle->SetOptFit(1);
305
306   TCanvas *c1=new TCanvas("c1","",0,0,700,850);
307
308   TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
309   p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); 
310   hp->SetFillColor(4);  hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd();
311
312   TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); 
313   p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
314   hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd();
315
316   TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
317   p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); 
318   hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd();
319
320   TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
321   p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
322   hmpt->SetXTitle("(%)"); 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
334   TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
335   line1->Draw("same");
336   TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
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   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   gBenchmark->Stop("AliTPCComparison2");
367   gBenchmark->Show("AliTPCComparison2");
368
369
370   return 0;
371 }
372
373
374 Int_t good_tracks(GoodTrackTPC *gt, Int_t max, Int_t firstev, Int_t eventn) {
375   //eventn  - number of events in file
376
377   TFile *file=TFile::Open("rfio:galice.root");
378   if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
379   //  delete gAlice; gAlice = 0;
380   if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
381     cerr<<"gAlice have not been found on galice.root !\n";
382     exit(5);
383   }
384    
385
386   AliTPC *TPC=(AliTPC*)gAlice->GetDetector("TPC");
387   Int_t ver = TPC->IsVersion(); 
388   cerr<<"TPC version "<<ver<<" has been found !\n";
389
390   AliTPCParam *digp=(AliTPCParam*)file->Get("75x40_100x60");
391   if (!digp) { cerr<<"TPC parameters have not been found !\n"; exit(6); }
392   TPC->SetParam(digp);
393
394   Int_t nrow_up=digp->GetNRowUp();
395   Int_t nrows=digp->GetNRowLow()+nrow_up;
396   Int_t zero=digp->GetZeroSup();
397   Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
398   Int_t good_number=Int_t(0.4*nrows);
399
400   Int_t nt=0;  //reset counter
401   char treeName[100];  //declare event identifier
402        
403   for (Int_t event=firstev;event<eventn;event++){
404     Int_t np=gAlice->GetEvent(event);
405     Int_t nopr = FindPrimaries(np);
406     cout<<"function good_tracks -- event "<<event<<", Charged primaries:"<<nopr<<"\n";
407
408     Int_t *good=new Int_t[np];
409     for (Int_t ii=0; ii<np; ii++) good[ii]=0;
410      
411      
412     Int_t sectors_by_rows=0;
413     TTree *TD=0;
414     AliSimDigits da, *digits=&da;
415     Int_t *count=0;
416     switch (ver) {
417     case 1:
418       {
419         TFile *cf=TFile::Open("AliTPCclusters.root");
420         if (!cf->IsOpen()){cerr<<"Can't open AliTPCclusters.root !\n";exit(5);}
421         AliTPCClustersArray *ca=new AliTPCClustersArray;
422         ca->Setup(digp);
423         ca->SetClusterType("AliTPCcluster");
424         ca->ConnectTree("TreeC_TPC_0");
425         Int_t nrows=Int_t(ca->GetTree()->GetEntries());
426         for (Int_t n=0; n<nrows; n++) {
427           AliSegmentID *s=ca->LoadEntry(n);
428           Int_t sec,row;
429           digp->AdjustSectorRow(s->GetID(),sec,row);
430           AliTPCClustersRow &clrow = *ca->GetRow(sec,row);
431           Int_t ncl=clrow.GetArray()->GetEntriesFast();
432           while (ncl--) {
433             AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
434             Int_t lab=c->GetLabel(0);
435             if (lab<0) continue; //noise cluster
436             lab=TMath::Abs(lab);
437             if (sec>=digp->GetNInnerSector())
438               if (row==nrow_up-1    ) good[lab]|=0x1000;
439             if (sec>=digp->GetNInnerSector())
440               if (row==nrow_up-1-gap) good[lab]|=0x800;
441             good[lab]++;
442           }
443           ca->ClearRow(sec,row);
444         }
445         delete ca;
446         cf->Close();
447       }
448       break;
449     case 2:
450       {
451         sprintf(treeName,"TreeD_75x40_100x60_%d",event);  
452         TD=(TTree*)gDirectory->Get(treeName);
453         TD->GetBranch("Segment")->SetAddress(&digits);
454         count = new Int_t[np];
455         Int_t i;
456         for (i=0; i<np; i++) count[i]=0;
457         sectors_by_rows=(Int_t)TD->GetEntries();
458         for (i=0; i<sectors_by_rows; i++) {
459           if (!TD->GetEvent(i)) continue;
460           Int_t sec,row;
461           digp->AdjustSectorRow(digits->GetID(),sec,row);
462           digits->First();
463           do { //Many thanks to J.Chudoba who noticed this
464             Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
465             Short_t dig = digits->GetDigit(it,ip);
466             Int_t idx0=digits->GetTrackID(it,ip,0); 
467             Int_t idx1=digits->GetTrackID(it,ip,1);
468             Int_t idx2=digits->GetTrackID(it,ip,2);
469             if (idx0>=0 && dig>=zero && idx0<np ) count[idx0]+=1;   //background events - higher track ID's
470             if (idx1>=0 && dig>=zero && idx1<np ) count[idx1]+=1;
471             if (idx2>=0 && dig>=zero && idx2<np ) count[idx2]+=1;
472           } while (digits->Next());
473           for (Int_t j=0; j<np; j++) {
474             if (count[j]>1) {
475               if (sec>=digp->GetNInnerSector())
476                 if (row==nrow_up-1    ) good[j]|=0x4000;
477               if (sec>=digp->GetNInnerSector())
478                 if (row==nrow_up-1-gap) good[j]|=0x1000;
479
480               if (sec>=digp->GetNInnerSector())
481                 if (row==nrow_up-1-shift) good[j]|=0x2000;
482               if (sec>=digp->GetNInnerSector())
483                 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
484               good[j]++;
485             }
486             count[j]=0;
487           }
488
489         }
490         delete[] count;
491         delete TD; //Thanks to Mariana Bondila
492       }
493       break;
494     default:
495       cerr<<"Invalid TPC version !\n";
496       file->Close();
497       exit(7);
498     }
499   
500     /** select tracks which are "good" enough **/
501     for (Int_t i=0; i<np; i++) {
502       if ((good[i]&0x5000) != 0x5000)
503         if ((good[i]&0x2800) != 0x2800) continue;
504       if ((good[i]&0x7FF ) < good_number) continue;
505
506       TParticle *p = (TParticle*)gAlice->Particle(i);
507
508       if (p->Pt()<0.100) continue;
509       if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
510       //  THIS IS GOOD IF YOU WANT SECONDARIES
511         Int_t j=p->GetFirstMother();
512         if (j>=0) {
513         TParticle *pp = (TParticle*)gAlice->Particle(j);
514         //if (!(pp->TestBit(kPrimaryCharged))) continue; //only one decay is allowed
515         }
516       
517         //if(!(p->TestBit(kPrimaryCharged)))continue; // only primaries
518   
519       gt[nt].fEventN=event;
520       gt[nt].lab=i;
521       gt[nt].code=p->GetPdgCode();
522       gt[nt].px=0.; gt[nt].py=0.; gt[nt].pz=0.;
523       gt[nt].x=0.; gt[nt].z=0.; gt[nt].z=0.;
524       nt++;
525       if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
526       cerr<<np-i<<"        \r";    
527     }
528
529     /** check if there is also information at the entrance of the TPC **/
530     TTree *TH=gAlice->TreeH(); np=(Int_t)TH->GetEntries();
531     cout<<endl;
532     for (Int_t i=0; i<np; i++) {
533       TPC->ResetHits();
534       TH->GetEvent(i);
535       AliTPChit *phit = (AliTPChit*)TPC->FirstHit(-1);
536       for ( ; phit; phit=(AliTPChit*)TPC->NextHit() ) {
537         if (phit->fQ !=0. ) continue;
538
539         Double_t px=phit->X(), py=phit->Y(), pz=phit->Z();
540
541         if ((phit=(AliTPChit*)TPC->NextHit())==0) break;
542         if (phit->fQ != 0.) continue;
543
544         Double_t x=phit->X(), y=phit->Y(), z=phit->Z();
545         if (TMath::Sqrt(x*x+y*y)>90.) continue;
546
547         Int_t j, lab=phit->Track();
548         for (j=0; j<nt; j++) {if (gt[j].fEventN==event && gt[j].lab==lab) break;}
549         if (j==nt) continue;         
550
551         // (px,py,pz) - in global coordinate system, (x,y,z) - in local !
552         gt[j].px=px; gt[j].py=py; gt[j].pz=pz;
553         Float_t cs,sn; digp->AdjustCosSin(phit->fSector,cs,sn);
554         gt[j].x = x*cs + y*sn;
555         gt[j].y =-x*sn + y*cs;
556         gt[j].z = z;
557       }
558       cerr<<np-i<<"        \r";
559     }
560     cout<<endl;
561     delete[] good;
562   }             // ///         loop on events
563   delete gAlice; gAlice=0;
564
565   file->Close();
566   gBenchmark->Stop("AliTPCComparison2");
567   gBenchmark->Show("AliTPCComparison2");   
568   return nt;
569
570 }
571
572 Int_t FindPrimaries(Int_t nparticles){
573
574   // cuts:
575   Double_t vertcut = 0.001;
576   Double_t decacut = 3.;
577   Double_t timecut = 0.;
578   Int_t nprch1=0;
579   TParticle * part = gAlice->Particle(0);
580   Double_t xori = part->Vx();
581   Double_t yori = part->Vy();
582   Double_t zori = part->Vz();
583   for(Int_t iprim = 0; iprim<nparticles; iprim++){   //loop on  tracks
584
585     part = gAlice->Particle(iprim);
586     char *xxx=strstr(part->GetName(),"XXX");
587     if(xxx)continue;
588
589     TParticlePDG *ppdg = part->GetPDG();
590     if(TMath::Abs(ppdg->Charge())!=3)continue;  // only charged (no quarks)
591
592     Double_t dist=TMath::Sqrt((part->Vx()-xori)*(part->Vx()-xori)+(part->Vy()-yori)*(part->Vy()-yori)+(part->Vz()-zori)*(part->Vz()-zori));
593     if(dist>vertcut)continue;  // cut on the vertex
594
595     if(part->T()>timecut)continue;
596
597     Double_t ptot=TMath::Sqrt(part->Px()*part->Px()+part->Py()*part->Py()+part->Pz()*part->Pz());
598     if(ptot==(TMath::Abs(part->Pz())))continue; // no beam particles
599
600     Bool_t prmch = kTRUE;   // candidate primary track
601     Int_t fidau=part->GetFirstDaughter();  // cut on daughters
602     Int_t lasdau=0;
603     Int_t ndau=0;
604     if(fidau>=0){
605       lasdau=part->GetLastDaughter();
606       ndau=lasdau-fidau+1;
607     }
608     if(ndau>0){
609       for(Int_t j=fidau;j<=lasdau;j++){
610         TParticle *dau=gAlice->Particle(j);
611         Double_t distd=TMath::Sqrt((dau->Vx()-xori)*(dau->Vx()-xori)+(dau->Vy()-yori)*(dau->Vy()-yori)+(dau->Vz()-zori)*(dau->Vz()-zori));
612         if(distd<decacut)prmch=kFALSE;  // eliminate if the decay is near the vertex
613       }
614     }
615
616     if(prmch){
617       nprch1++;
618       part->SetBit(kPrimaryCharged);
619     }
620   }
621
622   return nprch1;
623 }
624
625
626
627
628
629