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