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