]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/tracking-ca/AliHLTTPCCAPerformance.cxx
f06d70c15aa2f095f458ded584f40af7f9a1d6d6
[u/mrichter/AliRoot.git] / HLT / TPCLib / tracking-ca / AliHLTTPCCAPerformance.cxx
1 // $Id$
2 // **************************************************************************
3 // This file is property of and copyright by the ALICE HLT Project          * 
4 // ALICE Experiment at CERN, All rights reserved.                           *
5 //                                                                          *
6 // Primary Authors: Sergey Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de> *
7 //                  Ivan Kisel <kisel@kip.uni-heidelberg.de>                *
8 //                  for The ALICE HLT Project.                              *
9 //                                                                          *
10 // Permission to use, copy, modify and distribute this software and its     *
11 // documentation strictly for non-commercial purposes is hereby granted     *
12 // without fee, provided that the above copyright notice appears in all     *
13 // copies and that both the copyright notice and this permission notice     *
14 // appear in the supporting documentation. The authors make no claims       *
15 // about the suitability of this software for any purpose. It is            *
16 // provided "as is" without express or implied warranty.                    *
17 //                                                                          *
18 //***************************************************************************
19
20
21 #include "AliHLTTPCCAPerformance.h"
22 #include "AliHLTTPCCAGBHit.h"
23 #include "AliHLTTPCCAMCTrack.h"
24 #include "AliHLTTPCCAMCPoint.h"
25 #include "AliHLTTPCCAOutTrack.h"
26 #include "AliHLTTPCCAGBTrack.h"
27 #include "AliHLTTPCCAGBTracker.h"
28 #include "AliHLTTPCCATracker.h"
29 #include "AliHLTTPCCATracklet.h"
30
31
32 #include "TMath.h"
33 #include "TROOT.h"
34 #include "Riostream.h"
35 #include "TFile.h"
36 #include "TH1.h"
37 #include "TH2.h"
38 #include "TProfile.h"
39
40 AliHLTTPCCAPerformance &AliHLTTPCCAPerformance::Instance()
41 {
42   // reference to static object
43   static AliHLTTPCCAPerformance gAliHLTTPCCAPerformance;
44   return gAliHLTTPCCAPerformance;
45 }
46
47 AliHLTTPCCAPerformance::AliHLTTPCCAPerformance()
48   : 
49   fTracker(0),
50   fHitLabels(0),
51   fNHits(0),                   
52   fMCTracks(0),  
53   fNMCTracks(0),               
54   fMCPoints(0),  
55   fNMCPoints(0),               
56   fDoClusterPulls(0),         
57   fStatNEvents(0),
58   fStatTime(0),
59   fStatSeedNRecTot(0),
60   fStatSeedNRecOut(0),
61   fStatSeedNGhost(0),
62   fStatSeedNMCAll(0),
63   fStatSeedNRecAll(0),
64   fStatSeedNClonesAll(0),
65   fStatSeedNMCRef(0), 
66   fStatSeedNRecRef(0), 
67   fStatSeedNClonesRef(0),  
68   fStatCandNRecTot(0),
69   fStatCandNRecOut(0),
70   fStatCandNGhost(0),
71   fStatCandNMCAll(0),
72   fStatCandNRecAll(0),
73   fStatCandNClonesAll(0),
74   fStatCandNMCRef(0), 
75   fStatCandNRecRef(0), 
76   fStatCandNClonesRef(0),  
77   fStatNRecTot(0),
78   fStatNRecOut(0),
79   fStatNGhost(0),
80   fStatNMCAll(0),
81   fStatNRecAll(0),
82   fStatNClonesAll(0),
83   fStatNMCRef(0), 
84   fStatNRecRef(0), 
85   fStatNClonesRef(0),  
86   fStatGBNRecTot(0), 
87   fStatGBNRecOut(0), 
88   fStatGBNGhost(0),
89   fStatGBNMCAll(0),
90   fStatGBNRecAll(0),
91   fStatGBNClonesAll(0),
92   fStatGBNMCRef(0), 
93   fStatGBNRecRef(0), 
94   fStatGBNClonesRef(0),  
95   fHistoDir(0), 
96   fhResY(0),       
97   fhResZ(0),       
98   fhResSinPhi(0),  
99   fhResDzDs(0),    
100   fhResPt(0),      
101   fhPullY(0),      
102   fhPullZ(0),      
103   fhPullSinPhi(0), 
104   fhPullDzDs(0),
105   fhPullQPt(0), 
106   fhPullYS(0),      
107   fhPullZT(0),      
108   fhHitErrY(0), 
109   fhHitErrZ(0),
110   fhHitResY(0),
111   fhHitResZ(0),
112   fhHitPullY(0),
113   fhHitPullZ(0), 
114   fhHitShared(0),
115   fhHitResY1(0),
116   fhHitResZ1(0),
117   fhHitPullY1(0),
118   fhHitPullZ1(0),
119   fhCellPurity(0),
120   fhCellNHits(0),
121   fhCellPurityVsN(0), 
122   fhCellPurityVsPt(0),
123   fhEffVsP(0), 
124   fhSeedEffVsP(0),
125   fhCandEffVsP(0),
126   fhGBEffVsP(0),
127   fhGBEffVsPt(0),
128   fhNeighQuality(0),
129   fhNeighEff(0),
130   fhNeighQualityVsPt(0),
131   fhNeighEffVsPt(0),
132   fhNeighDy(0),
133   fhNeighDz(0),
134   fhNeighChi(0),
135   fhNeighDyVsPt(0),
136   fhNeighDzVsPt(0),
137   fhNeighChiVsPt(0), 
138   fhNeighNCombVsArea(0),
139   fhNHitsPerSeed (0),
140   fhNHitsPerTrackCand(0),
141   fhTrackLengthRef(0),
142   fhRefRecoX(0),
143   fhRefRecoY(0),
144   fhRefRecoZ(0),
145   fhRefRecoP(0),
146   fhRefRecoPt(0),
147   fhRefRecoAngleY(0),
148   fhRefRecoAngleZ(0),
149   fhRefRecoNHits(0),
150   fhRefNotRecoX(0),
151   fhRefNotRecoY(0),
152   fhRefNotRecoZ(0),
153   fhRefNotRecoP(0),
154   fhRefNotRecoPt(0),
155   fhRefNotRecoAngleY(0),
156   fhRefNotRecoAngleZ(0),
157   fhRefNotRecoNHits(0)
158 {
159   //* constructor
160 }
161
162
163 AliHLTTPCCAPerformance::AliHLTTPCCAPerformance(const AliHLTTPCCAPerformance&)
164   :
165   fTracker(0),
166   fHitLabels(0),
167   fNHits(0),                   
168   fMCTracks(0),  
169   fNMCTracks(0),               
170   fMCPoints(0),  
171   fNMCPoints(0),               
172   fDoClusterPulls(0),         
173   fStatNEvents(0),
174   fStatTime(0),
175   fStatSeedNRecTot(0),
176   fStatSeedNRecOut(0),
177   fStatSeedNGhost(0),
178   fStatSeedNMCAll(0),
179   fStatSeedNRecAll(0),
180   fStatSeedNClonesAll(0),
181   fStatSeedNMCRef(0), 
182   fStatSeedNRecRef(0), 
183   fStatSeedNClonesRef(0),  
184   fStatCandNRecTot(0),
185   fStatCandNRecOut(0),
186   fStatCandNGhost(0),
187   fStatCandNMCAll(0),
188   fStatCandNRecAll(0),
189   fStatCandNClonesAll(0),
190   fStatCandNMCRef(0), 
191   fStatCandNRecRef(0), 
192   fStatCandNClonesRef(0),  
193   fStatNRecTot(0),
194   fStatNRecOut(0),
195   fStatNGhost(0),
196   fStatNMCAll(0),
197   fStatNRecAll(0),
198   fStatNClonesAll(0),
199   fStatNMCRef(0), 
200   fStatNRecRef(0), 
201   fStatNClonesRef(0),  
202   fStatGBNRecTot(0), 
203   fStatGBNRecOut(0), 
204   fStatGBNGhost(0),
205   fStatGBNMCAll(0),
206   fStatGBNRecAll(0),
207   fStatGBNClonesAll(0),
208   fStatGBNMCRef(0), 
209   fStatGBNRecRef(0), 
210   fStatGBNClonesRef(0),  
211   fHistoDir(0), 
212   fhResY(0),       
213   fhResZ(0),       
214   fhResSinPhi(0),  
215   fhResDzDs(0),    
216   fhResPt(0),      
217   fhPullY(0),      
218   fhPullZ(0),      
219   fhPullSinPhi(0), 
220   fhPullDzDs(0),
221   fhPullQPt(0), 
222   fhPullYS(0),      
223   fhPullZT(0),      
224   fhHitErrY(0), 
225   fhHitErrZ(0),
226   fhHitResY(0),
227   fhHitResZ(0),
228   fhHitPullY(0),
229   fhHitPullZ(0),  
230   fhHitShared(0),
231   fhHitResY1(0),
232   fhHitResZ1(0),
233   fhHitPullY1(0),
234   fhHitPullZ1(0),
235   fhCellPurity(0),
236   fhCellNHits(0),
237   fhCellPurityVsN(0), 
238   fhCellPurityVsPt(0),
239   fhEffVsP(0), 
240   fhSeedEffVsP(0),
241   fhCandEffVsP(0),
242   fhGBEffVsP(0),
243   fhGBEffVsPt(0),
244   fhNeighQuality(0),
245   fhNeighEff(0),
246   fhNeighQualityVsPt(0),
247   fhNeighEffVsPt(0),
248   fhNeighDy(0),
249   fhNeighDz(0),
250   fhNeighChi(0),
251   fhNeighDyVsPt(0),
252   fhNeighDzVsPt(0),
253   fhNeighChiVsPt(0), 
254   fhNeighNCombVsArea(0),
255   fhNHitsPerSeed (0),
256   fhNHitsPerTrackCand(0),
257   fhTrackLengthRef(0),
258   fhRefRecoX(0),
259   fhRefRecoY(0),
260   fhRefRecoZ(0),
261   fhRefRecoP(0),
262   fhRefRecoPt(0),
263   fhRefRecoAngleY(0),
264   fhRefRecoAngleZ(0),
265   fhRefRecoNHits(0),
266   fhRefNotRecoX(0),
267   fhRefNotRecoY(0),
268   fhRefNotRecoZ(0),
269   fhRefNotRecoP(0),
270   fhRefNotRecoPt(0),
271   fhRefNotRecoAngleY(0),
272   fhRefNotRecoAngleZ(0),
273   fhRefNotRecoNHits(0)
274 {
275   //* dummy
276 }
277
278 const AliHLTTPCCAPerformance &AliHLTTPCCAPerformance::operator=(const AliHLTTPCCAPerformance&) const
279 {
280   //* dummy
281   return *this;
282 }
283
284 AliHLTTPCCAPerformance::~AliHLTTPCCAPerformance()
285 {
286   //* destructor
287   StartEvent();
288 }
289
290 void AliHLTTPCCAPerformance::SetTracker( AliHLTTPCCAGBTracker * const Tracker )
291 {
292   //* set pointer to HLT CA Global tracker
293   fTracker = Tracker;
294 }
295
296 void AliHLTTPCCAPerformance::StartEvent()
297 {
298   //* clean up arrays
299   if( !fHistoDir )  CreateHistos();
300   if( fHitLabels ) delete[] fHitLabels;
301   fHitLabels = 0;
302   fNHits = 0;
303   if( fMCTracks ) delete[] fMCTracks;
304   fMCTracks = 0;
305   fNMCTracks = 0;
306   if( fMCPoints ) delete[] fMCPoints;
307   fMCPoints = 0;
308   fNMCPoints = 0;
309 }
310
311 void AliHLTTPCCAPerformance::SetNHits( Int_t NHits )
312 {
313   //* set number of hits
314   if( fHitLabels ) delete[] fHitLabels;
315   fHitLabels = 0;
316   fHitLabels = new AliHLTTPCCAHitLabel[ NHits ];
317   fNHits = NHits;
318 }  
319
320 void AliHLTTPCCAPerformance::SetNMCTracks( Int_t NumberOfMCTracks )
321 {
322   //* set number of MC tracks
323   if( fMCTracks ) delete[] fMCTracks;
324   fMCTracks = 0;
325   fMCTracks = new AliHLTTPCCAMCTrack[ NumberOfMCTracks ];
326   fNMCTracks = NumberOfMCTracks;
327 }  
328
329 void AliHLTTPCCAPerformance::SetNMCPoints( Int_t NMCPoints )
330 {
331   //* set number of MC points
332   if( fMCPoints ) delete[] fMCPoints;
333   fMCPoints = 0;
334   fMCPoints = new AliHLTTPCCAMCPoint[ NMCPoints ];
335   fNMCPoints = 0;
336 }  
337
338 void AliHLTTPCCAPerformance::ReadHitLabel( Int_t HitID, 
339                                            Int_t lab0, Int_t lab1, Int_t lab2 )
340 {
341   //* read the hit labels
342   AliHLTTPCCAHitLabel hit;
343   hit.fLab[0] = lab0;
344   hit.fLab[1] = lab1;
345   hit.fLab[2] = lab2;
346   fHitLabels[HitID] = hit;
347 }
348
349 void AliHLTTPCCAPerformance::ReadMCTrack( Int_t index, const TParticle *part )
350 {
351   //* read mc track to the local array
352   fMCTracks[index] = AliHLTTPCCAMCTrack(part);
353 }
354
355 void AliHLTTPCCAPerformance::ReadMCTPCTrack( Int_t index, Float_t X, Float_t Y, Float_t Z, 
356                                              Float_t Px, Float_t Py, Float_t Pz )
357 {
358   //* read mc track parameters at TPC
359   fMCTracks[index].SetTPCPar(X,Y,Z,Px,Py,Pz);
360 }
361
362 void AliHLTTPCCAPerformance::ReadMCPoint( Int_t TrackID, Float_t X, Float_t Y, Float_t Z, Float_t Time, Int_t iSlice )
363 {
364   //* read mc point to the local array
365   AliHLTTPCCAMCPoint &p = fMCPoints[fNMCPoints];
366   p.SetTrackID( TrackID );
367   p.SetX( X );
368   p.SetY( Y );
369   p.SetZ( Z );
370   p.SetTime( Time );
371   p.SetISlice( iSlice );
372   Float_t sx, sy, sz;
373   fTracker->Slices()[iSlice].Param().Global2Slice( X, Y, Z, &sx, &sy, &sz );
374   p.SetSx(sx);
375   p.SetSy(sy); 
376   p.SetSz(sz); 
377   if( X*X + Y*Y>10.) fNMCPoints++;
378 }
379
380 void AliHLTTPCCAPerformance::CreateHistos()
381 {
382   //* create performance histogramms
383   TDirectory *curdir = gDirectory;
384   fHistoDir = gROOT->mkdir("HLTTPCCATrackerPerformance");
385   fHistoDir->cd();
386
387   gDirectory->mkdir("Links");
388   gDirectory->cd("Links");  
389
390   fhLinkEff[0] = new TProfile("fhLinkEffPrimRef", "fhLinkEffPrimRef vs row", 156, 2., 158.);
391   fhLinkEff[1] = new TProfile("fhLinkEffPrimExt", "fhLinkEffPrimExt vs row", 156, 2., 158.);
392   fhLinkEff[2] = new TProfile("fhLinkEffSecRef", "fhLinkEffSecRef vs row", 156, 2., 158.);
393   fhLinkEff[3] = new TProfile("fhLinkEffSecExt", "fhLinkEffSecExt vs row", 156, 2., 158.);
394   fhLinkAreaY[0] = new TH1D("fhLinkAreaYPrimRef","fhLinkAreaYPrimRef",100,0,10);
395   fhLinkAreaZ[0] = new TH1D("fhLinkAreaZPrimRef","fhLinkAreaZPrimRef",100,0,10);
396   fhLinkAreaY[1] = new TH1D("fhLinkAreaYPrimExt","fhLinkAreaYPrimExt",100,0,10);
397   fhLinkAreaZ[1] = new TH1D("fhLinkAreaZPrimExt","fhLinkAreaZPrimExt",100,0,10);
398   fhLinkAreaY[2] = new TH1D("fhLinkAreaYSecRef","fhLinkAreaYSecRef",100,0,10);
399   fhLinkAreaZ[2] = new TH1D("fhLinkAreaZSecRef","fhLinkAreaZSecRef",100,0,10);
400   fhLinkAreaY[3] = new TH1D("fhLinkAreaYSecExt","fhLinkAreaYSecExt",100,0,10);
401   fhLinkAreaZ[3] = new TH1D("fhLinkAreaZSecExt","fhLinkAreaZSecExt",100,0,10);
402   fhLinkChiRight[0] = new TH1D("fhLinkChiRightPrimRef","fhLinkChiRightPrimRef",100,0,10);
403   fhLinkChiRight[1] = new TH1D("fhLinkChiRightPrimExt","fhLinkChiRightPrimExt",100,0,10);
404   fhLinkChiRight[2] = new TH1D("fhLinkChiRightSecRef","fhLinkChiRightSecRef",100,0,10);
405   fhLinkChiRight[3] = new TH1D("fhLinkChiRightSecExt","fhLinkChiRightSecExt",100,0,10);
406   fhLinkChiWrong[0] = new TH1D("fhLinkChiWrongPrimRef","fhLinkChiWrongPrimRef",100,0,10);
407   fhLinkChiWrong[1] = new TH1D("fhLinkChiWrongPrimExt","fhLinkChiWrongPrimExt",100,0,10);
408   fhLinkChiWrong[2] = new TH1D("fhLinkChiWrongSecRef","fhLinkChiWrongSecRef",100,0,10);
409   fhLinkChiWrong[3] = new TH1D("fhLinkChiWrongSecExt","fhLinkChiWrongSecExt",100,0,10);
410
411   gDirectory->cd(".."); 
412
413   gDirectory->mkdir("Neighbours");
414   gDirectory->cd("Neighbours");  
415   
416   fhNeighQuality = new TProfile("NeighQuality", "Neighbours Quality vs row", 160, 0., 160.); 
417   fhNeighEff = new TProfile("NeighEff", "Neighbours Efficiency vs row", 160, 0., 160.); 
418   fhNeighQualityVsPt = new TProfile("NeighQualityVsPt", "Neighbours Quality vs Pt", 100, 0., 5.); 
419   fhNeighEffVsPt = new TProfile("NeighEffVsPt", "Neighbours Efficiency vs Pt", 100, 0., 5.); 
420   fhNeighDy = new TH1D("NeighDy","Neighbours dy",100,-10,10);
421   fhNeighDz =  new TH1D("NeighDz","Neighbours dz",100,-10,10);
422   fhNeighChi = new TH1D("NeighChi","Neighbours chi",100,0,20);
423
424   fhNeighDyVsPt = new TH2D("NeighDyVsPt","NeighDyVsPt", 100,0,5, 100, -20,20);
425   fhNeighDzVsPt = new TH2D("NeighDzVsPt","NeighDzVsPt", 100,0,5, 100, -20,20);
426   fhNeighChiVsPt = new TH2D("NeighChiVsPt","NeighChiVsPt", 100,0,5, 100, 0,40);
427   fhNeighNCombVsArea = new TH2D("NeighNCombVsArea","NeighNCombVsArea", 15,0,3, 40, 0,40);
428
429   gDirectory->cd(".."); 
430
431   gDirectory->mkdir("Tracklets");
432   gDirectory->cd("Tracklets");
433
434   fhNHitsPerSeed = new TH1D("NHitsPerSeed","NHitsPerSeed", 160,0,160);
435   fhSeedEffVsP = new TProfile("fhSeedEffVsP", "Track Seed Eff vs P", 100, 0., 5.);
436
437   gDirectory->cd(".."); 
438
439   gDirectory->mkdir("TrackCandidates");
440   gDirectory->cd("TrackCandidates");
441
442   fhNHitsPerTrackCand = new TH1D("NHitsPerTrackCand","NHitsPerTrackCand", 160,0,160);
443   fhCandEffVsP = new TProfile("fhCandEffVsP", "Track Candidate Eff vs P", 100, 0., 5.);
444
445   gDirectory->cd(".."); 
446   
447   gDirectory->mkdir("Tracks");
448   gDirectory->cd("Tracks");  
449
450   fhTrackLengthRef = new TH1D("TrackLengthRef", "TrackLengthRef", 100,0,1);
451
452   fhRefRecoX = new TH1D("fhRefRecoX","fhRefRecoX",100,0,200.);
453   fhRefRecoY = new TH1D("fhRefRecoY","fhRefRecoY",100,-200,200.);
454   fhRefRecoZ = new TH1D("fhRefRecoZ","fhRefRecoZ",100,-250,250.);
455
456
457   fhRefRecoP = new TH1D("fhRefRecoP","fhRefRecoP",100,0,10.);
458   fhRefRecoPt = new TH1D("fhRefRecoPt","fhRefRecoPt",100,0,10.);
459   fhRefRecoAngleY = new TH1D("fhRefRecoAngleY","fhRefRecoAngleY",100,-180.,180.);
460   fhRefRecoAngleZ = new TH1D("fhRefRecoAngleZ","fhRefRecoAngleZ",100,-180.,180);
461   fhRefRecoNHits = new TH1D("fhRefRecoNHits","fhRefRecoNHits",100,0.,200);
462
463   fhRefNotRecoX = new TH1D("fhRefNotRecoX","fhRefNotRecoX",100,0,200.);
464   fhRefNotRecoY = new TH1D("fhRefNotRecoY","fhRefNotRecoY",100,-200,200.);
465   fhRefNotRecoZ = new TH1D("fhRefNotRecoZ","fhRefNotRecoZ",100,-250,250.);
466
467
468   fhRefNotRecoP = new TH1D("fhRefNotRecoP","fhRefNotRecoP",100,0,10.);
469   fhRefNotRecoPt = new TH1D("fhRefNotRecoPt","fhRefNotRecoPt",100,0,10.);
470   fhRefNotRecoAngleY = new TH1D("fhRefNotRecoAngleY","fhRefNotRecoAngleY",100,-180.,180.);
471   fhRefNotRecoAngleZ = new TH1D("fhRefNotRecoAngleZ","fhRefNotRecoAngleZ",100,-180.,180);
472   fhRefNotRecoNHits = new TH1D("fhRefNotRecoNHits","fhRefNotRecoNHits",100,0.,200);
473   
474   gDirectory->cd(".."); 
475
476   gDirectory->mkdir("TrackFit");
477   gDirectory->cd("TrackFit");  
478   
479   fhResY = new TH1D("resY", "track Y resoltion [cm]", 30, -.5, .5);
480   fhResZ = new TH1D("resZ", "track Z resoltion [cm]", 30, -.5, .5);
481   fhResSinPhi = new TH1D("resSinPhi", "track SinPhi resoltion ", 30, -.03, .03);
482   fhResDzDs = new TH1D("resDzDs", "track DzDs resoltion ", 30, -.01, .01);
483   fhResPt = new TH1D("resPt", "track telative Pt resoltion", 30, -.2, .2);
484   fhPullY = new TH1D("pullY", "track Y pull", 30, -10., 10.);
485   fhPullZ = new TH1D("pullZ", "track Z pull", 30, -10., 10.);
486   fhPullSinPhi = new TH1D("pullSinPhi", "track SinPhi pull", 30, -10., 10.);
487   fhPullDzDs = new TH1D("pullDzDs", "track DzDs pull", 30, -10., 10.);
488   fhPullQPt = new TH1D("pullQPt", "track Q/Pt pull", 30, -10., 10.);
489   fhPullYS = new TH1D("pullYS", "track Y+SinPhi chi deviation", 100, 0., 30.);
490   fhPullZT = new TH1D("pullZT", "track Z+DzDs chi deviation ", 100, 0., 30.);
491
492   gDirectory->cd("..");  
493
494   fhEffVsP = new TProfile("EffVsP", "Eff vs P", 100, 0., 5.);
495   fhGBEffVsP = new TProfile("GBEffVsP", "Global tracker: Eff vs P", 100, 0., 5.);
496   fhGBEffVsPt = new TProfile("GBEffVsPt", "Global tracker: Eff vs Pt", 100, 0.2, 5.);
497
498   gDirectory->mkdir("Clusters");
499   gDirectory->cd("Clusters");  
500
501   fhHitShared = new TProfile("fhHitSharedf", "fhHitShared vs row", 160, 0., 160.);
502
503   fhHitResY = new TH1D("resHitY", "Y cluster resoltion [cm]", 100, -2., 2.);
504   fhHitResZ = new TH1D("resHitZ", "Z cluster resoltion [cm]", 100, -2., 2.);
505   fhHitPullY = new TH1D("pullHitY", "Y cluster pull", 50, -10., 10.);
506   fhHitPullZ = new TH1D("pullHitZ", "Z cluster pull", 50, -10., 10.);
507
508   fhHitResY1 = new TH1D("resHitY1", "Y cluster resoltion [cm]", 100, -2., 2.);
509   fhHitResZ1 = new TH1D("resHitZ1", "Z cluster resoltion [cm]", 100, -2., 2.);
510   fhHitPullY1 = new TH1D("pullHitY1", "Y cluster pull", 50, -10., 10.);
511   fhHitPullZ1 = new TH1D("pullHitZ1", "Z cluster pull", 50, -10., 10.);
512
513   fhHitErrY = new TH1D("HitErrY", "Y cluster error [cm]", 100, 0., 1.);
514   fhHitErrZ = new TH1D("HitErrZ", "Z cluster error [cm]", 100, 0., 1.);
515
516   gDirectory->cd("..");  
517
518   gDirectory->mkdir("Cells");
519   gDirectory->cd("Cells");  
520   fhCellPurity = new TH1D("CellPurity", "Cell Purity", 100, -0.1, 1.1);
521   fhCellNHits = new TH1D("CellNHits", "Cell NHits", 40, 0., 40.);
522   fhCellPurityVsN = new TProfile("CellPurityVsN", "Cell purity Vs N hits", 40, 2., 42.);
523   fhCellPurityVsPt = new TProfile("CellPurityVsPt", "Cell purity Vs Pt", 100, 0., 5.);
524   gDirectory->cd("..");  
525
526   curdir->cd();  
527 }
528
529 void AliHLTTPCCAPerformance::WriteDir2Current( TObject *obj )
530 {
531   //* recursive function to copy the directory 'obj' to the current one
532   if( !obj->IsFolder() ) obj->Write();
533   else{
534     TDirectory *cur = gDirectory;
535     TDirectory *sub = cur->mkdir(obj->GetName());
536     sub->cd();
537     TList *listSub = ((TDirectory*)obj)->GetList();
538     TIter it(listSub);
539     while( TObject *obj1=it() ) WriteDir2Current(obj1);
540     cur->cd();
541   }
542 }
543
544 void AliHLTTPCCAPerformance::WriteHistos()
545 {
546   //* write histograms to the file
547   TDirectory *curr = gDirectory;
548   // Open output file and write histograms
549   TFile* outfile = new TFile("HLTTPCCATrackerPerformance.root","RECREATE");
550   outfile->cd();
551   WriteDir2Current(fHistoDir);
552   outfile->Close();
553   curr->cd();
554 }
555
556
557 void AliHLTTPCCAPerformance::LinkPerformance( Int_t iSlice )
558
559   // Efficiency and quality of the found neighbours
560   std::cout<<"Link performance..."<<std::endl;
561   if( !fTracker ) return;  
562   AliHLTTPCCATracker &slice = fTracker->Slices()[iSlice];
563
564   const Int_t *id = slice.HitInputIDs();
565
566   Int_t mcType[fNMCTracks];
567
568   for (Int_t imc=0; imc<fNMCTracks; imc++){   
569     if( fMCTracks[imc].P()<.2 ){  mcType[imc] = -1; continue; }
570     Float_t x = fMCTracks[imc].Par()[0];
571     Float_t y = fMCTracks[imc].Par()[1];
572     //Float_t z = fMCTracks[imc].Par()[2];
573     if( x*x + y*y < 100. ){
574       if( fMCTracks[imc].P()>=1 ) mcType[imc] = 0;
575       else mcType[imc] = 1;
576     } else {
577       if( fMCTracks[imc].P()>=1 ) mcType[imc] = 2;
578       else mcType[imc] = 3;
579     }
580   }
581
582   struct AliHLTTPCCAMCHits{
583     Int_t fNHits;
584     Int_t fID[30];
585   }; 
586   AliHLTTPCCAMCHits *mcGbHitsUp = new AliHLTTPCCAMCHits[fNMCTracks];
587   AliHLTTPCCAMCHits *mcGbHitsDn = new AliHLTTPCCAMCHits[fNMCTracks];
588
589   for( Int_t iRow=2; iRow<slice.Param().NRows()-2; iRow++ ){
590     
591     const AliHLTTPCCARow &row = slice.Row(iRow);
592     const AliHLTTPCCARow &rowUp = slice.Row(iRow+2);
593     const AliHLTTPCCARow &rowDn = slice.Row(iRow-2);    
594     
595     Int_t gbHits  [row.NHits()];
596     Int_t gbHitsUp[rowUp.NHits()];
597     Int_t gbHitsDn[rowDn.NHits()];
598
599     for (Int_t ih = 0; ih<row.NHits()  ; ih++) gbHits  [ih] = fTracker->FirstSliceHit()[iSlice]+id[row  .FirstHit()+ih];    
600     for (Int_t ih = 0; ih<rowUp.NHits(); ih++) gbHitsUp[ih] = fTracker->FirstSliceHit()[iSlice]+id[rowUp.FirstHit()+ih];
601     for (Int_t ih = 0; ih<rowDn.NHits(); ih++) gbHitsDn[ih] = fTracker->FirstSliceHit()[iSlice]+id[rowDn.FirstHit()+ih];
602
603     for (Int_t imc=0; imc<fNMCTracks; imc++){ 
604       mcGbHitsUp[imc].fNHits = 0;
605       mcGbHitsDn[imc].fNHits = 0;
606     }
607   
608     for (Int_t ih = 0; ih<rowUp.NHits(); ih++){      
609       AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[gbHitsUp[ih]].ID()];
610       for( Int_t il=0; il<3; il++ ){
611         Int_t imc = l.fLab[il];
612         if( imc <0 ) break;
613         Int_t &nmc = mcGbHitsUp[imc].fNHits;
614         if( nmc>=30 ) continue;
615         mcGbHitsUp[imc].fID[nmc] = gbHitsUp[ih];
616         nmc++;
617       }
618     }
619
620     for (Int_t ih = 0; ih<rowDn.NHits(); ih++){
621       AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[gbHitsDn[ih]].ID()];    
622       for( Int_t il=0; il<3; il++ ){
623         Int_t imc = l.fLab[il];
624         if( imc <0 ) break;
625         Int_t &nmc = mcGbHitsDn[imc].fNHits;
626         if( nmc>=30 ) continue;
627         mcGbHitsDn[imc].fID[nmc] = gbHitsDn[ih];
628         nmc++;
629       }
630     }
631
632     Short_t *linkUp = ((Short_t*)(slice.RowData() + row.FullOffset())) + row.FullLinkOffset();
633     Short_t *linkDn = linkUp + row.NHits();
634
635     //Float_t dxUp = rowUp.X() - row.X();
636     //Float_t dxDn = row.X() - rowDn.X();
637     Float_t tUp = rowUp.X()/row.X();
638     Float_t tDn = rowDn.X()/row.X();
639
640     for( Int_t ih=0; ih<row.NHits(); ih++ ){      
641       
642       Int_t up = linkUp[ih];
643       Int_t dn = linkDn[ih];
644       
645       AliHLTTPCCAGBHit &h = fTracker->Hits()[gbHits[ih]];
646       AliHLTTPCCAHitLabel &l = fHitLabels[h.ID()];    
647
648       Int_t isMC = -1;
649       Int_t mcFound=-1;
650
651       Float_t yUp = h.Y()*tUp, zUp = h.Z()*tUp;
652       Float_t yDn = h.Y()*tDn, zDn = h.Z()*tDn;
653
654       for( Int_t il=0; il<3; il++ ){
655         Int_t imc = l.fLab[il];
656         if( imc<0 ) break;
657
658         Bool_t isMcUp=0, isMcDn=0;
659
660         Float_t dyMin=1.e8, dzMin=1.e8;
661         for( Int_t i=0; i<mcGbHitsUp[imc].fNHits; i++ ){
662           AliHLTTPCCAGBHit &h1 = fTracker->Hits()[mcGbHitsUp[imc].fID[i]];
663           Float_t dy = TMath::Abs(h1.Y()-yUp);
664           Float_t dz = TMath::Abs(h1.Z()-zUp);
665           if( dy*dy+dz*dz < dyMin*dyMin + dzMin*dzMin ){
666             dyMin = dy;
667             dzMin = dz;     
668           }
669         }
670
671         if( mcType[imc]>=0 && mcGbHitsUp[imc].fNHits>=0 ){
672           fhLinkAreaY[mcType[imc]]->Fill(dyMin);
673           fhLinkAreaZ[mcType[imc]]->Fill(dzMin);
674         }
675         if( dyMin*dyMin+dzMin*dzMin<100. ) isMcUp = 1;
676  
677         dyMin=1.e8; 
678         dzMin=1.e8;
679         for( Int_t i=0; i<mcGbHitsDn[imc].fNHits; i++ ){
680           AliHLTTPCCAGBHit &h1 = fTracker->Hits()[mcGbHitsDn[imc].fID[i]];
681           Float_t dy = TMath::Abs(h1.Y()-yDn);
682           Float_t dz = TMath::Abs(h1.Z()-zDn);
683           if( dy*dy+dz*dz < dyMin*dyMin + dzMin*dzMin ){
684             dyMin = dy;
685             dzMin = dz;     
686           }
687         }
688
689         if( mcType[imc]>=0 && mcGbHitsDn[imc].fNHits>=0 ){
690           fhLinkAreaY[mcType[imc]]->Fill(dyMin);
691           fhLinkAreaZ[mcType[imc]]->Fill(dzMin);
692         }
693         if( dyMin*dyMin+dzMin*dzMin<100. ) isMcDn = 1;
694  
695         if( !isMcUp || !isMcDn ) continue;
696         isMC = imc;
697
698         Bool_t found=0;
699         if( up>=0 && dn>=0 ){
700           //std::cout<<"row, ih, mc, up, dn = "<<iRow<<" "<<ih<<" "<<imc<<" "<<up<<" "<<dn<<std::endl;
701           AliHLTTPCCAGBHit &hUp = fTracker->Hits()[gbHitsUp[up]];
702           AliHLTTPCCAGBHit &hDn = fTracker->Hits()[gbHitsDn[dn]];
703           AliHLTTPCCAHitLabel &lUp = fHitLabels[hUp.ID()];
704           AliHLTTPCCAHitLabel &lDn = fHitLabels[hDn.ID()];
705           Bool_t foundUp=0, foundDn=0;
706           for( Int_t jl=0; jl<3; jl++ ){
707             if( lUp.fLab[jl] == imc ) foundUp=1;
708             if( lDn.fLab[jl] == imc ) foundDn=1;
709             //std::cout<<"mc up, dn = "<<lUp.fLab[jl]<<" "<<lDn.fLab[jl]<<std::endl;
710           }
711           if( foundUp && foundDn ) found = 1;
712         }
713         if( found ){ mcFound = imc; break;}     
714       }
715   
716       if( mcFound>=0 ){
717         //std::cout<<" mc "<<mcFound<<" found"<<std::endl;
718         if( mcType[mcFound]>=0 ) fhLinkEff[mcType[mcFound]]->Fill(iRow, 1 );    
719       }else if( isMC>=0 ){
720         //std::cout<<" mc "<<isMC<<" not found"<<std::endl;
721         if( mcType[isMC]>=0 ) fhLinkEff[mcType[isMC]]->Fill(iRow, 0 );                  
722       }
723
724     } // ih    
725   } // iRow
726   delete[] mcGbHitsUp;
727   delete[] mcGbHitsDn;
728 }
729
730
731 void AliHLTTPCCAPerformance::SliceTrackletPerformance( Int_t iSlice, Bool_t PrintFlag )
732
733   //* calculate slice tracker performance
734   if( !fTracker ) return;
735   
736   Int_t nRecTot = 0, nGhost=0, nRecOut=0;
737   Int_t nMCAll = 0, nRecAll=0, nClonesAll=0;
738   Int_t nMCRef = 0, nRecRef=0, nClonesRef=0;
739   AliHLTTPCCATracker &slice = fTracker->Slices()[iSlice];
740
741   Int_t firstSliceHit = fTracker->FirstSliceHit()[iSlice];
742   Int_t endSliceHit = fTracker->NHits();
743   if( iSlice<fTracker->NSlices()-1 ) endSliceHit = fTracker->FirstSliceHit()[iSlice+1];    
744
745   // Select reconstructable MC tracks
746
747   {
748     for (Int_t imc=0; imc<fNMCTracks; imc++) fMCTracks[imc].SetNHits( 0 );
749     
750     for( Int_t ih=firstSliceHit; ih<endSliceHit; ih++){
751       Int_t id = fTracker->Hits()[ih].ID();
752       if( id<0 || id>=fNHits ) break;
753       AliHLTTPCCAHitLabel &l = fHitLabels[id];
754       if( l.fLab[0]>=0 ) fMCTracks[l.fLab[0]].SetNHits(fMCTracks[l.fLab[0]].NHits()+1);
755       if( l.fLab[1]>=0 ) fMCTracks[l.fLab[1]].SetNHits(fMCTracks[l.fLab[1]].NHits()+1);
756       if( l.fLab[2]>=0 ) fMCTracks[l.fLab[2]].SetNHits(fMCTracks[l.fLab[2]].NHits()+1);
757     }
758     
759     for (Int_t imc=0; imc<fNMCTracks; imc++) {          
760       AliHLTTPCCAMCTrack &mc = fMCTracks[imc];
761       mc.SetSet( 0 );
762       mc.SetNReconstructed( 0 );
763       mc.SetNTurns( 1 );
764       if( mc.NHits() >=  30 && mc.P()>=.05 ){
765         mc.SetSet( 1 );
766         nMCAll++;
767         if( mc.NHits() >=  30 && mc.P()>=1. ){
768           mc.SetSet( 2 );
769           nMCRef++;
770         }
771       }
772     }
773   }
774
775
776   Int_t traN = *slice.NTracklets();
777   Int_t *traLabels = 0; 
778   Double_t *traPurity = 0;
779
780   traLabels = new Int_t[traN];
781   traPurity = new Double_t[traN];
782   {
783     for (Int_t itr=0; itr<traN; itr++) {
784       traLabels[itr]=-1;
785       traPurity[itr]= 0;
786       
787       Int_t hits[1600];
788       Int_t nHits=0;
789
790       {      
791         Int_t id = slice.TrackletStartHits()[itr];
792         Int_t iRow = AliHLTTPCCATracker::ID2IRow(id);
793         Int_t ih =  AliHLTTPCCATracker::ID2IHit(id);
794       
795         while( ih>=0 ){
796           const AliHLTTPCCARow &row = slice.Row(iRow);
797           hits[nHits] = firstSliceHit + slice.HitInputIDs()[row.FirstHit()+ih];
798           nHits++;      
799           Short_t *linkUp = ((Short_t*)(slice.RowData() + row.FullOffset())) + row.FullLinkOffset();
800           ih = linkUp[ih];
801           iRow++;
802         }
803       }
804
805       if( nHits<5 ) continue;
806         
807       Int_t lb[1600*3];
808       Int_t nla=0;
809
810       for( Int_t ih=0; ih<nHits; ih++ ){
811         AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[hits[ih]].ID()];   
812         if(l.fLab[0]>=0 ) lb[nla++]= l.fLab[0];
813         if(l.fLab[1]>=0 ) lb[nla++]= l.fLab[1];
814         if(l.fLab[2]>=0 ) lb[nla++]= l.fLab[2];
815       }
816       
817       sort( lb, lb+nla );
818       Int_t labmax = -1, labcur=-1, lmax = 0, lcurr=0;
819       for( Int_t i=0; i<nla; i++ ){
820         if( lb[i]!=labcur ){
821           if( labcur>=0 && lmax<lcurr ){
822             lmax = lcurr;
823             labmax = labcur;
824           }
825           labcur = lb[i];
826           lcurr = 0;
827         }
828         lcurr++;
829       }
830       if( labcur>=0 && lmax<lcurr ){
831         lmax = lcurr;
832         labmax = labcur;
833       }
834       lmax = 0;
835       for( Int_t ih=0; ih<nHits; ih++ ){
836         AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[hits[ih]].ID()];   
837         if( l.fLab[0] == labmax || l.fLab[1] == labmax || l.fLab[2] == labmax 
838             ) lmax++;
839       }
840       traLabels[itr] = labmax;
841       traPurity[itr] = ( (nHits>0) ?double(lmax)/double(nHits) :0 );
842     }
843   }
844
845   nRecTot+= traN;
846   
847   for(Int_t itr=0; itr<traN; itr++){      
848     if( traPurity[itr]<.9 || traLabels[itr]<0 || traLabels[itr]>=fNMCTracks){
849       nGhost++;
850       continue;
851     }
852     
853     AliHLTTPCCAMCTrack &mc = fMCTracks[traLabels[itr]]; 
854     mc.SetNReconstructed( mc.NReconstructed()+1 );
855     if( mc.Set()== 0 ) nRecOut++;
856     else{
857       if( mc.NReconstructed()==1 ) nRecAll++;
858       else if(mc.NReconstructed() > mc.NTurns() ) nClonesAll++;
859       if( mc.Set()==2 ){
860         if( mc.NReconstructed()==1 ) nRecRef++;
861         else if(mc.NReconstructed() > mc.NTurns() ) nClonesRef++;
862       }
863     }      
864   }
865
866   for (Int_t ipart=0; ipart<fNMCTracks; ipart++) {              
867     AliHLTTPCCAMCTrack &mc = fMCTracks[ipart];
868     if( mc.Set()>0 ) fhSeedEffVsP->Fill(mc.P(), ( mc.NReconstructed()>0 ?1 :0));
869   }  
870
871   if( traLabels ) delete[] traLabels;
872   if( traPurity ) delete[] traPurity;
873
874   fStatSeedNRecTot += nRecTot;
875   fStatSeedNRecOut += nRecOut;
876   fStatSeedNGhost  += nGhost;
877   fStatSeedNMCAll  += nMCAll;
878   fStatSeedNRecAll  += nRecAll;
879   fStatSeedNClonesAll  += nClonesAll;
880   fStatSeedNMCRef  += nMCRef;
881   fStatSeedNRecRef  += nRecRef;
882   fStatSeedNClonesRef  += nClonesRef;
883
884   if( nMCAll ==0 ) return;
885
886   if( PrintFlag ){
887     cout<<"Track seed performance for slice "<<iSlice<<" : "<<endl;
888     cout<<" N tracks : "
889         <<nMCAll<<" mc all, "
890         <<nMCRef<<" mc ref, "
891         <<nRecTot<<" rec total, "
892         <<nRecAll<<" rec all, "
893         <<nClonesAll<<" clones all, "
894         <<nRecRef<<" rec ref, "
895         <<nClonesRef<<" clones ref, "
896         <<nRecOut<<" out, "
897         <<nGhost<<" ghost"<<endl;
898   
899     Int_t nRecExtr = nRecAll - nRecRef;
900     Int_t nMCExtr = nMCAll - nMCRef;
901     Int_t nClonesExtr = nClonesAll - nClonesRef;
902   
903     Double_t dRecTot = (nRecTot>0 ) ? nRecTot :1;
904     Double_t dMCAll = (nMCAll>0 ) ? nMCAll :1;
905     Double_t dMCRef = (nMCRef>0 ) ? nMCRef :1;
906     Double_t dMCExtr = (nMCExtr>0 ) ? nMCExtr :1;
907     Double_t dRecAll = (nRecAll+nClonesAll>0 ) ? nRecAll+nClonesAll :1;
908     Double_t dRecRef = (nRecRef+nClonesRef>0 ) ? nRecRef+nClonesRef :1;
909     Double_t dRecExtr = (nRecExtr+nClonesExtr>0 ) ? nRecExtr+nClonesExtr :1;
910     
911     cout<<" EffRef = ";
912     if( nMCRef>0 ) cout<<nRecRef/dMCRef; else cout<<"_";
913     cout<<", CloneRef = ";
914     if( nRecRef >0 ) cout << nClonesRef/dRecRef; else cout<<"_";
915     cout<<endl;
916     cout<<" EffExtra = ";
917     if( nMCExtr>0 ) cout << nRecExtr/dMCExtr; else cout<<"_";
918     cout <<", CloneExtra = ";
919     if( nRecExtr>0 ) cout << nClonesExtr/dRecExtr; else cout<<"_";
920     cout<<endl;
921     cout<<" EffAll = ";
922     if( nMCAll>0 ) cout<<nRecAll/dMCAll; else cout<<"_";
923     cout <<", CloneAll = ";
924     if( nRecAll>0 ) cout << nClonesAll/dRecAll; else cout<<"_";
925     cout <<endl;
926     cout<<" Out = ";
927     if( nRecTot>0 ) cout <<nRecOut/dRecTot; else cout<<"_";
928     cout <<", Ghost = ";
929     if( nRecTot>0 ) cout<<nGhost/dRecTot; else cout<<"_";
930     cout<<endl;
931   }
932 }
933
934
935
936
937 void AliHLTTPCCAPerformance::SliceTrackCandPerformance( Int_t iSlice, Bool_t PrintFlag )
938
939   //* calculate slice tracker performance
940   if( !fTracker ) return;
941   
942   Int_t nRecTot = 0, nGhost=0, nRecOut=0;
943   Int_t nMCAll = 0, nRecAll=0, nClonesAll=0;
944   Int_t nMCRef = 0, nRecRef=0, nClonesRef=0;
945   AliHLTTPCCATracker &slice = fTracker->Slices()[iSlice];
946
947   Int_t firstSliceHit = fTracker->FirstSliceHit()[iSlice];
948   Int_t endSliceHit = fTracker->NHits();
949   if( iSlice<fTracker->NSlices()-1 ) endSliceHit = fTracker->FirstSliceHit()[iSlice+1];    
950
951   // Select reconstructable MC tracks
952
953   {
954     for (Int_t imc=0; imc<fNMCTracks; imc++) fMCTracks[imc].SetNHits( 0 );
955     
956     for( Int_t ih=firstSliceHit; ih<endSliceHit; ih++){
957       Int_t id = fTracker->Hits()[ih].ID();
958       if( id<0 || id>=fNHits ) break;
959       AliHLTTPCCAHitLabel &l = fHitLabels[id];
960       if( l.fLab[0]>=0 ) fMCTracks[l.fLab[0]].SetNHits(fMCTracks[l.fLab[0]].NHits()+1);
961       if( l.fLab[1]>=0 ) fMCTracks[l.fLab[1]].SetNHits(fMCTracks[l.fLab[1]].NHits()+1);
962       if( l.fLab[2]>=0 ) fMCTracks[l.fLab[2]].SetNHits(fMCTracks[l.fLab[2]].NHits()+1);
963     }
964     
965     for (Int_t imc=0; imc<fNMCTracks; imc++) {          
966       AliHLTTPCCAMCTrack &mc = fMCTracks[imc];
967       mc.SetSet( 0 );
968       mc.SetNReconstructed( 0 );
969       mc.SetNTurns( 1 );
970       if( mc.NHits() >=  30 && mc.P()>=.05 ){
971         mc.SetSet( 1 );
972         nMCAll++;
973         if( mc.NHits() >=  30 && mc.P()>=1. ){
974           mc.SetSet( 2 );
975           nMCRef++;
976         }
977       }
978     }
979   }
980
981   Int_t traN = *slice.NTracklets();
982   Int_t *traLabels = 0; 
983   Double_t *traPurity = 0;
984   traLabels = new Int_t[traN];
985   traPurity = new Double_t[traN];
986   {
987     for (Int_t itr=0; itr<traN; itr++) {
988       traLabels[itr]=-1;
989       traPurity[itr]= 0;
990
991       AliHLTTPCCATracklet &t = slice.Tracklets()[itr];
992         
993       Int_t nHits = t.NHits();
994       if( nHits<10 ) continue;
995       Int_t firstRow = t.FirstRow();
996       Int_t lastRow = t.LastRow();
997       nHits=0;
998
999       Int_t lb[1600*3];
1000       Int_t nla=0;
1001
1002       for( Int_t irow=firstRow; irow<=lastRow; irow++ ){
1003         Int_t ih = t.RowHit(irow);
1004         if( ih<0 ) continue;
1005         Int_t index = firstSliceHit + slice.HitInputIDs()[slice.Row(irow).FirstHit()+ih];
1006         AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[index].ID()];      
1007         if(l.fLab[0]>=0 ) lb[nla++]= l.fLab[0];
1008         if(l.fLab[1]>=0 ) lb[nla++]= l.fLab[1];
1009         if(l.fLab[2]>=0 ) lb[nla++]= l.fLab[2];
1010         nHits++;     
1011       }
1012       if( nHits<10 ) continue;      
1013       
1014       sort( lb, lb+nla );
1015       Int_t labmax = -1, labcur=-1, lmax = 0, lcurr=0;
1016       for( Int_t i=0; i<nla; i++ ){
1017         if( lb[i]!=labcur ){
1018           if( labcur>=0 && lmax<lcurr ){
1019             lmax = lcurr;
1020             labmax = labcur;
1021           }
1022           labcur = lb[i];
1023           lcurr = 0;
1024         }
1025         lcurr++;
1026       }
1027       if( labcur>=0 && lmax<lcurr ){
1028         lmax = lcurr;
1029         labmax = labcur;
1030       }
1031       lmax = 0;
1032       for( Int_t irow=firstRow; irow<=lastRow; irow++ ){
1033         Int_t ih = t.RowHit(irow);
1034         if( ih<0 ) continue;
1035         Int_t index = firstSliceHit + slice.HitInputIDs()[slice.Row(irow).FirstHit()+ih];
1036         AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[index].ID()];      
1037         if( l.fLab[0] == labmax || l.fLab[1] == labmax || l.fLab[2] == labmax 
1038             ) lmax++;
1039       }
1040       traLabels[itr] = labmax;
1041       traPurity[itr] = ( (nHits>0) ?double(lmax)/double(nHits) :0 );
1042     }
1043   }
1044   
1045   nRecTot+= traN;
1046   
1047   for(Int_t itr=0; itr<traN; itr++){      
1048     if( traPurity[itr]<.9 || traLabels[itr]<0 || traLabels[itr]>=fNMCTracks){
1049       nGhost++;
1050       continue;
1051     }
1052     
1053     AliHLTTPCCAMCTrack &mc = fMCTracks[traLabels[itr]]; 
1054     mc.SetNReconstructed( mc.NReconstructed()+1);
1055     if( mc.Set()== 0 ) nRecOut++;
1056     else{
1057       if( mc.NReconstructed()==1 ) nRecAll++;
1058       else if(mc.NReconstructed() > mc.NTurns() ) nClonesAll++;
1059       if( mc.Set()==2 ){
1060         if( mc.NReconstructed()==1 ) nRecRef++;
1061         else if(mc.NReconstructed() > mc.NTurns() ) nClonesRef++;
1062       }
1063     }      
1064   }
1065
1066   for (Int_t ipart=0; ipart<fNMCTracks; ipart++) {              
1067     AliHLTTPCCAMCTrack &mc = fMCTracks[ipart];
1068     if( mc.Set()>0 ) fhCandEffVsP->Fill(mc.P(), ( mc.NReconstructed()>0 ?1 :0));
1069   }  
1070
1071   if( traLabels ) delete[] traLabels;
1072   if( traPurity ) delete[] traPurity;
1073
1074   fStatCandNRecTot += nRecTot;
1075   fStatCandNRecOut += nRecOut;
1076   fStatCandNGhost  += nGhost;
1077   fStatCandNMCAll  += nMCAll;
1078   fStatCandNRecAll  += nRecAll;
1079   fStatCandNClonesAll  += nClonesAll;
1080   fStatCandNMCRef  += nMCRef;
1081   fStatCandNRecRef  += nRecRef;
1082   fStatCandNClonesRef  += nClonesRef;
1083
1084   if( nMCAll ==0 ) return;
1085
1086   if( PrintFlag ){
1087     cout<<"Track candidate performance for slice "<<iSlice<<" : "<<endl;
1088     cout<<" N tracks : "
1089         <<nMCAll<<" mc all, "
1090         <<nMCRef<<" mc ref, "
1091         <<nRecTot<<" rec total, "
1092         <<nRecAll<<" rec all, "
1093         <<nClonesAll<<" clones all, "
1094         <<nRecRef<<" rec ref, "
1095         <<nClonesRef<<" clones ref, "
1096         <<nRecOut<<" out, "
1097         <<nGhost<<" ghost"<<endl;
1098   
1099     Int_t nRecExtr = nRecAll - nRecRef;
1100     Int_t nMCExtr = nMCAll - nMCRef;
1101     Int_t nClonesExtr = nClonesAll - nClonesRef;
1102   
1103     Double_t dRecTot = (nRecTot>0 ) ? nRecTot :1;
1104     Double_t dMCAll = (nMCAll>0 ) ? nMCAll :1;
1105     Double_t dMCRef = (nMCRef>0 ) ? nMCRef :1;
1106     Double_t dMCExtr = (nMCExtr>0 ) ? nMCExtr :1;
1107     Double_t dRecAll = (nRecAll+nClonesAll>0 ) ? nRecAll+nClonesAll :1;
1108     Double_t dRecRef = (nRecRef+nClonesRef>0 ) ? nRecRef+nClonesRef :1;
1109     Double_t dRecExtr = (nRecExtr+nClonesExtr>0 ) ? nRecExtr+nClonesExtr :1;
1110     
1111     cout<<" EffRef = ";
1112     if( nMCRef>0 ) cout<<nRecRef/dMCRef; else cout<<"_";
1113     cout<<", CloneRef = ";
1114     if( nRecRef >0 ) cout << nClonesRef/dRecRef; else cout<<"_";
1115     cout<<endl;
1116     cout<<" EffExtra = ";
1117     if( nMCExtr>0 ) cout << nRecExtr/dMCExtr; else cout<<"_";
1118     cout <<", CloneExtra = ";
1119     if( nRecExtr>0 ) cout << nClonesExtr/dRecExtr; else cout<<"_";
1120     cout<<endl;
1121     cout<<" EffAll = ";
1122     if( nMCAll>0 ) cout<<nRecAll/dMCAll; else cout<<"_";
1123     cout <<", CloneAll = ";
1124     if( nRecAll>0 ) cout << nClonesAll/dRecAll; else cout<<"_";
1125     cout <<endl;
1126     cout<<" Out = ";
1127     if( nRecTot>0 ) cout <<nRecOut/dRecTot; else cout<<"_";
1128     cout <<", Ghost = ";
1129     if( nRecTot>0 ) cout<<nGhost/dRecTot; else cout<<"_";
1130     cout<<endl;
1131   }
1132 }
1133
1134
1135
1136 void AliHLTTPCCAPerformance::SlicePerformance( Int_t iSlice, Bool_t PrintFlag )
1137
1138   //* calculate slice tracker performance
1139   if( !fTracker ) return;
1140   
1141   Int_t nRecTot = 0, nGhost=0, nRecOut=0;
1142   Int_t nMCAll = 0, nRecAll=0, nClonesAll=0;
1143   Int_t nMCRef = 0, nRecRef=0, nClonesRef=0;
1144   AliHLTTPCCATracker &slice = fTracker->Slices()[iSlice];
1145
1146   Int_t firstSliceHit = 0;
1147   for( ; firstSliceHit<fTracker->NHits(); firstSliceHit++){
1148     if( fTracker->Hits()[firstSliceHit].ISlice()==iSlice ) break;
1149   }
1150   Int_t endSliceHit = firstSliceHit;  
1151
1152   for( ; endSliceHit<fTracker->NHits(); endSliceHit++){
1153     if( fTracker->Hits()[endSliceHit].ISlice()!=iSlice ) break;
1154   }
1155     
1156
1157   // Select reconstructable MC tracks
1158
1159   {
1160     for (Int_t imc=0; imc<fNMCTracks; imc++) fMCTracks[imc].SetNHits( 0 );
1161         
1162     for( Int_t ih=firstSliceHit; ih<endSliceHit; ih++){
1163       Int_t id = fTracker->Hits()[ih].ID();
1164       if( id<0 || id>fNHits ) break;
1165       AliHLTTPCCAHitLabel &l = fHitLabels[id];
1166       if( l.fLab[0]>=0 ) fMCTracks[l.fLab[0]].SetNHits(fMCTracks[l.fLab[0]].NHits()+1);
1167       if( l.fLab[1]>=0 ) fMCTracks[l.fLab[1]].SetNHits(fMCTracks[l.fLab[1]].NHits()+1);
1168       if( l.fLab[2]>=0 ) fMCTracks[l.fLab[2]].SetNHits(fMCTracks[l.fLab[2]].NHits()+1);
1169     }
1170     
1171     for (Int_t imc=0; imc<fNMCTracks; imc++) {          
1172       AliHLTTPCCAMCTrack &mc = fMCTracks[imc];
1173       mc.SetSet( 0 );
1174       mc.SetNReconstructed( 0 );
1175       mc.SetNTurns( 1 );
1176       if( mc.NHits() >=  30 && mc.P()>=.05 ){
1177         mc.SetSet(1);
1178         nMCAll++;
1179         if( mc.NHits() >=  30 && mc.P()>=1. ){
1180           mc.SetSet( 2 );
1181           nMCRef++;
1182         }
1183       }
1184     }
1185   }
1186
1187   Int_t traN = *slice.NOutTracks();
1188   Int_t *traLabels = 0; 
1189   Double_t *traPurity = 0;
1190   traLabels = new Int_t[traN];
1191   traPurity = new Double_t[traN];
1192   {
1193     for (Int_t itr=0; itr<traN; itr++) {
1194       traLabels[itr]=-1;
1195       traPurity[itr]= 0;
1196       AliHLTTPCCAOutTrack &tCA = slice.OutTracks()[itr];
1197       Int_t nhits = tCA.NHits();
1198       Int_t *lb = new Int_t[nhits*3];
1199       Int_t nla=0;
1200       //cout<<"\nHit labels:"<<endl;
1201       for( Int_t ihit=0; ihit<nhits; ihit++){
1202         Int_t index = firstSliceHit + slice.OutTrackHits()[tCA.FirstHitRef()+ihit];
1203         AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[index].ID()];
1204         //cout<<l.fLab[0]<<" "<<l.fLab[1]<<" "<<l.fLab[2]<<endl;
1205         if(l.fLab[0]>=0 ) lb[nla++]= l.fLab[0];
1206         if(l.fLab[1]>=0 ) lb[nla++]= l.fLab[1];
1207         if(l.fLab[2]>=0 ) lb[nla++]= l.fLab[2];
1208       }
1209       sort( lb, lb+nla );
1210       Int_t labmax = -1, labcur=-1, lmax = 0, lcurr=0;
1211       for( Int_t i=0; i<nla; i++ ){
1212         if( lb[i]!=labcur ){
1213           if( labcur>=0 && lmax<lcurr ){
1214             lmax = lcurr;
1215             labmax = labcur;
1216           }
1217           labcur = lb[i];
1218           lcurr = 0;
1219         }
1220         lcurr++;
1221       }
1222       if( labcur>=0 && lmax<lcurr ){
1223         lmax = lcurr;
1224         labmax = labcur;
1225       }
1226       lmax = 0;
1227       for( Int_t ihit=0; ihit<nhits; ihit++){
1228         Int_t index = firstSliceHit + slice.OutTrackHits()[tCA.FirstHitRef()+ihit];
1229         AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[index].ID()];
1230         if( l.fLab[0] == labmax || l.fLab[1] == labmax || l.fLab[2] == labmax 
1231             ) lmax++;
1232       }
1233       traLabels[itr] = labmax;
1234       traPurity[itr] = ( (nhits>0) ?double(lmax)/double(nhits) :0 );
1235       //cout<<"perf track "<<itr<<": "<<nhits<<" "<<labmax<<" "<<traPurity[itr]<<endl;
1236       if( lb ) delete[] lb;
1237     }
1238   }
1239
1240   nRecTot+= traN;
1241
1242   for(Int_t itr=0; itr<traN; itr++){      
1243     if( traPurity[itr]<.9 || traLabels[itr]<0 || traLabels[itr]>=fNMCTracks){
1244       nGhost++;
1245       continue;
1246     }
1247
1248     AliHLTTPCCAMCTrack &mc = fMCTracks[traLabels[itr]]; 
1249     mc.SetNReconstructed(mc.NReconstructed()+1);
1250     if( mc.Set()== 0 ) nRecOut++;
1251     else{
1252       if( mc.NReconstructed()==1 ) nRecAll++;
1253       else if(mc.NReconstructed() > mc.NTurns() ) nClonesAll++;
1254       if( mc.Set()==2 ){
1255         if( mc.NReconstructed()==1 ) nRecRef++;
1256         else if(mc.NReconstructed() > mc.NTurns() ) nClonesRef++;
1257       }
1258     }      
1259   }
1260
1261   for (Int_t ipart=0; ipart<fNMCTracks; ipart++) {              
1262     AliHLTTPCCAMCTrack &mc = fMCTracks[ipart];
1263     if( mc.Set()>0 ) fhEffVsP->Fill(mc.P(), ( mc.NReconstructed()>0 ?1 :0));
1264   }  
1265
1266   if( traLabels ) delete[] traLabels;
1267   if( traPurity ) delete[] traPurity;
1268
1269   fStatNRecTot += nRecTot;
1270   fStatNRecOut += nRecOut;
1271   fStatNGhost  += nGhost;
1272   fStatNMCAll  += nMCAll;
1273   fStatNRecAll  += nRecAll;
1274   fStatNClonesAll  += nClonesAll;
1275   fStatNMCRef  += nMCRef;
1276   fStatNRecRef  += nRecRef;
1277   fStatNClonesRef  += nClonesRef;
1278
1279   if( nMCAll ==0 ) return;
1280
1281   if( PrintFlag ){
1282     cout<<"Performance for slice "<<iSlice<<" : "<<endl;
1283     cout<<" N tracks : "
1284         <<nMCAll<<" mc all, "
1285         <<nMCRef<<" mc ref, "
1286         <<nRecTot<<" rec total, "
1287         <<nRecAll<<" rec all, "
1288         <<nClonesAll<<" clones all, "
1289         <<nRecRef<<" rec ref, "
1290         <<nClonesRef<<" clones ref, "
1291         <<nRecOut<<" out, "
1292         <<nGhost<<" ghost"<<endl;
1293   
1294     Int_t nRecExtr = nRecAll - nRecRef;
1295     Int_t nMCExtr = nMCAll - nMCRef;
1296     Int_t nClonesExtr = nClonesAll - nClonesRef;
1297   
1298     Double_t dRecTot = (nRecTot>0 ) ? nRecTot :1;
1299     Double_t dMCAll = (nMCAll>0 ) ? nMCAll :1;
1300     Double_t dMCRef = (nMCRef>0 ) ? nMCRef :1;
1301     Double_t dMCExtr = (nMCExtr>0 ) ? nMCExtr :1;
1302     Double_t dRecAll = (nRecAll+nClonesAll>0 ) ? nRecAll+nClonesAll :1;
1303     Double_t dRecRef = (nRecRef+nClonesRef>0 ) ? nRecRef+nClonesRef :1;
1304     Double_t dRecExtr = (nRecExtr+nClonesExtr>0 ) ? nRecExtr+nClonesExtr :1;
1305     
1306     cout<<" EffRef = ";
1307     if( nMCRef>0 ) cout<<nRecRef/dMCRef; else cout<<"_";
1308     cout<<", CloneRef = ";
1309     if( nRecRef >0 ) cout << nClonesRef/dRecRef; else cout<<"_";
1310     cout<<endl;
1311     cout<<" EffExtra = ";
1312     if( nMCExtr>0 ) cout << nRecExtr/dMCExtr; else cout<<"_";
1313     cout <<", CloneExtra = ";
1314     if( nRecExtr>0 ) cout << nClonesExtr/dRecExtr; else cout<<"_";
1315     cout<<endl;
1316     cout<<" EffAll = ";
1317     if( nMCAll>0 ) cout<<nRecAll/dMCAll; else cout<<"_";
1318     cout <<", CloneAll = ";
1319     if( nRecAll>0 ) cout << nClonesAll/dRecAll; else cout<<"_";
1320     cout <<endl;
1321     cout<<" Out = ";
1322     if( nRecTot>0 ) cout <<nRecOut/dRecTot; else cout<<"_";
1323     cout <<", Ghost = ";
1324     if( nRecTot>0 ) cout<<nGhost/dRecTot; else cout<<"_";
1325     cout<<endl;
1326   }
1327 }
1328
1329
1330 void AliHLTTPCCAPerformance::Performance( fstream *StatFile )
1331
1332   // main routine for performance calculation  
1333   //SG!!!
1334   /*  
1335   fStatNEvents=0;
1336     fStatNRecTot=0;
1337     fStatNRecOut=0;
1338     fStatNGhost=0;
1339     fStatNMCAll=0;
1340     fStatNRecAll=0;
1341     fStatNClonesAll=0;
1342     fStatNMCRef=0;
1343     fStatNRecRef=0;
1344     fStatNClonesRef=0;
1345   */
1346   fStatNEvents++;  
1347   for( Int_t islice=0; islice<fTracker->NSlices(); islice++){ 
1348     //SliceTrackletPerformance( islice,0 );
1349     SliceTrackCandPerformance( islice, 0 );
1350     SlicePerformance(islice,0);
1351   }
1352
1353
1354   // global tracker performance
1355   {
1356       if( !fTracker ) return;
1357
1358       Int_t nRecTot = 0, nGhost=0, nRecOut=0;
1359       Int_t nMCAll = 0, nRecAll=0, nClonesAll=0;
1360       Int_t nMCRef = 0, nRecRef=0, nClonesRef=0;
1361
1362       // Select reconstructable MC tracks
1363    
1364       {
1365         for (Int_t imc=0; imc<fNMCTracks; imc++) fMCTracks[imc].SetNHits( 0 );
1366           
1367         for( Int_t ih=0; ih<fNHits; ih++){
1368           AliHLTTPCCAHitLabel &l = fHitLabels[ih];
1369           if( l.fLab[0]>=0 ) fMCTracks[l.fLab[0]].SetNHits(fMCTracks[l.fLab[0]].NHits()+1);
1370           if( l.fLab[1]>=0 ) fMCTracks[l.fLab[1]].SetNHits(fMCTracks[l.fLab[1]].NHits()+1);
1371           if( l.fLab[2]>=0 ) fMCTracks[l.fLab[2]].SetNHits(fMCTracks[l.fLab[2]].NHits()+1);
1372         }
1373     
1374         for (Int_t imc=0; imc<fNMCTracks; imc++) {              
1375           AliHLTTPCCAMCTrack &mc = fMCTracks[imc];
1376           mc.SetSet( 0 );
1377           mc.SetNReconstructed( 0 );
1378           mc.SetNTurns( 1 );
1379           if( mc.NHits() >=  50 && mc.P()>=.05 ){
1380             mc.SetSet(1);
1381             nMCAll++;
1382             if( mc.P()>=1. ){
1383               mc.SetSet(2);
1384               nMCRef++;
1385             }
1386           }
1387         }
1388       }
1389
1390       Int_t traN = fTracker->NTracks();
1391       Int_t *traLabels = 0;
1392       Double_t *traPurity = 0;
1393       traLabels = new Int_t[traN];
1394       traPurity = new Double_t[traN];
1395       {
1396         for (Int_t itr=0; itr<traN; itr++) {
1397           traLabels[itr]=-1;
1398           traPurity[itr]= 0;
1399           AliHLTTPCCAGBTrack &tCA = fTracker->Tracks()[itr];
1400           Int_t nhits = tCA.NHits();
1401           Int_t *lb = new Int_t[nhits*3];
1402           Int_t nla=0;
1403           for( Int_t ihit=0; ihit<nhits; ihit++){
1404             Int_t index = fTracker->TrackHits()[tCA.FirstHitRef()+ihit];
1405             AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[index].ID()];
1406             if(l.fLab[0]>=0 ) lb[nla++]= l.fLab[0];
1407             if(l.fLab[1]>=0 ) lb[nla++]= l.fLab[1];
1408             if(l.fLab[2]>=0 ) lb[nla++]= l.fLab[2];
1409           }
1410           sort( lb, lb+nla );
1411           Int_t labmax = -1, labcur=-1, lmax = 0, lcurr=0;
1412           for( Int_t i=0; i<nla; i++ ){
1413             if( lb[i]!=labcur ){
1414               if( labcur>=0 && lmax<lcurr ){
1415                 lmax = lcurr;
1416                 labmax = labcur;
1417               }
1418               labcur = lb[i];
1419               lcurr = 0;
1420             }
1421             lcurr++;
1422           }
1423           if( labcur>=0 && lmax<lcurr ){
1424             lmax = lcurr;
1425             labmax = labcur;
1426           }
1427           lmax = 0;
1428           for( Int_t ihit=0; ihit<nhits; ihit++){
1429             Int_t index = fTracker->TrackHits()[tCA.FirstHitRef()+ihit];
1430             AliHLTTPCCAHitLabel &l = fHitLabels[fTracker->Hits()[index].ID()];
1431             if( l.fLab[0] == labmax || l.fLab[1] == labmax || l.fLab[2] == labmax 
1432                 ) lmax++;
1433           }
1434           traLabels[itr] = labmax;
1435           traPurity[itr] = ( (nhits>0) ?double(lmax)/double(nhits) :0 );
1436           if( lb ) delete[] lb;
1437         }
1438       }
1439
1440       nRecTot+= traN;
1441       for(Int_t itr=0; itr<traN; itr++){      
1442         if( traPurity[itr]<.9 || traLabels[itr]<0 || traLabels[itr]>=fNMCTracks){
1443           nGhost++;
1444           continue;
1445         }
1446         AliHLTTPCCAGBTrack &tCA = fTracker->Tracks()[itr];
1447         AliHLTTPCCAMCTrack &mc = fMCTracks[traLabels[itr]];     
1448         
1449         mc.SetNReconstructed(mc.NReconstructed()+1);
1450         if( mc.Set()== 0 ) nRecOut++;
1451         else{
1452           if( mc.NReconstructed()==1 ) nRecAll++;
1453           else if(mc.NReconstructed() > mc.NTurns() ) nClonesAll++;
1454           if( mc.Set()==2 ){
1455             if( mc.NReconstructed()==1 ) nRecRef++;
1456             else if(mc.NReconstructed() > mc.NTurns() ) nClonesRef++;
1457             fhTrackLengthRef->Fill( tCA.NHits()/((Double_t) mc.NHits()));
1458           }
1459         }      
1460
1461         // track resolutions
1462         while( mc.Set()==2 && TMath::Abs(mc.TPCPar()[0]) + TMath::Abs(mc.TPCPar()[1])>1 ){
1463           if( traPurity[itr]<.90 ) break;
1464           AliHLTTPCCAGBTrack &t = fTracker->Tracks()[itr];
1465           AliHLTTPCCATrackParam p = t.Param();
1466           Double_t cosA = TMath::Cos( t.Alpha() );
1467           Double_t sinA = TMath::Sin( t.Alpha() );
1468           Double_t mcX =  mc.TPCPar()[0]*cosA + mc.TPCPar()[1]*sinA;
1469           Double_t mcY = -mc.TPCPar()[0]*sinA + mc.TPCPar()[1]*cosA;
1470           Double_t mcZ =  mc.TPCPar()[2];
1471           Double_t mcEx =  mc.TPCPar()[3]*cosA + mc.TPCPar()[4]*sinA;
1472           Double_t mcEy = -mc.TPCPar()[3]*sinA + mc.TPCPar()[4]*cosA;
1473           Double_t mcEz =  mc.TPCPar()[5];
1474           Double_t mcEt = TMath::Sqrt(mcEx*mcEx + mcEy*mcEy);
1475           if( TMath::Abs(mcEt)<1.e-4 ) break;
1476           Double_t mcSinPhi = mcEy / mcEt;
1477           Double_t mcDzDs   = mcEz / mcEt;
1478           Double_t mcQPt = mc.TPCPar()[6]/ mcEt;
1479           if( TMath::Abs(mcQPt)<1.e-4 ) break;
1480           Double_t mcPt = 1./TMath::Abs(mcQPt);
1481           if( mcPt<1. ) break;
1482           if( t.NHits() <  50 ) break;
1483           Double_t bz = fTracker->Slices()[0].Param().Bz();
1484           if( !p.TransportToXWithMaterial( mcX, bz ) ) break;
1485           if( p.GetCosPhi()*mcEx < 0 ){ // change direction
1486             mcSinPhi = -mcSinPhi;
1487             mcDzDs = -mcDzDs;
1488             mcQPt = -mcQPt;
1489           }
1490           const Double_t kCLight = 0.000299792458;  
1491           Double_t k2QPt = 100;
1492           if( TMath::Abs(bz)>1.e-4 ) k2QPt= 1./(bz*kCLight);
1493           Double_t qPt = p.GetKappa()*k2QPt;
1494           Double_t pt = 100;
1495           if( TMath::Abs(qPt) >1.e-4 ) pt = 1./TMath::Abs(qPt);
1496           
1497           fhResY->Fill( p.GetY() - mcY ); 
1498           fhResZ->Fill( p.GetZ() - mcZ );
1499           fhResSinPhi->Fill( p.GetSinPhi() - mcSinPhi );
1500           fhResDzDs->Fill( p.GetDzDs() - mcDzDs );
1501           fhResPt->Fill( ( pt - mcPt )/mcPt );
1502
1503           if( p.GetErr2Y()>0 ) fhPullY->Fill( (p.GetY() - mcY)/TMath::Sqrt(p.GetErr2Y()) ); 
1504           if( p.GetErr2Z()>0 ) fhPullZ->Fill( (p.GetZ() - mcZ)/TMath::Sqrt(p.GetErr2Z()) ); 
1505
1506           if( p.GetErr2SinPhi()>0 ) fhPullSinPhi->Fill( (p.GetSinPhi() - mcSinPhi)/TMath::Sqrt(p.GetErr2SinPhi()) ); 
1507           if( p.GetErr2DzDs()>0 ) fhPullDzDs->Fill( (p.DzDs() - mcDzDs)/TMath::Sqrt(p.GetErr2DzDs()) ); 
1508           if( p.GetErr2Kappa()>0 ) fhPullQPt->Fill( (qPt - mcQPt)/TMath::Sqrt(p.GetErr2Kappa()*k2QPt*k2QPt) ); 
1509           fhPullYS->Fill( TMath::Sqrt(fTracker->GetChi2( p.GetY(), p.GetSinPhi(), p.GetCov()[0], p.GetCov()[3], p.GetCov()[5], 
1510                                              mcY, mcSinPhi, 0,0,0 )));
1511           fhPullZT->Fill( TMath::Sqrt(fTracker->GetChi2( p.GetZ(), p.GetDzDs(), p.GetCov()[2], p.GetCov()[7], p.GetCov()[9], 
1512                           mcZ, mcDzDs, 0,0,0 ) ));
1513
1514           break;
1515         }
1516       }
1517
1518       for (Int_t ipart=0; ipart<fNMCTracks; ipart++) {          
1519         AliHLTTPCCAMCTrack &mc = fMCTracks[ipart];
1520         if( mc.Set()>0 ) fhGBEffVsP->Fill(mc.P(), ( mc.NReconstructed()>0 ?1 :0));
1521         if( mc.Set()>0 ) fhGBEffVsPt->Fill(mc.Pt(), ( mc.NReconstructed()>0 ?1 :0));
1522         if( mc.Set()==2 ){ 
1523           const Double_t *p = mc.TPCPar();
1524           Double_t r = TMath::Sqrt(p[0]*p[0] + p[1]*p[1]);
1525           Double_t cosA = p[0]/r;
1526           Double_t sinA = p[1]/r;
1527
1528
1529           Double_t phipos = TMath::Pi()+TMath::ATan2(-p[1], -p[0]);
1530           Double_t alpha =  TMath::Pi()*(20*((((Int_t)(phipos*180/TMath::Pi()))/20))+10)/180.;
1531           cosA = TMath::Cos(alpha);
1532           sinA = TMath::Sin(alpha);
1533
1534           Double_t mcX =  p[0]*cosA + p[1]*sinA;
1535           Double_t mcY = -p[0]*sinA + p[1]*cosA;
1536           Double_t mcZ =  p[2];
1537           Double_t mcEx =  p[3]*cosA + p[4]*sinA;
1538           Double_t mcEy = -p[3]*sinA + p[4]*cosA;
1539           Double_t mcEz =  p[5];
1540           //Double_t mcEt = TMath::Sqrt(mcEx*mcEx + mcEy*mcEy);
1541           Double_t angleY = TMath::ATan2(mcEy, mcEx)*180./TMath::Pi();
1542           Double_t angleZ = TMath::ATan2(mcEz, mcEx)*180./TMath::Pi();
1543
1544           if( mc.NReconstructed()>0 ){
1545             fhRefRecoX->Fill(mcX);
1546             fhRefRecoY->Fill(mcY);
1547             fhRefRecoZ->Fill(mcZ);
1548             fhRefRecoP->Fill(mc.P());
1549             fhRefRecoPt->Fill(mc.Pt());
1550             fhRefRecoAngleY->Fill(angleY);
1551             fhRefRecoAngleZ->Fill(angleZ);
1552             fhRefRecoNHits->Fill(mc.NHits());
1553           } else {
1554             fhRefNotRecoX->Fill(mcX);
1555             fhRefNotRecoY->Fill(mcY);
1556             fhRefNotRecoZ->Fill(mcZ);
1557             fhRefNotRecoP->Fill(mc.P());
1558             fhRefNotRecoPt->Fill(mc.Pt());
1559             fhRefNotRecoAngleY->Fill(angleY);
1560             fhRefNotRecoAngleZ->Fill(angleZ);
1561             fhRefNotRecoNHits->Fill(mc.NHits());
1562           }
1563         }
1564       }
1565
1566       if( traLabels ) delete[] traLabels;
1567       if( traPurity ) delete[] traPurity;
1568
1569       fStatGBNRecTot += nRecTot;
1570       fStatGBNRecOut += nRecOut;
1571       fStatGBNGhost  += nGhost;
1572       fStatGBNMCAll  += nMCAll;
1573       fStatGBNRecAll  += nRecAll;
1574       fStatGBNClonesAll  += nClonesAll;
1575       fStatGBNMCRef  += nMCRef;
1576       fStatGBNRecRef  += nRecRef;
1577       fStatGBNClonesRef  += nClonesRef;
1578   }
1579
1580   
1581   // distribution of cluster errors
1582
1583   {    
1584     Int_t nHits = fTracker->NHits();
1585     for( Int_t ih=0; ih<nHits; ih++ ){
1586       AliHLTTPCCAGBHit &hit = fTracker->Hits()[ih];
1587       AliHLTTPCCAHitLabel &l = fHitLabels[hit.ID()];
1588       fhHitErrY->Fill(hit.ErrY());
1589       fhHitErrZ->Fill(hit.ErrZ());
1590       Int_t nmc=0;
1591       for( Int_t il=0; il<3; il++ ) if( l.fLab[il]>=0 ) nmc++;
1592       if( nmc==1 ) fhHitShared->Fill(hit.IRow(),0);
1593       else if( nmc>1 ) fhHitShared->Fill(hit.IRow(),1);
1594     }
1595   }
1596
1597   // cluster pulls
1598
1599   if( fDoClusterPulls && fNMCPoints>0 ) {
1600
1601     {
1602       for (Int_t ipart=0; ipart<fNMCTracks; ipart++) {          
1603         AliHLTTPCCAMCTrack &mc = fMCTracks[ipart];
1604         mc.SetNMCPoints(0);
1605       }
1606       sort(fMCPoints, fMCPoints+fNMCPoints, AliHLTTPCCAMCPoint::Compare );
1607       
1608       for( Int_t ip=0; ip<fNMCPoints; ip++ ){
1609         AliHLTTPCCAMCPoint &p = fMCPoints[ip];
1610         AliHLTTPCCAMCTrack &t = fMCTracks[p.TrackID()];
1611         if( t.NMCPoints()==0 ) t.SetFirstMCPointID( ip );
1612         t.SetNMCPoints(t.NMCPoints()+1);
1613       }
1614     }
1615
1616     for( Int_t ih=0; ih<fNHits; ih++ ){
1617
1618       AliHLTTPCCAGBHit &hit = fTracker->Hits()[ih];
1619       AliHLTTPCCAHitLabel &l = fHitLabels[ih];
1620
1621       if( l.fLab[0]<0 || l.fLab[0]>=fNMCTracks
1622           || l.fLab[1]>=0 || l.fLab[2]>=0       ) continue;
1623
1624       Int_t lab = l.fLab[0];
1625
1626       AliHLTTPCCAMCTrack &track = fMCTracks[lab];
1627       //if( track.Pt()<1. ) continue;
1628       Int_t ip1=-1, ip2=-1;
1629       Double_t d1 = 1.e20, d2=1.e20;
1630       for( Int_t ip=0; ip<track.NMCPoints(); ip++ ){
1631         AliHLTTPCCAMCPoint &p = fMCPoints[track.FirstMCPointID() + ip];
1632         if( p.ISlice() != hit.ISlice() ) continue;        
1633         Double_t dx = p.Sx()-hit.X();
1634         Double_t dy = p.Sy()-hit.Y();
1635         Double_t dz = p.Sz()-hit.Z();
1636         Double_t d = dx*dx + dy*dy + dz*dz;
1637         if( p.Sx()< hit.X() ){
1638           if( d<d1 ){
1639             d1 = d;
1640             ip1 = ip;
1641           }
1642         }else{
1643           if( d<d2 ){
1644             d2 = d;
1645             ip2 = ip;
1646           }
1647         }
1648       }
1649
1650       if( ip1<0 || ip2<0 ) continue;
1651
1652       AliHLTTPCCAMCPoint &p1 = fMCPoints[track.FirstMCPointID() + ip1];
1653       AliHLTTPCCAMCPoint &p2 = fMCPoints[track.FirstMCPointID() + ip2];
1654       Double_t dx = p2.Sx() - p1.Sx();
1655       Double_t dy = p2.Sy() - p1.Sy();
1656       Double_t dz = p2.Sz() - p1.Sz();
1657       if( TMath::Abs(dx)>1.e-8 && TMath::Abs(p1.Sx()-hit.X())<2. && TMath::Abs(p2.Sx()-hit.X())<2.  ){
1658         Double_t sx = hit.X();
1659         Double_t sy = p1.Sy() + dy/dx*(sx-p1.Sx());
1660         Double_t sz = p1.Sz() + dz/dx*(sx-p1.Sx());
1661         
1662         Float_t errY, errZ;
1663         {
1664           AliHLTTPCCATrackParam t;
1665           t.SetZ( sz );
1666           t.SetSinPhi( dy/TMath::Sqrt(dx*dx+dy*dy) );
1667           t.SetCosPhi( dx/TMath::Sqrt(dx*dx+dy*dy) );
1668           t.SetDzDs( dz/TMath::Sqrt(dx*dx+dy*dy) );
1669           fTracker->GetErrors2(hit,t,errY, errZ );
1670           errY = TMath::Sqrt(errY);
1671           errZ = TMath::Sqrt(errZ);
1672         } 
1673              
1674         fhHitResY->Fill((hit.Y()-sy));
1675         fhHitResZ->Fill((hit.Z()-sz));
1676         fhHitPullY->Fill((hit.Y()-sy)/errY);
1677         fhHitPullZ->Fill((hit.Z()-sz)/errZ);
1678         if( track.Pt()>=1. ){
1679           fhHitResY1->Fill((hit.Y()-sy));
1680           fhHitResZ1->Fill((hit.Z()-sz));
1681           fhHitPullY1->Fill((hit.Y()-sy)/errY);
1682           fhHitPullZ1->Fill((hit.Z()-sz)/errZ);
1683         }
1684       }
1685     }   
1686   }
1687
1688   {
1689     cout<<"\nSlice Track Seed performance: \n"<<endl;
1690     cout<<" N tracks : "
1691         <<fStatNMCAll/fStatNEvents<<" mc all, "
1692         <<fStatSeedNMCRef/fStatNEvents<<" mc ref, "
1693         <<fStatSeedNRecTot/fStatNEvents<<" rec total, "
1694         <<fStatSeedNRecAll/fStatNEvents<<" rec all, "
1695         <<fStatSeedNClonesAll/fStatNEvents<<" clones all, "
1696         <<fStatSeedNRecRef/fStatNEvents<<" rec ref, "
1697         <<fStatSeedNClonesRef/fStatNEvents<<" clones ref, "
1698         <<fStatSeedNRecOut/fStatNEvents<<" out, "
1699         <<fStatSeedNGhost/fStatNEvents<<" ghost"<<endl;
1700   
1701     Int_t nRecExtr = fStatSeedNRecAll - fStatSeedNRecRef;
1702     Int_t nMCExtr = fStatNMCAll - fStatNMCRef;
1703     Int_t nClonesExtr = fStatSeedNClonesAll - fStatSeedNClonesRef;
1704     
1705     Double_t dRecTot = (fStatSeedNRecTot>0 ) ? fStatSeedNRecTot :1;
1706     Double_t dMCAll = (fStatNMCAll>0 ) ? fStatNMCAll :1;
1707     Double_t dMCRef = (fStatNMCRef>0 ) ? fStatNMCRef :1;
1708     Double_t dMCExtr = (nMCExtr>0 ) ? nMCExtr :1;
1709     Double_t dRecAll = (fStatSeedNRecAll+fStatSeedNClonesAll>0 ) ? fStatSeedNRecAll+fStatSeedNClonesAll :1;
1710     Double_t dRecRef = (fStatSeedNRecRef+fStatSeedNClonesRef>0 ) ? fStatSeedNRecRef+fStatSeedNClonesRef :1;
1711     Double_t dRecExtr = (nRecExtr+nClonesExtr>0 ) ? nRecExtr+nClonesExtr :1;
1712     
1713     cout<<" EffRef = "<< fStatSeedNRecRef/dMCRef
1714         <<", CloneRef = " << fStatSeedNClonesRef/dRecRef <<endl;
1715     cout<<" EffExtra = "<<nRecExtr/dMCExtr
1716         <<", CloneExtra = " << nClonesExtr/dRecExtr<<endl;
1717     cout<<" EffAll = "<<fStatSeedNRecAll/dMCAll
1718         <<", CloneAll = " << fStatSeedNClonesAll/dRecAll<<endl;
1719     cout<<" Out = "<<fStatSeedNRecOut/dRecTot
1720         <<", Ghost = "<<fStatSeedNGhost/dRecTot<<endl;
1721   }
1722
1723   {
1724     cout<<"\nSlice Track candidate performance: \n"<<endl;
1725     cout<<" N tracks : "
1726         <<fStatNMCAll/fStatNEvents<<" mc all, "
1727         <<fStatCandNMCRef/fStatNEvents<<" mc ref, "
1728         <<fStatCandNRecTot/fStatNEvents<<" rec total, "
1729         <<fStatCandNRecAll/fStatNEvents<<" rec all, "
1730         <<fStatCandNClonesAll/fStatNEvents<<" clones all, "
1731         <<fStatCandNRecRef/fStatNEvents<<" rec ref, "
1732         <<fStatCandNClonesRef/fStatNEvents<<" clones ref, "
1733         <<fStatCandNRecOut/fStatNEvents<<" out, "
1734         <<fStatCandNGhost/fStatNEvents<<" ghost"<<endl;
1735   
1736     Int_t nRecExtr = fStatCandNRecAll - fStatCandNRecRef;
1737     Int_t nMCExtr = fStatNMCAll - fStatNMCRef;
1738     Int_t nClonesExtr = fStatCandNClonesAll - fStatCandNClonesRef;
1739     
1740     Double_t dRecTot = (fStatCandNRecTot>0 ) ? fStatCandNRecTot :1;
1741     Double_t dMCAll = (fStatNMCAll>0 ) ? fStatNMCAll :1;
1742     Double_t dMCRef = (fStatNMCRef>0 ) ? fStatNMCRef :1;
1743     Double_t dMCExtr = (nMCExtr>0 ) ? nMCExtr :1;
1744     Double_t dRecAll = (fStatCandNRecAll+fStatCandNClonesAll>0 ) ? fStatCandNRecAll+fStatCandNClonesAll :1;
1745     Double_t dRecRef = (fStatCandNRecRef+fStatCandNClonesRef>0 ) ? fStatCandNRecRef+fStatCandNClonesRef :1;
1746     Double_t dRecExtr = (nRecExtr+nClonesExtr>0 ) ? nRecExtr+nClonesExtr :1;
1747     
1748     cout<<" EffRef = "<< fStatCandNRecRef/dMCRef
1749         <<", CloneRef = " << fStatCandNClonesRef/dRecRef <<endl;
1750     cout<<" EffExtra = "<<nRecExtr/dMCExtr
1751         <<", CloneExtra = " << nClonesExtr/dRecExtr<<endl;
1752     cout<<" EffAll = "<<fStatCandNRecAll/dMCAll
1753         <<", CloneAll = " << fStatCandNClonesAll/dRecAll<<endl;
1754     cout<<" Out = "<<fStatCandNRecOut/dRecTot
1755         <<", Ghost = "<<fStatCandNGhost/dRecTot<<endl;
1756   }
1757
1758   {
1759     cout<<"\nSlice tracker performance: \n"<<endl;
1760     cout<<" N tracks : "
1761         <<fStatNMCAll/fStatNEvents<<" mc all, "
1762         <<fStatNMCRef/fStatNEvents<<" mc ref, "
1763         <<fStatNRecTot/fStatNEvents<<" rec total, "
1764         <<fStatNRecAll/fStatNEvents<<" rec all, "
1765         <<fStatNClonesAll/fStatNEvents<<" clones all, "
1766         <<fStatNRecRef/fStatNEvents<<" rec ref, "
1767         <<fStatNClonesRef/fStatNEvents<<" clones ref, "
1768         <<fStatNRecOut/fStatNEvents<<" out, "
1769         <<fStatNGhost/fStatNEvents<<" ghost"<<endl;
1770   
1771     Int_t nRecExtr = fStatNRecAll - fStatNRecRef;
1772     Int_t nMCExtr = fStatNMCAll - fStatNMCRef;
1773     Int_t nClonesExtr = fStatNClonesAll - fStatNClonesRef;
1774     
1775     Double_t dRecTot = (fStatNRecTot>0 ) ? fStatNRecTot :1;
1776     Double_t dMCAll = (fStatNMCAll>0 ) ? fStatNMCAll :1;
1777     Double_t dMCRef = (fStatNMCRef>0 ) ? fStatNMCRef :1;
1778     Double_t dMCExtr = (nMCExtr>0 ) ? nMCExtr :1;
1779     Double_t dRecAll = (fStatNRecAll+fStatNClonesAll>0 ) ? fStatNRecAll+fStatNClonesAll :1;
1780     Double_t dRecRef = (fStatNRecRef+fStatNClonesRef>0 ) ? fStatNRecRef+fStatNClonesRef :1;
1781     Double_t dRecExtr = (nRecExtr+nClonesExtr>0 ) ? nRecExtr+nClonesExtr :1;
1782     
1783     cout<<" EffRef = "<< fStatNRecRef/dMCRef
1784         <<", CloneRef = " << fStatNClonesRef/dRecRef <<endl;
1785     cout<<" EffExtra = "<<nRecExtr/dMCExtr
1786         <<", CloneExtra = " << nClonesExtr/dRecExtr<<endl;
1787     cout<<" EffAll = "<<fStatNRecAll/dMCAll
1788         <<", CloneAll = " << fStatNClonesAll/dRecAll<<endl;
1789     cout<<" Out = "<<fStatNRecOut/dRecTot
1790         <<", Ghost = "<<fStatNGhost/dRecTot<<endl;
1791     cout<<" Time = "<<fTracker->StatTime(0)/fTracker->StatNEvents()*1.e3<<" msec/event "<<endl;
1792     cout<<" Local timers = "
1793         <<fTracker->StatTime(1)/fTracker->StatNEvents()*1.e3<<" "
1794         <<fTracker->StatTime(2)/fTracker->StatNEvents()*1.e3<<" "
1795         <<fTracker->StatTime(3)/fTracker->StatNEvents()*1.e3<<" "
1796         <<fTracker->StatTime(4)/fTracker->StatNEvents()*1.e3<<" "
1797         <<fTracker->StatTime(5)/fTracker->StatNEvents()*1.e3<<" "
1798         <<fTracker->StatTime(6)/fTracker->StatNEvents()*1.e3<<" "
1799         <<fTracker->StatTime(7)/fTracker->StatNEvents()*1.e3<<" "
1800         <<fTracker->StatTime(8)/fTracker->StatNEvents()*1.e3<<" "
1801         <<" msec/event "<<endl;
1802   }
1803
1804
1805   {
1806     cout<<"\nGlobal tracker performance for "<<fStatNEvents<<" events: \n"<<endl;
1807     cout<<" N tracks : "
1808         <<fStatGBNMCAll<<" mc all, "
1809         <<fStatGBNMCRef<<" mc ref, "
1810         <<fStatGBNRecTot<<" rec total, "
1811         <<fStatGBNRecAll<<" rec all, "
1812         <<fStatGBNClonesAll<<" clones all, "
1813         <<fStatGBNRecRef<<" rec ref, "
1814         <<fStatGBNClonesRef<<" clones ref, "
1815         <<fStatGBNRecOut<<" out, "
1816         <<fStatGBNGhost<<" ghost"<<endl;
1817      cout<<" N tracks average : "
1818         <<fStatGBNMCAll/fStatNEvents<<" mc all, "
1819         <<fStatGBNMCRef/fStatNEvents<<" mc ref, "
1820         <<fStatGBNRecTot/fStatNEvents<<" rec total, "
1821         <<fStatGBNRecAll/fStatNEvents<<" rec all, "
1822         <<fStatGBNClonesAll/fStatNEvents<<" clones all, "
1823         <<fStatGBNRecRef/fStatNEvents<<" rec ref, "
1824         <<fStatGBNClonesRef/fStatNEvents<<" clones ref, "
1825         <<fStatGBNRecOut/fStatNEvents<<" out, "
1826         <<fStatGBNGhost/fStatNEvents<<" ghost"<<endl;
1827   
1828     Int_t nRecExtr = fStatGBNRecAll - fStatGBNRecRef;
1829     Int_t nMCExtr = fStatGBNMCAll - fStatGBNMCRef;
1830     Int_t nClonesExtr = fStatGBNClonesAll - fStatGBNClonesRef;
1831     
1832     Double_t dRecTot = (fStatGBNRecTot>0 ) ? fStatGBNRecTot :1;
1833     Double_t dMCAll = (fStatGBNMCAll>0 ) ? fStatGBNMCAll :1;
1834     Double_t dMCRef = (fStatGBNMCRef>0 ) ? fStatGBNMCRef :1;
1835     Double_t dMCExtr = (nMCExtr>0 ) ? nMCExtr :1;
1836     Double_t dRecAll = (fStatGBNRecAll+fStatGBNClonesAll>0 ) ? fStatGBNRecAll+fStatGBNClonesAll :1;
1837     Double_t dRecRef = (fStatGBNRecRef+fStatGBNClonesRef>0 ) ? fStatGBNRecRef+fStatGBNClonesRef :1;
1838     Double_t dRecExtr = (nRecExtr+nClonesExtr>0 ) ? nRecExtr+nClonesExtr :1;
1839     
1840     cout<<" EffRef = "<< fStatGBNRecRef/dMCRef
1841         <<", CloneRef = " << fStatGBNClonesRef/dRecRef <<endl;
1842     cout<<" EffExtra = "<<nRecExtr/dMCExtr
1843         <<", CloneExtra = " << nClonesExtr/dRecExtr<<endl;
1844     cout<<" EffAll = "<<fStatGBNRecAll/dMCAll
1845         <<", CloneAll = " << fStatGBNClonesAll/dRecAll<<endl;
1846     cout<<" Out = "<<fStatGBNRecOut/dRecTot
1847         <<", Ghost = "<<fStatGBNGhost/dRecTot<<endl;
1848     cout<<" Time = "<<( fTracker->StatTime(0)+fTracker->StatTime(9) )/fTracker->StatNEvents()*1.e3<<" msec/event "<<endl;
1849     cout<<" Local timers: "<<endl;
1850     cout<<" slice tracker "<<fTracker->StatTime(0)/fTracker->StatNEvents()*1.e3<<": "
1851         <<fTracker->StatTime(1)/fTracker->StatNEvents()*1.e3<<" "
1852         <<fTracker->StatTime(2)/fTracker->StatNEvents()*1.e3<<" "
1853         <<fTracker->StatTime(3)/fTracker->StatNEvents()*1.e3<<" "
1854         <<fTracker->StatTime(4)/fTracker->StatNEvents()*1.e3<<" "
1855         <<fTracker->StatTime(5)/fTracker->StatNEvents()*1.e3<<"["
1856         <<fTracker->StatTime(6)/fTracker->StatNEvents()*1.e3<<"/"
1857         <<fTracker->StatTime(7)/fTracker->StatNEvents()*1.e3<<"] "
1858         <<fTracker->StatTime(8)/fTracker->StatNEvents()*1.e3
1859         <<" msec/event "<<endl;
1860     cout<<" GB merger "<<fTracker->StatTime(9)/fTracker->StatNEvents()*1.e3<<": "
1861         <<fTracker->StatTime(10)/fTracker->StatNEvents()*1.e3<<", "
1862         <<fTracker->StatTime(11)/fTracker->StatNEvents()*1.e3<<", "
1863         <<fTracker->StatTime(12)/fTracker->StatNEvents()*1.e3<<" "
1864         <<" msec/event "<<endl;
1865
1866     if( StatFile && StatFile->is_open() ){
1867       fstream &out = *StatFile;
1868
1869       //out<<"\nGlobal tracker performance for "<<fStatNEvents<<" events: \n"<<endl;
1870       //out<<" N tracks : "
1871       //<<fStatGBNMCAll/fStatNEvents<<" mc all, "
1872       //<<fStatGBNMCRef/fStatNEvents<<" mc ref, "
1873       // <<fStatGBNRecTot/fStatNEvents<<" rec total, "
1874       // <<fStatGBNRecAll/fStatNEvents<<" rec all, "
1875       // <<fStatGBNClonesAll/fStatNEvents<<" clones all, "
1876       // <<fStatGBNRecRef/fStatNEvents<<" rec ref, "
1877       // <<fStatGBNClonesRef/fStatNEvents<<" clones ref, "
1878       // <<fStatGBNRecOut/fStatNEvents<<" out, "
1879       // <<fStatGBNGhost/fStatNEvents<<" ghost"<<endl;
1880       fStatTime+=fTracker->SliceTrackerTime();
1881       double timeHz=0;
1882       if( fStatTime>1.e-4 ) timeHz = 1./fStatTime*fStatNEvents;
1883
1884       out<<"<table border>"<<endl;
1885       out<<"<tr>"<<endl; 
1886       out<<"<td>      </td> <td align=center> RefSet </td> <td align=center> AllSet </td> <td align=center> ExtraSet </td>"<<endl;
1887       out<<"</tr>"<<endl;
1888       out<<"<tr>"<<endl;
1889       out<<"<td>Efficiency</td> <td align=center>"<<fStatGBNRecRef/dMCRef
1890          <<"</td> <td align=center>"<<fStatGBNRecAll/dMCAll
1891          <<"</td> <td align=center>"<<nRecExtr/dMCExtr
1892          <<"</td>"<<endl;
1893       out<<"</tr>"<<endl;
1894       out<<"<tr> "<<endl;
1895       out<<"<td>Clone</td>      <td align=center>"<<fStatGBNClonesRef/dRecRef
1896          <<"</td> <td align=center>"<<fStatGBNClonesAll/dRecAll
1897          <<"</td> <td align=center>"<<nClonesExtr/dRecExtr
1898          <<"</td>"<<endl;
1899       out<<"</tr>"<<endl;
1900       out<<"<tr> "<<endl;
1901       out<<"<td>Ghost</td>      <td colspan=3 align=center>"<<fStatGBNGhost/dRecTot
1902          <<"</td>"<<endl;
1903       out<<"</tr>"<<endl;
1904       out<<"<tr> "<<endl;
1905       out<<"<td>Time</td>      <td colspan=3 align=center>"<<timeHz 
1906          <<" ev/s</td>"<<endl;
1907       out<<"</tr>"<<endl;
1908       out<<"<tr> "<<endl;
1909       out<<"<td>N Events</td>      <td colspan=3 align=center>"<<fStatNEvents
1910          <<"</td>"<<endl;
1911       out<<"</tr>"<<endl;
1912       out<<"</table>"<<endl;
1913     }
1914
1915   }
1916
1917   WriteHistos();
1918
1919
1920 }
1921
1922 void AliHLTTPCCAPerformance::WriteMCEvent( ostream &out ) const
1923 {
1924   // write MC information to the file
1925   out<<fNMCTracks<<endl;
1926   for( Int_t it=0; it<fNMCTracks; it++ ){
1927     AliHLTTPCCAMCTrack &t = fMCTracks[it];
1928     out<<it<<" ";
1929     out<<t.PDG()<<endl;
1930     for( Int_t i=0; i<7; i++ ) out<<t.Par()[i]<<" ";
1931     out<<endl<<"    ";
1932     for( Int_t i=0; i<7; i++ ) out<<t.TPCPar()[i]<<" ";
1933     out<<endl<<"    ";
1934     out<< t.P()<<" ";
1935     out<< t.Pt()<<" ";
1936     out<< t.NMCPoints()<<" ";
1937     out<< t.FirstMCPointID()<<" ";
1938     out<< t.NHits()<<" ";
1939     out<< t.NReconstructed()<<" ";
1940     out<< t.Set()<<" ";
1941     out<< t.NTurns()<<endl;
1942   }
1943
1944   out<<fNHits<<endl;
1945   for( Int_t ih=0; ih<fNHits; ih++ ){
1946     AliHLTTPCCAHitLabel &l = fHitLabels[ih];
1947     out<<l.fLab[0]<<" "<<l.fLab[1]<<" "<<l.fLab[2]<<endl;
1948   }
1949 }
1950
1951 void AliHLTTPCCAPerformance::WriteMCPoints( ostream &out ) const
1952 {  
1953   // write Mc points to the file
1954   out<<fNMCPoints<<endl;
1955   for( Int_t ip=0; ip<fNMCPoints; ip++ ){
1956     AliHLTTPCCAMCPoint &p = fMCPoints[ip];
1957     out<< p.X()<<" ";
1958     out<< p.Y()<<" ";
1959     out<< p.Z()<<" ";
1960     out<< p.Sx()<<" ";
1961     out<< p.Sy()<<" ";
1962     out<< p.Sz()<<" ";
1963     out<< p.Time()<<" ";
1964     out<< p.ISlice()<<" ";
1965     out<< p.TrackID()<<endl;
1966   }
1967 }
1968
1969 void AliHLTTPCCAPerformance::ReadMCEvent( istream &in )
1970 {
1971   // read mc info from the file
1972   StartEvent();
1973   if( fMCTracks ) delete[] fMCTracks;
1974   fMCTracks = 0;
1975   fNMCTracks = 0;
1976   if( fHitLabels ) delete[] fHitLabels;
1977   fHitLabels = 0;
1978   fNHits = 0;
1979   if( fMCPoints ) delete[] fMCPoints;
1980   fMCPoints = 0;
1981   fNMCPoints = 0;
1982
1983   in>>fNMCTracks;
1984   fMCTracks = new AliHLTTPCCAMCTrack[fNMCTracks];
1985   for( Int_t it=0; it<fNMCTracks; it++ ){
1986     AliHLTTPCCAMCTrack &t = fMCTracks[it];
1987     Int_t j;
1988     Float_t f;
1989     in>>j;
1990     in>> j; t.SetPDG(j);
1991     for( Int_t i=0; i<7; i++ ){ in>>f; t.SetPar(i,f);}
1992     for( Int_t i=0; i<7; i++ ){ in>>f; t.SetTPCPar(i,f);}
1993     in>> f; t.SetP(f);
1994     in>> f; t.SetPt(f);
1995     in>> j; t.SetNHits(j);
1996     in>> j; t.SetNMCPoints(j);
1997     in>> j; t.SetFirstMCPointID(j);
1998     in>> j; t.SetNReconstructed(j);
1999     in>> j; t.SetSet(j);
2000     in>> j; t.SetNTurns(j);
2001   }
2002   
2003   in>>fNHits;
2004   fHitLabels = new AliHLTTPCCAHitLabel[fNHits];
2005   for( Int_t ih=0; ih<fNHits; ih++ ){
2006     AliHLTTPCCAHitLabel &l = fHitLabels[ih];
2007     in>>l.fLab[0]>>l.fLab[1]>>l.fLab[2];
2008   }
2009 }
2010
2011 void AliHLTTPCCAPerformance::ReadMCPoints( istream &in )
2012 {
2013   // read mc points from the file
2014   if( fMCPoints ) delete[] fMCPoints;
2015   fMCPoints = 0;
2016   fNMCPoints = 0;
2017   
2018   in>>fNMCPoints;
2019   fMCPoints = new AliHLTTPCCAMCPoint[fNMCPoints];
2020   for( Int_t ip=0; ip<fNMCPoints; ip++ ){
2021     AliHLTTPCCAMCPoint &p = fMCPoints[ip];
2022     Float_t f;
2023     Int_t i;
2024     in>> f;
2025     p.SetX( f );
2026     in >> f;
2027     p.SetY( f );
2028     in >> f;
2029     p.SetZ( f );
2030     in >> f;
2031     p.SetSx( f );
2032     in >> f;
2033     p.SetSy( f );
2034     in >> f;
2035     p.SetSz( f );
2036     in >> f;
2037     p.SetTime( f );
2038     in >> i;
2039     p.SetISlice( i );
2040     in >> i;
2041     p.SetTrackID( i );
2042   }
2043 }