]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/src/AliL3Evaluate.cxx
minor changes, bugfix for FastSim
[u/mrichter/AliRoot.git] / HLT / src / AliL3Evaluate.cxx
1 //Author:        Anders Strand Vestbo
2 //Last Modified: 5.01.2001
3
4 #include <stdio.h>
5 #include <TFile.h>
6 #include <TH1.h>
7 #include <TParticle.h>
8 #include <TTree.h>
9 #include <TClonesArray.h>
10
11 #include "AliRun.h"
12 #include "AliSimDigits.h"
13 #include "AliTPC.h"
14 #include "AliTPCClustersArray.h"
15 #include "AliTPCClustersRow.h"
16 #include "AliTPCcluster.h"
17 #include "AliTPCParam.h"
18
19 #include "AliL3Transform.h"
20 #include "AliL3SpacePointData.h"
21 #include "AliL3Track.h"
22 #include "AliL3FileHandler.h"
23 #include "AliL3TrackArray.h"
24 #include "AliL3Evaluate.h"
25 #include "AliL3Logging.h"
26
27 //AliL3Evaluate
28 //Class for tracking evaluation.
29
30 ClassImp(AliL3Evaluate)
31
32 AliL3Evaluate::AliL3Evaluate()
33 {
34   fMCFile = NULL;
35   fTracks = NULL;
36   fMCclusterfile = NULL;
37 }
38
39 AliL3Evaluate::AliL3Evaluate(Char_t *mcfile,Int_t *slice)
40 {
41   //Normal constructor. Input are the rootfile containing the 
42   //original MC information of the simulated event. 
43
44   fMCFile = new TFile(mcfile,"READ");
45   if(!fMCFile->IsOpen())
46     {
47       LOG(AliL3Log::kError,"AliL3Evaluation::AliL3Evaluation","File Open")
48         <<"Inputfile "<<mcfile<<" does not exist"<<ENDLOG;
49       return;
50     }
51   
52   fParam = (AliTPCParam*)fMCFile->Get("75x40_100x60");
53   fTransform = new AliL3Transform();
54   
55   fMinSlice = slice[0];
56   fMaxSlice = slice[1];
57
58 }
59
60 AliL3Evaluate::~AliL3Evaluate()
61 {
62   if(fDigitsTree) 
63     fDigitsTree->Delete();
64   
65
66   /*if(fMCFile) 
67     {
68       fMCFile->Close();
69       delete fMCFile;
70     }
71   */
72   if(fTransform)
73     delete fTransform;
74   if(fTracks)
75     delete fTracks;
76   if(fPtRes)
77     delete fPtRes;
78   if(fNGoodTracksPt)
79     delete fNGoodTracksPt;
80   if(fNFoundTracksPt)
81     delete fNFoundTracksPt;
82   if(fNFakeTracksPt)
83     delete fNFakeTracksPt;
84   if(fTrackEffPt)
85     delete fTrackEffPt;
86   if(fFakeTrackEffPt)
87     delete fFakeTrackEffPt;
88   if(fNGoodTracksEta)
89     delete fNGoodTracksEta;
90   if(fNFoundTracksEta)
91     delete fNFoundTracksEta;
92   if(fNFakeTracksEta)
93     delete fNFakeTracksEta;
94   if(fTrackEffEta)
95     delete fTrackEffEta;
96   if(fFakeTrackEffEta)
97     delete fFakeTrackEffEta;
98 }
99
100 void AliL3Evaluate::Setup(Char_t *trackfile,Char_t *path)
101 {
102   //Read in the hit and track information from produced files.
103   
104   Char_t fname[256];
105   AliL3FileHandler *clusterfile[36][5];
106   for(Int_t s=fMinSlice; s<=fMaxSlice; s++)
107     {
108       for(Int_t p=0; p<5; p++)
109         {
110           fClusters[s][p] = 0;
111           clusterfile[s][p] = new AliL3FileHandler();
112           sprintf(fname,"%s/points_%d_%d.raw",path,s,p);
113           if(!clusterfile[s][p]->SetBinaryInput(fname))
114             {
115               LOG(AliL3Log::kError,"AliL3Evaluation::Setup","File Open")
116                 <<"Inputfile "<<fname<<" does not exist"<<ENDLOG; 
117               delete clusterfile[s][p];
118               clusterfile[s][p] = 0; 
119               continue;
120             }
121           fClusters[s][p] = (AliL3SpacePointData*)clusterfile[s][p]->Allocate();
122           clusterfile[s][p]->Binary2Memory(fNcl[s][p],fClusters[s][p]);
123           clusterfile[s][p]->CloseBinaryInput();
124         }
125     }
126   
127   
128   AliL3FileHandler *tfile = new AliL3FileHandler();
129   if(!tfile->SetBinaryInput(trackfile)){
130     LOG(AliL3Log::kError,"AliL3Evaluation::Setup","File Open")
131     <<"Inputfile "<<trackfile<<" does not exist"<<ENDLOG; 
132     return;
133   }
134   fTracks = new AliL3TrackArray();
135   tfile->Binary2TrackArray(fTracks);
136   tfile->CloseBinaryInput();
137   delete tfile;
138 }
139
140 void AliL3Evaluate::SetupSlow(Char_t *trackfile,Char_t *path)
141 {
142   //Setup for using the slow simulator.
143   
144   fIsSlow = true;
145   Setup(trackfile,path);
146   
147   if(!SetDigitsTree())
148     LOG(AliL3Log::kError,"AliL3Evaluation::SetupSlow","Digits Tree")
149     <<"Error setting up digits tree"<<ENDLOG;
150   if(!SetMCParticleArray())
151     LOG(AliL3Log::kError,"AliL3Evaluation::Setup","Particle array")
152     <<"Error setting up particle array"<<ENDLOG;
153 }
154
155 void AliL3Evaluate::SetupFast(Char_t *trackfile,Char_t *mcClusterfile,Char_t *path)
156 {
157   //Setup for using the fast simulator.
158
159   fIsSlow = false;
160   
161   fMCclusterfile = new TFile(mcClusterfile);
162   if(!fMCclusterfile->IsOpen())
163     LOG(AliL3Log::kError,"AliL3Evaluation::SetupFast","File Open")
164       <<"Inputfile "<<mcClusterfile<<" does not exist"<<ENDLOG; 
165
166   Setup(trackfile,path);
167
168   if(!SetMCParticleArray())
169     LOG(AliL3Log::kError,"AliL3Evaluation::SetupFast","Particle array")
170       <<"Error setting up particle array"<<ENDLOG;
171 }
172
173 void AliL3Evaluate::CreateHistos(Int_t nbin,Int_t xlow,Int_t xup)
174 {
175   //Create the histograms 
176   
177   fPtRes = new TH1F("fPtRes","Relative Pt resolution",30,-10.,10.); 
178   fNGoodTracksPt = new TH1F("fNGoodTracksPt","Good tracks vs pt",nbin,xlow,xup);    
179   fNFoundTracksPt = new TH1F("fNFoundTracksPt","Found tracks vs pt",nbin,xlow,xup);
180   fNFakeTracksPt = new TH1F("fNFakeTracksPt","Fake tracks vs pt",nbin,xlow,xup);
181   fTrackEffPt = new TH1F("fTrackEffPt","Tracking efficiency vs pt",nbin,xlow,xup);
182   fFakeTrackEffPt = new TH1F("fFakeTrackEffPt","Efficiency for fake tracks vs pt",nbin,xlow,xup);
183   
184   fNGoodTracksEta = new TH1F("fNGoodTracksEta","Good tracks vs eta",20,-50,50);
185   fNFoundTracksEta = new TH1F("fNFoundTracksEta","Found tracks vs eta",20,-50,50);
186   fNFakeTracksEta = new TH1F("fNFakeTracksEta","Fake tracks vs eta",20,-50,50);
187   fTrackEffEta = new TH1F("fTrackEffEta","Tracking efficienct vs eta",20,-50,50);
188   fFakeTrackEffEta = new TH1F("fFakeTrackEffEta","Efficiency for fake tracks vs eta",20,-50,50);
189
190 }
191
192 Bool_t AliL3Evaluate::SetDigitsTree()
193 {
194
195   fMCFile->cd();
196   fDigitsTree = (TTree*)fMCFile->Get("TreeD_75x40_100x60");
197   if(!fDigitsTree) return false;
198   fDigitsTree->GetBranch("Segment")->SetAddress(&fDigits);
199   for(Int_t i=0; i<fDigitsTree->GetEntries(); i++)
200     {
201       if(!fDigitsTree->GetEvent(i)) continue;
202       Int_t se,ro,slice,slicerow;
203       fParam->AdjustSectorRow(fDigits->GetID(),se,ro);
204       fTransform->Sector2Slice(slice,slicerow,se,ro);
205       fRowid[slice][slicerow] = i;
206     }
207   
208   return true;
209 }
210
211 Bool_t AliL3Evaluate::SetMCParticleArray()
212 {
213   fMCFile->cd();
214   AliRun *gAlice = (AliRun*)fMCFile->Get("gAlice");
215   if (!gAlice) 
216     {
217       LOG(AliL3Log::kError,"AliL3Evaluate::SetParticleArray","gAlice")
218         <<"AliRun object non existing on file"<<ENDLOG;
219       return false;
220     }
221   
222   gAlice->GetEvent(0);
223   fParticles=gAlice->Particles(); 
224   return true;
225 }
226
227
228 TObjArray *AliL3Evaluate::DefineGoodTracks(Int_t slice,Int_t *padrow,Int_t good_number,Int_t *particle_id)
229 {
230   //Loop over MC particles, and mark the good ones
231   //(which the tracker should find...)
232
233   Int_t np=fParticles->GetEntriesFast();
234   AliTPC *TPC = (AliTPC*)gAlice->GetDetector("TPC");
235   TPC->SetParam(fParam);
236   Int_t ver = TPC->IsVersion();
237   LOG(AliL3Log::kInformational,"AliL3Evaluate::DefineGoodTracks","TPC version")
238     <<"TPC version "<<ver<<" found on file"<<ENDLOG;
239   
240   //Int_t nrow_up=TPC->GetParam()->GetNRowUp();
241   //Int_t nrows=TPC->GetParam()->GetNRowLow()+nrow_up;
242   Int_t zero=TPC->GetParam()->GetZeroSup();
243   
244   //Int_t number_of_rows = padrow[1] - padrow[0] + 1;
245   
246   Int_t *good = new Int_t[np];
247   for(Int_t ii=0; ii<np; ii++)
248     good[ii] = 0;
249
250   if(ver==1)
251     {
252       if(fIsSlow)
253         LOG(AliL3Log::kError,"AliL3Evaluate::DefineGoodTracks","TPC version")
254           <<"TPC version "<<ver<<" does not match."<<ENDLOG;
255       fMCclusterfile->cd();
256       AliTPCClustersArray carray;
257       carray.Setup(fParam);
258       carray.SetClusterType("AliTPCcluster");
259       Bool_t clusterok = carray.ConnectTree("Segment Tree");
260       if(!clusterok) 
261         LOG(AliL3Log::kError,"AliL3Evaluate::DefineGoodTracks","Cluster Array")
262           <<"Error loading clusters from rootfile"<<ENDLOG;
263             
264       for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++)  
265         {
266           Int_t sec,row,sl,lr;
267           AliSegmentID *s = carray.LoadEntry(i);
268           fParam->AdjustSectorRow(s->GetID(),sec,row);
269           fTransform->Sector2Slice(sl,lr,sec,row);
270           
271           if(sl != slice) {carray.ClearRow(sec,row); continue;}
272           if(lr < padrow[0]) {carray.ClearRow(sec,row); continue;}
273           if(lr > padrow[1]) {carray.ClearRow(sec,row); continue;}
274           AliTPCClustersRow *cRow = carray.GetRow(sec,row);
275           for(Int_t j=0; j<cRow->GetArray()->GetEntriesFast(); j++)
276             {
277               AliTPCcluster *cluster=(AliTPCcluster*)(*cRow)[j];
278               Int_t lab=cluster->GetLabel(0);
279               if(lab<0) continue;
280               lab=TMath::Abs(lab);
281               good[lab]++;
282             }
283           if(carray.GetRow(sec,row)) 
284             carray.ClearRow(sec,row);
285         }
286     }
287   else if(ver==2)
288     {
289       if(!fIsSlow)
290         LOG(AliL3Log::kError,"AliL3Evaluate::DefineGoodTracks","TPC version")
291           <<"TPC version "<<ver<<" does not match."<<ENDLOG;
292       Int_t *count = new Int_t[np]; //np number of particles.
293       Int_t i;
294       for (i=0; i<np; i++) count[i]=0;
295       for (i=padrow[0]; i<=padrow[1]; i++) {
296         Int_t index = fRowid[slice][i];
297         if (!fDigitsTree->GetEvent(index)) continue;
298         Int_t sec,row;
299         fParam->AdjustSectorRow(fDigits->GetID(),sec,row);
300         fDigits->First();
301         while (fDigits->Next()) {
302           Int_t it=fDigits->CurrentRow(), ip=fDigits->CurrentColumn();
303           Short_t dig = fDigits->GetDigit(it,ip);
304           Int_t idx0=fDigits->GetTrackID(it,ip,0); 
305           Int_t idx1=fDigits->GetTrackID(it,ip,1);
306           Int_t idx2=fDigits->GetTrackID(it,ip,2);
307           if (idx0>=0 && dig>=zero) count[idx0]+=1;
308           if (idx1>=0 && dig>=zero) count[idx1]+=1;
309           if (idx2>=0 && dig>=zero) count[idx2]+=1;
310         }
311         for (Int_t j=0; j<np; j++) {
312           if (count[j]>1) {//at least two digits at this padrow 
313             good[j]++;
314           }
315           count[j]=0;
316         }
317       }
318       delete[] count;
319     }
320   else 
321     {
322       LOG(AliL3Log::kError,"AliL3Evaluation::FillEffHistos","TPC version")
323         <<"No valid TPC version found"<<ENDLOG;
324       return 0;
325     }
326   
327   Float_t torad=TMath::Pi()/180;
328   
329   Float_t phi_min = slice*20 - 10;
330   Float_t phi_max = slice*20 + 10;
331   TObjArray *good_part = new TObjArray();
332   
333   for(Int_t i=0; i<fParticles->GetEntriesFast(); i++)
334    {
335      TParticle *p = (TParticle*)fParticles->UncheckedAt(i);
336      if(p->GetFirstMother()>0) continue; //secondary particle
337      if(good[i] < good_number) {continue;}
338      
339      Double_t ptg=p->Pt(),pxg=p->Px(),pyg=p->Py(),pzg=p->Pz();
340      Double_t phi_part = TMath::ATan2(pyg,pxg);
341      if (phi_part < 0) phi_part += 2*TMath::Pi();
342      
343
344      if(phi_part < phi_min*torad || phi_part > phi_max*torad) {continue;}
345      if(ptg<0.100) continue;
346      if(fabs(pzg/ptg)>0.999) {continue;}
347      Int_t entries = good_part->GetEntriesFast();
348      good_part->AddLast(p);
349      particle_id[entries] = i;
350    }
351   delete [] good;
352   LOG(AliL3Log::kInformational,"AliL3Evaluate::DefineGoodTracks","NPart")
353   <<AliL3Log::kDec<<"Found "<<good_part->GetEntriesFast()
354   <<" good Particles (Tracks) out of "<<fParticles->GetEntriesFast()<<ENDLOG;
355
356   return good_part;
357 }
358
359 void AliL3Evaluate::EvaluatePatch(Int_t slice,Int_t patch,Int_t min_points,Int_t good_number)
360 {
361   //Make efficiency plots for tracking on patch level (before any merging).
362   
363   Int_t row[5][2] = {{ 0, 45},{46,77},{78,109},{110,141},{142,173}};
364   Int_t *particle_id = new Int_t[fParticles->GetEntriesFast()];
365   TObjArray *good_particles = DefineGoodTracks(slice,row[patch],good_number,particle_id);
366   SetMinPoints(min_points);
367   AssignIDs();
368   CreateHistos();
369   FillEffHistos(good_particles,particle_id);
370   CalcEffHistos();
371   delete good_particles;
372   delete [] particle_id;
373 }
374
375 void AliL3Evaluate::EvaluateSlice(Int_t slice,Int_t min_points,Int_t good_number)
376 {
377   //Make efficiency plots for tracking on a slice (after merging).
378   //min_points = minimum points on track to be considered for evaluation
379   //good_number = minimum hits (padrows) produced by simulated track for consideration.
380
381   Int_t row[2] = {0,173};
382   Int_t *particle_id = new Int_t[fParticles->GetEntriesFast()];
383   TObjArray *good_particles = DefineGoodTracks(slice,row,good_number,particle_id);
384   SetMinPoints(min_points);
385   AssignIDs();
386   CreateHistos();
387   FillEffHistos(good_particles,particle_id);
388   CalcEffHistos();
389   delete good_particles;
390   delete [] particle_id;
391 }
392
393 void AliL3Evaluate::EvaluateGlobal(Int_t min_points,Int_t good_number)
394 {
395   //Make efficiency plots for tracking on several slices.
396   
397   Int_t row[2] = {0,173};
398   Int_t *particle_id = new Int_t[fParticles->GetEntriesFast()];
399   SetMinPoints(min_points);
400   AssignIDs();
401   CreateHistos();
402   for(Int_t slice=fMinSlice;slice<=fMaxSlice;slice++){
403     TObjArray *good_particles = DefineGoodTracks(slice,row,good_number,particle_id);
404     FillEffHistos(good_particles,particle_id);
405     delete good_particles;
406   }
407   CalcEffHistos();
408 }
409
410 void AliL3Evaluate::FillEffHistos(TObjArray *good_particles,Int_t *particle_id)
411 {  
412   //Fill the efficiency histograms.
413
414   for(Int_t i=0; i<good_particles->GetEntriesFast(); i++)
415     {
416       TParticle *p = (TParticle*)good_particles->UncheckedAt(i);
417       Double_t ptg=p->Pt(),pzg=p->Pz();
418       Float_t dipangle=TMath::ATan2(pzg,ptg)*180./TMath::Pi();
419       fNGoodTracksPt->Fill(ptg);
420       fNGoodTracksEta->Fill(dipangle);
421       Int_t found = 0;
422       for(Int_t k=0; k<fTracks->GetNTracks(); k++)
423         {
424           AliL3Track *track = fTracks->GetCheckedTrack(k);
425           if(!track) continue;
426           Int_t nHits = track->GetNumberOfPoints();
427           if(nHits < fMinPointsOnTrack) break;
428           
429           Int_t tracklabel;
430           tracklabel = track->GetMCid();
431           
432           if(TMath::Abs(tracklabel) != particle_id[i]) continue;
433           found=1;
434           if(tracklabel == particle_id[i]) {fNFoundTracksPt->Fill(ptg); fNFoundTracksEta->Fill(dipangle);}
435           else {fNFakeTracksPt->Fill(ptg); fNFakeTracksEta->Fill(dipangle);}
436           Float_t pt=track->GetPt();
437           fPtRes->Fill((pt-ptg)/ptg*100.);
438           break;
439           
440         }
441     }
442 }
443
444 void AliL3Evaluate::CalcEffHistos(){  
445   
446   Stat_t ngood=fNGoodTracksPt->GetEntries();
447   Stat_t nfound=fNFoundTracksPt->GetEntries();
448   Stat_t nfake=fNFakeTracksPt->GetEntries();
449   
450   LOG(AliL3Log::kInformational,"AliL3Evaluate::FillEffHistos","Efficiency")
451     <<AliL3Log::kDec<<"There was "<<ngood<<" generated good tracks"<<ENDLOG;
452   LOG(AliL3Log::kInformational,"AliL3Evaluate::FillEffHistos","Efficiency")
453     <<AliL3Log::kDec<<"Found "<<nfound<<" tracks"<<ENDLOG;
454   LOG(AliL3Log::kInformational,"AliL3Evaluate::FillEffHistos","Efficiency")
455     <<AliL3Log::kDec<<"Integral efficiency is about "<<nfound/ngood*100<<ENDLOG;
456   LOG(AliL3Log::kInformational,"AliL3Evaluate::FillEffHistos","Efficiency")
457     <<AliL3Log::kDec<<"Fake tracks relative is about "<<nfake/ngood*100<<ENDLOG;
458   
459   fNFoundTracksPt->Sumw2(); fNGoodTracksPt->Sumw2();
460   fTrackEffPt->Divide(fNFoundTracksPt,fNGoodTracksPt,1,1.,"b");
461   fFakeTrackEffPt->Divide(fNFakeTracksPt,fNGoodTracksPt,1,1.,"b");
462   fTrackEffPt->SetMaximum(1.4);
463   fTrackEffPt->SetXTitle("P_{T} [GeV]");
464   fTrackEffPt->SetLineWidth(2);
465   fFakeTrackEffPt->SetFillStyle(3013);
466   fTrackEffPt->SetLineColor(4);
467   fFakeTrackEffPt->SetFillColor(2);
468   
469   fNFoundTracksEta->Sumw2(); fNGoodTracksEta->Sumw2();
470   fTrackEffEta->Divide(fNFoundTracksEta,fNGoodTracksEta,1,1.,"b");
471   fFakeTrackEffEta->Divide(fNFakeTracksEta,fNGoodTracksEta,1,1.,"b");
472   fTrackEffEta->SetMaximum(1.4);
473   fTrackEffEta->SetXTitle("#lambda [degrees]");
474   fTrackEffEta->SetLineWidth(2);
475   fFakeTrackEffEta->SetFillStyle(3013);
476   fTrackEffEta->SetLineColor(4);
477   fFakeTrackEffEta->SetFillColor(2);
478        
479 }
480
481 void AliL3Evaluate::AssignIDs()
482 {
483   //Assign MC id to the tracks.
484   
485   UInt_t *index=0,tmp_ind=0;
486   Int_t *pID=0,npoints=0;
487   
488   if(!fIsSlow)
489     {
490       pID = GetFastIDs(tmp_ind,npoints);
491       index = (UInt_t*)tmp_ind;
492     }
493   
494   fTracks->QSort();
495   LOG(AliL3Log::kDebug,"AliL3Evaluate::AssignIDs","Track Loop")
496     <<"Assigning MC id to the found tracks...."<<ENDLOG;
497   for(Int_t i=0; i<fTracks->GetNTracks(); i++)
498     {
499       AliL3Track *track = (AliL3Track*)fTracks->GetCheckedTrack(i);
500       if(!track) 
501         {
502           LOG(AliL3Log::kWarning,"AliL3Evaluate::AssignIDs","Track Loop")
503             <<AliL3Log::kDec<<"No track in track array, index "<<i<<ENDLOG;
504           continue;
505         }
506       if(track->GetNumberOfPoints() < fMinPointsOnTrack) break;
507       Int_t tID;
508       if(!fIsSlow)
509         tID = GetMCTrackLabel(track,index,pID,npoints);
510       else
511         tID = GetMCTrackLabel(track);
512       track->SetMCid(tID);
513       printf("track %i id %d\n",i,tID);
514     }
515   
516   if(pID)   delete [] pID;
517   if(index) delete [] index;
518 }
519
520
521 struct S {Int_t lab; Int_t max;};
522 Int_t AliL3Evaluate::GetMCTrackLabel(AliL3Track *track,UInt_t *index,Int_t *pID,Int_t npoints) 
523 {
524   //Returns the MCtrackID of the belonging clusters.
525   //If MCLabel < 0, means that track is fake.
526   //Fake track means that more than 10 percent of clusters are assigned incorrectly.
527   
528   Int_t num_of_clusters = track->GetNumberOfPoints();
529   S *s=new S[num_of_clusters];
530   Int_t i;
531   for (i=0; i<num_of_clusters; i++) s[i].lab=s[i].max=0;
532   UInt_t *hitnum = track->GetHitNumbers();  
533   UInt_t id;
534
535   Int_t **trackID;
536     
537   if(fIsSlow)
538     trackID = GetClusterIDs(track);
539   else
540     trackID = GetClusterIDs(track,index,pID,npoints);
541   
542
543   Int_t lab=123456789;
544   for (i=0; i<num_of_clusters; i++) 
545     {
546       //Tricks to get the clusters belonging to this track:
547       id = hitnum[i];
548       Int_t slice = (id>>25) & 0x7f;
549       Int_t patch = (id>>22) & 0x7;
550       UInt_t pos = id&0x3fffff;       
551       
552       AliL3SpacePointData *points = fClusters[slice][patch];
553       
554       if(!points) 
555         continue;
556       
557       if(pos>=fNcl[slice][patch]) 
558         {
559           LOG(AliL3Log::kError,"AliL3Evaluate::GetMCTrackLabel","Clusterarray")
560             <<AliL3Log::kDec<<"ERROR"<<ENDLOG;
561           continue;
562         }
563       
564       //Get the label of the cluster:
565       //printf("label %d %d %d\n",trackID[i][0],trackID[i][1],trackID[i][2]);
566       lab=abs(trackID[i][0]);
567       Int_t j;
568       for (j=0; j<num_of_clusters; j++)
569         if (s[j].lab==lab || s[j].max==0) break;
570       s[j].lab=lab;
571       s[j].max++;
572     }
573   
574   Int_t max=0;
575   for (i=0; i<num_of_clusters; i++) 
576     if (s[i].max>max) {max=s[i].max; lab=s[i].lab;}
577   
578   delete[] s;
579   
580   for (i=0; i<num_of_clusters; i++) 
581     {
582       id = hitnum[i];
583       Int_t slice = (id>>25) & 0x7f;
584       Int_t patch = (id>>22) & 0x7;
585       UInt_t pos = id&0x3fffff;       
586       
587       AliL3SpacePointData *points = fClusters[slice][patch];
588       
589       if(!points) 
590         continue;
591       
592       if(pos>=fNcl[slice][patch]) 
593         {
594           LOG(AliL3Log::kError,"AliL3Evaluate::GetMCTrackLabel","Clusterarray")
595             <<AliL3Log::kDec<<"ERROR"<<ENDLOG;
596           continue;
597         }
598       
599       if (abs(trackID[i][1]) == lab || 
600           abs(trackID[i][2]) == lab ) max++;
601     }
602   
603   //check if more than 10% of the clusters are incorrectly assigned (fake track):
604   if (1.-Float_t(max)/num_of_clusters > 0.10) 
605     {
606       return -lab;
607     }
608   
609   delete [] trackID;
610   return lab;
611 }
612
613
614 Int_t **AliL3Evaluate::GetClusterIDs(AliL3Track *track,UInt_t *index,Int_t *pID,Int_t npoints)
615 {
616   //Return the MC information of all clusters belonging to track.
617   
618   Int_t num_of_clusters = track->GetNumberOfPoints();
619   Int_t **trackID = new Int_t*[num_of_clusters];
620   
621   UInt_t *hitnum = track->GetHitNumbers();  
622   UInt_t id;
623
624   
625   Float_t xyz[3];
626   Int_t sector,padrow;
627   for(Int_t i=0; i<num_of_clusters; i++)
628     {
629       id = hitnum[i];
630       Int_t slice = (id>>25) & 0x7f;
631       Int_t patch = (id>>22) & 0x7;
632       UInt_t pos = id&0x3fffff;       
633       
634       AliL3SpacePointData *points = fClusters[slice][patch];
635       
636       if(!points) 
637         {
638           LOG(AliL3Log::kError,"AliL3Evaluate::GetClusterIDs","Clusterarray")
639             <<"No points at slice "<<slice<<" patch "<<patch<<" pos "<<pos<<ENDLOG;
640           continue;
641         }
642       if(pos>=fNcl[slice][patch]) 
643         {
644           LOG(AliL3Log::kError,"AliL3Evaluate::GetClusterIDs","Clusterarray")
645             <<AliL3Log::kDec<<"ERROR"<<ENDLOG;
646           continue;
647         }
648       
649       xyz[0] = points[pos].fX;
650       xyz[1] = points[pos].fY;
651       xyz[2] = points[pos].fZ;
652       //sector = points[pos].fSector;
653       padrow = points[pos].fPadRow;
654       Int_t se,ro;
655       fTransform->Slice2Sector(slice,padrow,se,ro);
656       fTransform->Global2Raw(xyz,se,ro);
657       
658       if(fIsSlow)
659         {
660           Int_t p = fRowid[slice][padrow];
661           
662           if(!fDigitsTree->GetEvent(p)) 
663             LOG(AliL3Log::kError,"AliL3Evaluate::GetClusterIDs","Digits Tree")
664               <<"Error reading digits tree"<<ENDLOG;
665           
666           trackID[i] = new Int_t[3];
667           trackID[i][0] = fDigits->GetTrackID((Int_t)xyz[2],(Int_t)xyz[1],0);
668           trackID[i][1] = fDigits->GetTrackID((Int_t)xyz[2],(Int_t)xyz[1],1);
669           trackID[i][2] = fDigits->GetTrackID((Int_t)xyz[2],(Int_t)xyz[1],2);
670         }
671       else
672         {
673           Int_t tmp_pid=0;
674           for(Int_t ii=0; ii<npoints; ii++)
675             {
676               tmp_pid = pID[ii];
677               if(index[ii] == id) break;
678             }
679           trackID[i] = new Int_t[3];
680           trackID[i][0] = tmp_pid;
681           trackID[i][1] = -1;
682           trackID[i][2] = -1;
683         }
684     }
685   return trackID;
686 }
687
688 Int_t *AliL3Evaluate::GetFastIDs(UInt_t &tmp_ind,Int_t &npoints)
689 {
690   //Get the MC id of space points in case of using the fast simulator. 
691
692   FILE *infile = fopen("point_mc.dat","r");
693   if(!infile) return 0;
694   Int_t hitid,hitmc,i;
695   
696   for(i=0; ; i++)
697     if(fscanf(infile,"%d %d",&hitid,&hitmc)==EOF) break;
698   npoints = i;
699   rewind(infile);
700   Int_t *pID = new Int_t[npoints];
701   UInt_t *ind = new UInt_t[npoints];
702   tmp_ind = (UInt_t)ind;
703   
704   for(i=0; i<npoints; i++)
705     {
706       if(fscanf(infile,"%d %d",&hitid,&hitmc)==EOF) break;
707       pID[i] = hitmc;
708       ind[i] = hitid;
709     }
710   fclose(infile);
711
712   return pID;
713
714 }
715
716 void AliL3Evaluate::Write2File(Char_t *outputfile)
717 {
718   //Write histograms to file:
719   
720   TFile *of = new TFile(outputfile,"RECREATE");
721   if(!of->IsOpen())
722     {
723       LOG(AliL3Log::kError,"AliL3Evaluate::Write2File","File Open")
724         <<"Problems opening rootfile"<<ENDLOG;
725       return;
726     }
727   
728   of->cd();
729   fPtRes->Write();
730   fNGoodTracksPt->Write();
731   fNFoundTracksPt->Write();
732   fNFakeTracksPt->Write();
733   fTrackEffPt->Write();
734   fFakeTrackEffPt->Write();
735   fNGoodTracksEta->Write();
736   fNFoundTracksEta->Write();
737   fNFakeTracksEta->Write();
738   fTrackEffEta->Write();
739   fFakeTrackEffEta->Write();
740   
741   of->Close();
742   delete of;
743
744 }
745