]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/comp/AliHLTTPCCompModelAnalysis.cxx
cleaning memory on error condition and loop break
[u/mrichter/AliRoot.git] / HLT / TPCLib / comp / AliHLTTPCCompModelAnalysis.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: J. Wagner <jwagner@cern.ch>                           *
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 /** @file   AliHLTTPCCompModelAnalysis.cxx
20     @author J. Wagner jwagner@cern.ch
21     @date   17-11-2007
22     @brief  A processing analysis component for the HLT */
23
24 #if __GNUC__ >= 3
25 using namespace std;
26 #endif
27
28 #include "AliHLTTPCCompModelAnalysis.h"
29 #include "AliHLTTPCTransform.h"
30 #include "AliHLTTPCModelTrack.h"
31 #include "AliHLTTPCCompDataCompressorHelper.h"
32 #include "TFile.h"
33 #include "TH1.h"
34 #include <cerrno>
35
36 /** constructor **/
37 AliHLTTPCCompModelAnalysis::AliHLTTPCCompModelAnalysis(Bool_t modelanalysis, Bool_t trackanalysis, TString dumpfilename, TString graphfilename):
38   fModelAnalysis(modelanalysis),
39   fTrackAnalysis(trackanalysis),
40   fDumpFileName(dumpfilename),
41   fGraphFileName(graphfilename),
42   fFirstTrackArray(),
43   fSecondTrackArray(),
44   fFirstTrackList(NULL),
45   fSecondTrackList(NULL),
46   fFirstTrashTracks(0),
47   fSecondTrashTracks(0),
48   fTotalComparedTracks(0),
49   fMatchedFirstTrashTracks(0),
50   fMatchedSecondTrashTracks(0),
51   fFirstUnmatchedTracks(0),
52   fSecondUnmatchedTracks(0),
53   fToleranceDeviation(0.0),
54   fTrackListPointer(NULL),
55   fTotalDiscardedClusters(0),
56   fValuableDiscardedClusters(0),
57   fTrashTracks(0)
58 {
59   // see header file for class documentation
60   // or
61   // refer to README to build package
62   // or
63   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
64 }
65
66 /** destructor **/
67 AliHLTTPCCompModelAnalysis::~AliHLTTPCCompModelAnalysis()
68     {
69   // see header file for class documentation
70     for ( UInt_t slice=0; slice<36; slice++ )
71         for ( UInt_t patch=0; patch<6; patch++ )
72             {
73             if ( fDiscardedClusters[slice][patch]!=NULL )
74                 {
75                 delete [] fDiscardedClusters[slice][patch];
76                 fDiscardedClusters[slice][patch]=NULL;
77                 }
78             }
79     }
80
81 /** initialise arrays for tracks/ discarded clusters depending on model-flags **/
82 Int_t AliHLTTPCCompModelAnalysis::Init()
83 {
84   // see header file for class documentation
85
86   if(fTrackAnalysis) // track quantities are to be initialised
87     {
88       fFirstTrackArray.Reset();
89       fSecondTrackArray.Reset();
90
91       fFirstTrackList = NULL;
92       fSecondTrackList = NULL;
93
94       fFirstTrashTracks = 0;
95       fSecondTrashTracks = 0;
96
97       fTotalComparedTracks = 0;
98       fMatchedFirstTrashTracks = 0;
99       fMatchedSecondTrashTracks = 0;
100
101       fFirstUnmatchedTracks = 0;
102       fSecondUnmatchedTracks = 0;
103
104       fToleranceDeviation = 0.001;
105
106     }
107   
108   if(fModelAnalysis) // cluster array to be initialised
109     {
110       for ( UInt_t slice=0; slice<36; slice++ )
111         {
112           for ( UInt_t patch=0; patch<6; patch++ )
113             {
114               fDiscardedClusters[slice][patch] = NULL;
115             }
116         }
117
118       // initialise trash track list to store discarded tracks
119       fTrackListPointer = NULL;
120
121       // set all counters to zero:
122       fTrashTracks = 0;
123       fTotalDiscardedClusters = 0;
124       fValuableDiscardedClusters = 0;
125     }
126
127   return 0;
128 }
129
130 Int_t AliHLTTPCCompModelAnalysis::DisplayResults()
131 {
132   // see header file for class documentation
133   HLTInfo("--------------------DISPLAYING RESULTS---------------------");
134   // if model loss analysis, then display these results, else display track results
135   if (fModelAnalysis)
136     {
137       DisplayModelResults();
138     };
139
140   if (fTrackAnalysis)
141     {
142       DisplayTrackResults();
143     };
144
145   // error message: if no analysis flag is switched on
146   if( (!fModelAnalysis) && (!fTrackAnalysis) )
147     {
148       HLTError("Error! Display Results called without any analysis flag switched on.");
149       return 1;
150     }
151
152   return 0;
153 }
154
155 Int_t AliHLTTPCCompModelAnalysis::SetTracks(AliHLTTPCTrackletData* tracklets, Bool_t fillingfirsttracks)
156 {
157   // see header file for class documentation
158   // if fFillingFirstTrackArray is true (i.e. first input file in component is processed)
159   // first input track array is filled
160   // else the second track array is filled
161
162   if(fillingfirsttracks)
163     {
164       HLTDebug( "Reading %u tracks in first array", (unsigned) tracklets->fTrackletCnt );
165       
166       if(tracklets->fTrackletCnt == 0)
167         {
168           HLTError("Error! No tracklets to fill into first track array!");
169           return EINVAL;
170         }
171       else
172         {
173           fFirstTrackArray.FillTracks(tracklets->fTrackletCnt, tracklets->fTracklets );
174         }
175
176     }
177   else
178     {
179       if(tracklets->fTrackletCnt == 0)
180         {
181           HLTError("Error! No tracklets to fill into second track array!");
182           return EINVAL;
183         }
184       else
185         {
186           fSecondTrackArray.FillTracks(tracklets->fTrackletCnt, tracklets->fTracklets );
187         }
188       // read in tracks in second array (tracks from clusters after Vestbo in second tracking process)
189       HLTDebug( "Reading %u tracks in second array", (unsigned) tracklets->fTrackletCnt );
190     
191     }
192
193   return 0;
194 }
195
196 Int_t AliHLTTPCCompModelAnalysis::SetClusters(AliHLTTPCClusterData* clusters, UInt_t slice, UInt_t patch, Bool_t fillingfirstclusters)
197 {
198   // see header file for class documentation
199   if(fillingfirstclusters == 1)
200     {
201       if ( slice>=36 || patch>=6 )
202         return EINVAL;
203       if ( fOriginalClusters[slice][patch] )
204         return EBUSY;
205       fOriginalClusters[slice][patch] = clusters;
206
207       //HLTDebug( "Filling %u clusters in first array", (unsigned)clusters->fSpacePointCnt);
208     }
209   else
210     {
211       if ( slice>=36 || patch>=6 )
212         return EINVAL;
213       if ( fSecondaryClusters[slice][patch] )
214         return EBUSY;
215       fSecondaryClusters[slice][patch] = clusters;
216
217       //HLTDebug( "Filling %u clusters in second array", (unsigned)clusters->fSpacePointCnt);
218     }
219       return 0;
220 }
221
222 Int_t AliHLTTPCCompModelAnalysis::CompareTracks()
223 {
224   // see header file for class documentation
225   // define variables for number of tracks in track arrays:
226   Int_t firsttracks = fFirstTrackArray.GetNTracks();
227   Int_t secondtracks = fSecondTrackArray.GetNTracks();
228
229   // error checking: if second array has been filled or not:
230   if(firsttracks == 0)
231     {
232       HLTError("No tracks in first track array!");
233       return -EINVAL;
234     };
235  
236
237   if(secondtracks == 0)
238     {
239       HLTError("No tracks in second track array!");
240       return -EINVAL;
241     };
242
243   // take track from first tracking,
244   for(Int_t ii=0; ii < firsttracks; ii++)
245     {
246       // build track list for all tracks in first array
247       // FIXME: I can't find the cleanup of the linked list fFirstTrackList
248       AliHLTTPCTrackList* currenttrackentry = new AliHLTTPCTrackList;
249       if (!currenttrackentry) return -ENOMEM;
250       currenttrackentry->fTrack = *(fFirstTrackArray.GetCheckedTrack(ii));
251
252       // get its pythia information, 
253       currenttrackentry->fPythiatrack = GetComparableTrackPythiaInfo(currenttrackentry->fTrack);
254
255       currenttrackentry->fMatchingindicator = 0;
256
257       // put this element as first in list
258       currenttrackentry->fNext = fFirstTrackList;
259       fFirstTrackList = currenttrackentry;
260
261       // count tracks below 0.1GeV
262       if(currenttrackentry->fTrack.GetPt()<0.1)
263         {
264           ++fFirstTrashTracks;
265         }
266
267     }
268
269  // take track from second tracking,
270   for(Int_t ii=0; ii < secondtracks; ii++)
271     {
272       // build track list for all tracks in second array
273       // FIXME: I can't find the cleanup of the linked list fSecondTrackArray
274       AliHLTTPCTrackList* currenttrackentry = new AliHLTTPCTrackList;
275       if (!currenttrackentry) return -ENOMEM;
276       currenttrackentry->fTrack = *(fSecondTrackArray.GetCheckedTrack(ii));
277
278       // get its pythia information, 
279       currenttrackentry->fPythiatrack = GetComparableTrackPythiaInfo(currenttrackentry->fTrack);
280
281       // put this element as first in list
282       currenttrackentry->fNext = fSecondTrackList;
283       fSecondTrackList = currenttrackentry;
284
285       // count tracks below 0.1GeV
286       if(currenttrackentry->fTrack.GetPt()<0.1)
287         {
288           ++fSecondTrashTracks;
289         }
290
291     }
292
293   // search for matching track from secondary tracking
294   AliHLTTPCTrackList* firstmatchpointer = fFirstTrackList;
295  
296   while(firstmatchpointer != NULL)
297    {
298      AliHLTTPCTrackList* secondmatchpointer = fSecondTrackList;
299
300      while(secondmatchpointer != NULL)
301        {
302         
303          // compare paramters of the two tracks, 
304          // match only when coincidence >= 50% in fToleranceDeviation range!
305          if ((CompareTrackInfo(firstmatchpointer, secondmatchpointer) > 4))
306            { 
307     
308              if((CompareTrackInfo(firstmatchpointer, secondmatchpointer) > firstmatchpointer->fMatchingindicator))
309                {
310                  // look if current better matching track has already been matched before
311                  if((CompareTrackInfo(firstmatchpointer, secondmatchpointer) > secondmatchpointer->fMatchingindicator))
312                    {
313                      
314                      // set previously assigned matchingindicator (if there was one) of secondary track back to zero
315                      if(firstmatchpointer->fMatchingindicator > 0)
316                        {
317                          firstmatchpointer->fMatchingtrack->fMatchingindicator = 0;
318                          firstmatchpointer->fMatchingtrack->fMatchingtrack = NULL;
319                        }
320
321                      if(secondmatchpointer->fMatchingindicator > 0)
322                        {
323                          secondmatchpointer->fMatchingtrack->fMatchingindicator = 0;
324                          secondmatchpointer->fMatchingtrack->fMatchingtrack = NULL;
325                        }
326                      
327                      // compare according to tracks themselves (other possibility: compare pythiatracks - better!)
328                      secondmatchpointer->fMatchingindicator = CompareTrackInfo(firstmatchpointer, secondmatchpointer) ;
329                      firstmatchpointer->fMatchingindicator = CompareTrackInfo(firstmatchpointer, secondmatchpointer);
330                      
331                      // remember which track matches which
332                      secondmatchpointer->fMatchingtrack = firstmatchpointer;
333                      firstmatchpointer->fMatchingtrack = secondmatchpointer;
334
335                    } // end if compare > second matching indicator
336
337                } // end if compare > first matching indicator
338              
339            }// end if compare > 4
340                  
341          secondmatchpointer = secondmatchpointer->fNext;
342        }
343      
344      // go on with next original track
345      firstmatchpointer = firstmatchpointer->fNext;
346    }
347   
348   // count not matched tracks in first and second track list
349   AliHLTTPCTrackList* nomatchcounter = fFirstTrackList;
350
351   while(nomatchcounter != NULL)
352     {
353       if(nomatchcounter->fMatchingindicator == 0)
354         {
355           ++fFirstUnmatchedTracks;
356         }
357       else
358         {
359          ++fTotalComparedTracks;
360
361          // count matched trash tracks
362          if(nomatchcounter->fTrack.GetPt() < 0.1)
363            {
364              ++fMatchedFirstTrashTracks;
365            };
366         }
367       nomatchcounter = nomatchcounter->fNext;
368     }
369
370   nomatchcounter = fSecondTrackList;
371   while(nomatchcounter != NULL)
372     {
373       if(nomatchcounter->fMatchingindicator == 0)
374         {
375           ++fSecondUnmatchedTracks;
376         }
377       else
378         {
379           // count matched trash tracks
380          if(nomatchcounter->fTrack.GetPt() < 0.1)
381            {
382              ++fMatchedSecondTrashTracks;
383            };
384         }
385     
386       nomatchcounter = nomatchcounter->fNext;
387     }
388  
389   // consistency check: fFirstUnmatchedTracks + fTotalComparedTracks = # of tracks in first array
390   // ...and analogously for second array:
391   if(fFirstUnmatchedTracks + fTotalComparedTracks !=  firsttracks)
392     {
393       HLTWarning("Warning! Possible inconsistency in original track array: Number of compared and unmatched tracks not equal to total number of tracks!");
394     };
395   
396   if(fSecondUnmatchedTracks + fTotalComparedTracks !=  secondtracks)
397     {
398       HLTWarning("Warning! Possible inconsistency in second track array: Number of compared and unmatched tracks not equal to total number of tracks!");
399     };
400   
401   return 0;
402 }
403
404 Int_t AliHLTTPCCompModelAnalysis::CompareClusters(Bool_t relativedifferences)
405 {
406
407   Int_t totaloriginal = 0;
408   Int_t totalsecondary = 0;
409   Int_t usedclusters = 0;
410   Int_t comparedclusters = 0;
411   Int_t notcomparedclusters = 0;
412
413   // create graphs out of differences and leave loop
414   TFile* clustergraphrootfile = NULL;
415   if(!fGraphFileName.IsNull())
416     {
417       clustergraphrootfile = new TFile(fGraphFileName, "recreate");
418     }
419
420   // specifications of histograms
421   Double_t clusterdifffxmin, clusterdifffxmax;
422   Double_t clusterdifffymin, clusterdifffymax;
423   Double_t clusterdifffzmin, clusterdifffzmax;
424   Int_t clusterdifffxbins, clusterdifffybins, clusterdifffzbins;
425
426   Double_t clusterdifffsigmay2min, clusterdifffsigmay2max;
427   Double_t clusterdifffsigmaz2min, clusterdifffsigmaz2max;
428   Int_t clusterdifffsigmay2bins, clusterdifffsigmaz2bins;
429
430   if (!relativedifferences) // not tested yet!
431     {
432       clusterdifffxmin = -1;
433       clusterdifffxmax = +1;
434       clusterdifffxbins = (Int_t) ((clusterdifffxmax - clusterdifffxmin)/0.0001);
435   
436       clusterdifffymin = -1;
437       clusterdifffymax = +1;
438       clusterdifffybins = (Int_t) ((clusterdifffymax - clusterdifffymin)/0.0001);
439   
440       clusterdifffzmin = -1;
441       clusterdifffzmax = +1;
442       clusterdifffzbins = (Int_t) ((clusterdifffzmax - clusterdifffzmin)/0.0001);
443   
444       clusterdifffsigmay2min = -1;
445       clusterdifffsigmay2max = +1;
446       clusterdifffsigmay2bins = (Int_t) ((clusterdifffsigmay2max - clusterdifffsigmay2min)/0.0001);
447   
448       clusterdifffsigmaz2min = -1;
449       clusterdifffsigmaz2max = +1;
450       clusterdifffsigmaz2bins = (Int_t) ((clusterdifffsigmaz2max - clusterdifffsigmaz2min)/0.0001);
451     }
452   else
453     {
454       clusterdifffxmin = -1;
455       clusterdifffxmax = +1;
456       clusterdifffxbins = (Int_t) ((clusterdifffxmax - clusterdifffxmin)/0.0001);
457   
458       clusterdifffymin = -1;
459       clusterdifffymax = +1;
460       clusterdifffybins = (Int_t) ((clusterdifffymax - clusterdifffymin)/0.0001);
461   
462       clusterdifffzmin = -1;
463       clusterdifffzmax = +1;
464       clusterdifffzbins = (Int_t) ((clusterdifffzmax - clusterdifffzmin)/0.0001);
465   
466       clusterdifffsigmay2min = -1;
467       clusterdifffsigmay2max = +1;
468       clusterdifffsigmay2bins = (Int_t) ((clusterdifffsigmay2max - clusterdifffsigmay2min)/0.0001);
469   
470       clusterdifffsigmaz2min = -1;
471       clusterdifffsigmaz2max = +1;
472       clusterdifffsigmaz2bins = (Int_t) ((clusterdifffsigmaz2max - clusterdifffsigmaz2min)/0.0001);
473     }
474   
475   // intialise histogramms
476   TH1F* clusterfxhisto = new TH1F("Differences of x (original - secondary) clusters", "Differences of x (original - secondary) clusters", clusterdifffxbins, clusterdifffxmin, clusterdifffxmax);
477   TH1F* clusterfyhisto = new TH1F("Differences of y (original - secondary) clusters", "Differences of y (original - secondary) clusters", clusterdifffybins, clusterdifffymin, clusterdifffymax);
478   TH1F* clusterfzhisto = new TH1F("Differences of z (original - secondary) clusters", "Differences of z (original - secondary) clusters", clusterdifffzbins, clusterdifffzmin, clusterdifffzmax);
479   TH1F* clusterfsigmay2histo = new TH1F("Differences of sigmay2 (original - secondary) clusters", "Differences of sigmay2 (original - secondary) clusters", clusterdifffsigmay2bins, clusterdifffsigmay2min, clusterdifffsigmay2max);
480   TH1F* clusterfsigmaz2histo = new TH1F("Differences of sigmaz2 (original - secondary) clusters", "Differences of sigmaz2 (original - secondary) clusters", clusterdifffsigmaz2bins, clusterdifffsigmaz2min, clusterdifffsigmaz2max);
481   
482
483   // see headerfile for class documentation
484   // compare for each slice and patch the clusters in the original cluster array
485   // to the ones of the secondary cluster array
486   for(Int_t slicecntr = 0; slicecntr < 36; slicecntr++)
487     {
488       for (Int_t patchcntr = 0; patchcntr < 6; patchcntr++)
489         {
490           if(!fOriginalClusters[slicecntr][patchcntr])
491             {
492               // HLTDebug("No original clusters for slice %d patch %d", slicecntr, patchcntr);
493               continue;
494             }
495
496           if(!fSecondaryClusters[slicecntr][patchcntr])
497             {
498               //HLTDebug("No secondary clusters for slice %d patch %d", slicecntr, patchcntr);
499               continue;
500             }
501
502           for ( unsigned long ii=0; ii<fOriginalClusters[slicecntr][patchcntr]->fSpacePointCnt; ii++ )
503             {
504               ++totaloriginal;
505
506               // search matching secondary cluster by charge and padrow, 
507               // fill histograms if cluster has been used in tracking process
508               Int_t found = 0;
509              
510               for(unsigned long jj=0; jj<fSecondaryClusters[slicecntr][patchcntr]->fSpacePointCnt; jj++)
511                 {
512                   // if fTrackN != -1 -> draw histograms out of used clusters
513                   // if fTrackN == -1 -> draw histograms out of unused clusters
514                   if (fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fTrackN == -1)
515                     {
516                       if((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fCharge == fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fCharge) && (fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fPadRow == fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fPadRow) )
517                         {
518                           
519                           if (relativedifferences == 1)
520                             {
521                               // check whether first entries in cluster array are zero
522                               if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX == 0)
523                                 {
524                                   HLTWarning("Warning! x value of original cluster is zero, relative differences cannot be calculated!");
525                                 };
526
527                               if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY == 0)
528                                 {
529                                   HLTWarning("Warning! y value of original cluster is zero, relative differences cannot be calculated!");
530                                 };
531
532                               if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ == 0)
533                                 {
534                                   HLTWarning("Warning! z value of original cluster is zero, relative differences cannot be calculated!");
535                                 };
536
537                               if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2 == 0)
538                                 {
539                                   HLTWarning("Warning! sigmay2 value of original cluster is zero, relative differences cannot be calculated!");
540                                 };
541
542                               if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2 == 0)
543                                 {
544                                   HLTWarning("Warning! sigmaz2 value of original cluster is zero, relative differences cannot be calculated!");
545                                 };
546
547                               // fill relative differences in histograms
548                               clusterfxhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fX)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX,1);
549                               clusterfyhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fY)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY,1);
550                               clusterfzhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fZ)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ,1);
551                               clusterfsigmay2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaY2)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2,1);
552                               clusterfsigmaz2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaZ2)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2,1);
553                             }
554                           else
555                             {
556                               // fill absolute differences histograms
557                               clusterfxhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fX),1);
558                               clusterfyhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fY),1);
559                               clusterfzhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fZ),1);
560                               clusterfsigmay2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaY2),1);
561                               clusterfsigmaz2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaZ2),1);
562                             }
563                               
564                           found = 1;
565                           ++comparedclusters;
566                           break;
567                         }
568                     }
569                 }
570               
571               if(found == 0)
572                 {
573                   ++notcomparedclusters;
574                 }
575
576             }
577         }
578     }
579
580   // write graphs to rootfile
581  if(!fGraphFileName.IsNull())
582     {
583       clustergraphrootfile->WriteObject(clusterfxhisto, "clusterfxhistogram");
584       clustergraphrootfile->WriteObject(clusterfyhisto, "clusterfyhistogram");
585       clustergraphrootfile->WriteObject(clusterfzhisto, "clusterfzhistogram");
586       clustergraphrootfile->WriteObject(clusterfsigmay2histo, "clusterfsigmay2histogram");
587       clustergraphrootfile->WriteObject(clusterfsigmaz2histo, "clusterfsigmaz2histogram");
588       clustergraphrootfile->Close();
589     }
590   
591   // count clusters used for tracking
592   for (Int_t slicecount=0; slicecount<36; slicecount++)
593     {
594       for(Int_t patchcount=0; patchcount<6; patchcount++)
595         {
596
597           if(!fSecondaryClusters[slicecount][patchcount])
598             {
599               //HLTDebug("No secondary clusters for slice %d patch %d", slicecntr, patchcntr);
600               continue;
601             }
602
603           for(Int_t count=0; count < (Int_t) fSecondaryClusters[slicecount][patchcount]->fSpacePointCnt; count++)
604             {
605
606               ++totalsecondary;
607
608               if(fSecondaryClusters[slicecount][patchcount]->fSpacePoints[count].fTrackN != -1)
609                 {
610                   ++usedclusters;
611                 };
612             }
613         }
614     }
615   
616   // Display results of cluster analysis
617   HLTInfo("Number of original clusters: %d", totaloriginal);
618   HLTInfo("Number of 2ndary clusters: %d", totalsecondary);
619   HLTInfo("Number of 2ndary clusters used for tracking: %d", usedclusters);
620   HLTInfo("Number of compared (tracked) original clusters: %d", comparedclusters);
621   HLTInfo("Number of uncompared (tracked) original clusters: %d", notcomparedclusters);
622
623   //////////////////////////////////////////////////////
624   FILE* clusterfile = fopen("/afsuser/jwagner/TrackerTest_25092007/cosmics/fullanalysis/clusteranalysis.out", "a");
625
626   fprintf(clusterfile, "%d \t %d \t %d \t %d \t %d \t %d \n", totaloriginal, totalsecondary, usedclusters, comparedclusters, notcomparedclusters, totaloriginal-comparedclusters-notcomparedclusters);
627
628   fclose(clusterfile);
629   ////////////////////////////////////////////////////////
630
631   // consistency check
632   if(comparedclusters + notcomparedclusters != totaloriginal)
633     {
634       HLTWarning("Possible inconsistency: number of compared clusters %d + number of not compared clusters %d not equal to total number of original clusters %d", comparedclusters, notcomparedclusters, totaloriginal);
635     }
636
637   return 0;
638 }
639
640 AliHLTTPCTrack AliHLTTPCCompModelAnalysis::GetComparableTrackPythiaInfo(const AliHLTTPCTrack& comparabletrack) const
641 {
642   // see headerfile for class documentation
643
644   AliHLTTPCTrack pythiatrack = comparabletrack;
645
646   return pythiatrack;
647 }
648
649 Int_t AliHLTTPCCompModelAnalysis::MarkTrashTrack(AliHLTTPCTrack *lowpttrack)
650 {
651   // see header file for class documentation
652
653   // save track first in lowpttrack list (all lowpttracks are displayed in display function altogether)
654   AliHLTTPCTrackList* tracklistentry =  new AliHLTTPCTrackList;
655   tracklistentry->fTrack = *lowpttrack;
656   tracklistentry->fWronglydiscarded = GetTrashTrackPythiaInfo(lowpttrack);
657   tracklistentry->fMatchingindicator = 0; // not needed here, therefore initialised to zero
658   tracklistentry->fNext = fTrackListPointer;
659   tracklistentry->fMatchingtrack = NULL; // not needed here, therefore initialised to NULL
660   fTrackListPointer = tracklistentry;
661
662   ++fTrashTracks;
663
664   return 0;
665 }
666
667 Int_t AliHLTTPCCompModelAnalysis::MarkTrashCluster(AliHLTTPCClusterData *discardedcluster, UInt_t slice, UInt_t patch)
668 {
669   // see header file for class documentation
670
671   // get Pythia information of discarded cluster
672   Bool_t wronglydiscarded = GetClusterPythiaInfo(discardedcluster);
673
674   // if cluster has been discarded wrongly, save information
675   if(wronglydiscarded)
676     {
677       // store cluster
678       fDiscardedClusters[slice][patch] = discardedcluster;
679       // increase number of valuable discarded clusters
680       ++fValuableDiscardedClusters;
681     }
682       
683       ++fTotalDiscardedClusters;   
684       
685   return 0;
686 }
687
688 Bool_t AliHLTTPCCompModelAnalysis::GetTrashTrackPythiaInfo(const AliHLTTPCTrack* /*discardedtrack*/ ) const
689 {
690   // see header file for class documentation
691   // store information from pythia in current track list entry
692   // fTrackListPointer.pythiatrack = FillFromPythia...
693
694   return 0;
695 }
696
697 Bool_t AliHLTTPCCompModelAnalysis::GetClusterPythiaInfo(const AliHLTTPCClusterData* /*discardedcluster*/) const
698 {
699   // see header file for class documentation
700   // Pythia information can be
701   // either: cluster belongs to discarded track with pt < 0.1 Gev (--> cluster correctly discarded)
702   //     or: cluster is discarded and does not belong to any pythia track (--> correctly discarded)
703   //     or: cluster is discarded but belongs to pythia track (--> cluster WRONGLY discarded!!!)
704
705   return 0;
706 }
707  
708 Int_t AliHLTTPCCompModelAnalysis::CompareTrackInfo(AliHLTTPCTrackList* firsttracklistelement, AliHLTTPCTrackList* secondtracklistelement)
709 {
710   // see header file for class documentation
711   // calculate matching indicator accoring to the track information
712   // ++matchingindicator for every paramter that matches
713
714   Int_t currentmatchingindicator = 0;
715
716   // tolerance range of 1 percent deviation for each quantity
717
718   // compare start point (x,y,z)
719   if(abs((firsttracklistelement->fTrack.GetFirstPointX() - secondtracklistelement->fTrack.GetFirstPointX()))/firsttracklistelement->fTrack.GetFirstPointX() <= fToleranceDeviation)
720     ++currentmatchingindicator;
721
722   if(abs((firsttracklistelement->fTrack.GetFirstPointY() - secondtracklistelement->fTrack.GetFirstPointY()))/firsttracklistelement->fTrack.GetFirstPointY() <= fToleranceDeviation)
723     ++currentmatchingindicator;
724
725   if(abs((firsttracklistelement->fTrack.GetFirstPointZ() - secondtracklistelement->fTrack.GetFirstPointZ()))/firsttracklistelement->fTrack.GetFirstPointZ() <= fToleranceDeviation)
726     ++currentmatchingindicator;
727   
728   // compare end point
729   if(abs((firsttracklistelement->fTrack.GetLastPointX() - secondtracklistelement->fTrack.GetLastPointX()))/firsttracklistelement->fTrack.GetLastPointX() <= fToleranceDeviation)
730     ++currentmatchingindicator;
731
732   if(abs((firsttracklistelement->fTrack.GetLastPointY() - secondtracklistelement->fTrack.GetLastPointY()))/firsttracklistelement->fTrack.GetLastPointY() <= fToleranceDeviation)
733     ++currentmatchingindicator;
734
735   if(abs((firsttracklistelement->fTrack.GetLastPointZ() - secondtracklistelement->fTrack.GetLastPointZ()))/firsttracklistelement->fTrack.GetLastPointZ() <= fToleranceDeviation)
736     ++currentmatchingindicator;
737
738   // compare pt, psi, tgl
739   if(abs((firsttracklistelement->fTrack.GetPt() - secondtracklistelement->fTrack.GetPt()))/firsttracklistelement->fTrack.GetPt() <= fToleranceDeviation)
740     ++currentmatchingindicator;
741
742   if(abs((firsttracklistelement->fTrack.GetPsi() - secondtracklistelement->fTrack.GetPsi()))/firsttracklistelement->fTrack.GetPsi() <= fToleranceDeviation)
743     ++currentmatchingindicator;
744
745   if(abs((firsttracklistelement->fTrack.GetTgl() - secondtracklistelement->fTrack.GetTgl()))/firsttracklistelement->fTrack.GetTgl() <= fToleranceDeviation)
746     ++currentmatchingindicator;
747
748   // compare number of assigned cluster hits
749   if(abs((firsttracklistelement->fTrack.GetNHits() - secondtracklistelement->fTrack.GetNHits()))/firsttracklistelement->fTrack.GetNHits() <= fToleranceDeviation)
750     ++currentmatchingindicator;
751
752   return currentmatchingindicator;
753 }
754
755 Int_t AliHLTTPCCompModelAnalysis::ComparePythiaTrackInfo(AliHLTTPCTrackList* firsttracklistelement, AliHLTTPCTrackList* secondtracklistelement)
756 {
757   // see header file for class documentation
758   // calculate matching indicator accoring to the track information
759   // ++matchingindicator for every paramter that matches
760
761   Int_t currentmatchingindicator = 0;
762
763   // tolerance range of 1 percent deviation for each quantity
764
765   // compare start point (x,y,z)
766   if(firsttracklistelement->fPythiatrack.GetFirstPointX() == secondtracklistelement->fPythiatrack.GetFirstPointX())
767     ++currentmatchingindicator;
768
769   if(firsttracklistelement->fPythiatrack.GetFirstPointY() == secondtracklistelement->fPythiatrack.GetFirstPointY())
770     ++currentmatchingindicator;
771
772   if(firsttracklistelement->fPythiatrack.GetFirstPointZ() == secondtracklistelement->fPythiatrack.GetFirstPointZ())
773     ++currentmatchingindicator;
774   
775   // compare end point
776   if(firsttracklistelement->fPythiatrack.GetLastPointX() == secondtracklistelement->fPythiatrack.GetLastPointX())
777     ++currentmatchingindicator;
778
779   if(firsttracklistelement->fPythiatrack.GetLastPointY() == secondtracklistelement->fPythiatrack.GetLastPointY())
780     ++currentmatchingindicator;
781
782   if(firsttracklistelement->fPythiatrack.GetLastPointZ() == secondtracklistelement->fPythiatrack.GetLastPointZ())
783     ++currentmatchingindicator;
784
785   // compare pt, psi, tgl
786   if(firsttracklistelement->fPythiatrack.GetPt() == secondtracklistelement->fPythiatrack.GetPt())
787     ++currentmatchingindicator;
788
789   if(firsttracklistelement->fPythiatrack.GetPsi() == secondtracklistelement->fPythiatrack.GetPsi())
790     ++currentmatchingindicator;
791
792   if(firsttracklistelement->fPythiatrack.GetTgl() == secondtracklistelement->fPythiatrack.GetTgl())
793     ++currentmatchingindicator;
794
795   // compare number of assigned cluster hits
796   if(firsttracklistelement->fPythiatrack.GetNHits() == secondtracklistelement->fPythiatrack.GetNHits())
797     ++currentmatchingindicator;
798
799   return currentmatchingindicator;
800 }
801
802 Int_t AliHLTTPCCompModelAnalysis::CreateGraphs(Bool_t relativedifferences)
803 {
804   // see header file for class documentation
805   AliHLTTPCTrackList* tracklistpointer = fFirstTrackList;
806
807   AliHLTTPCTrackList* trackmatchingpointer;
808
809   // set up histogram ranges
810   Double_t difffirstxmin, difffirstxmax;
811   Double_t difffirstymin, difffirstymax;
812   Double_t difffirstzmin, difffirstzmax;
813   Int_t difffirstxbins, difffirstybins, difffirstzbins;
814
815   Double_t difflastxmin, difflastxmax;
816   Double_t difflastymin, difflastymax;
817   Double_t difflastzmin, difflastzmax;
818   Int_t difflastxbins, difflastybins, difflastzbins; 
819
820   Double_t diffptmin, diffptmax;
821   Double_t diffpsimin, diffpsimax;
822   Double_t difftglmin, difftglmax;
823   Int_t diffptbins, diffpsibins, difftglbins;
824
825   Double_t diffclustermin, diffclustermax;
826   Int_t diffclusterbins;
827  
828   // resolution histograms (currently not working since pterr = 0 for every track!)
829   Double_t diffpterrmin, diffpterrmax;
830   Double_t diffpsierrmin, diffpsierrmax;
831   Double_t difftglerrmin, difftglerrmax;
832   Int_t diffpterrbins, diffpsierrbins, difftglerrbins;
833
834   if(!relativedifferences)
835     {
836       difffirstxmin = -2;
837       difffirstxmax = +2;
838       difffirstxbins = (Int_t) ((difffirstxmax - difffirstxmin)/0.0001);
839
840       difffirstymin = -2;
841       difffirstymax = +2;
842       difffirstybins = (Int_t) ((difffirstymax - difffirstymin)/0.0001);
843
844       difffirstzmin = -2;
845       difffirstzmax = +2;
846       difffirstzbins = (Int_t) ((difffirstzmax - difffirstzmin)/0.0001);
847
848       difflastxmin = -2;
849       difflastxmax = +2;
850       difflastxbins = (Int_t) ((difflastxmax - difflastxmin)/0.0001);
851
852       difflastymin = -2;
853       difflastymax = +2;
854       difflastybins = (Int_t) ((difflastymax - difflastymin)/0.0001);
855
856       difflastzmin = -2;
857       difflastzmax = +2;
858       difflastzbins = (Int_t) ((difflastzmax - difflastzmin)/0.0001);
859
860       diffptmin = -1;
861       diffptmax = +1;
862       diffptbins = (Int_t) ((diffptmax - diffptmin)/0.0001);
863       
864       diffpsimin = -1;
865       diffpsimax = +1;
866       diffpsibins = (Int_t) ((diffpsimax - diffpsimin)/0.0001);
867
868       difftglmin = -1;
869       difftglmax = +1;
870       difftglbins = (Int_t) ((difftglmax - difftglmin)/0.0001);
871
872       diffclustermin = -50;
873       diffclustermax = +50;
874       diffclusterbins = (Int_t) ((diffclustermax - diffclustermin)/1);
875
876       //#if 0
877       diffpterrmin = -1;
878       diffpterrmax = 1;
879       diffpterrbins = (Int_t) ((diffpterrmax - diffpterrmin)/1);
880
881       diffpsierrmin = -1;
882       diffpsierrmax = 1;
883       diffpsierrbins = (Int_t) ((diffpsierrmax - diffpsierrmin)/1);
884
885       difftglerrmin = -1;
886       difftglerrmax = 1;
887       difftglerrbins = (Int_t) ((difftglerrmax - difftglerrmin)/1);
888       //#endif
889
890     }
891   else
892     {
893       difffirstxmin = -1;
894       difffirstxmax = +1;
895       difffirstxbins = (Int_t) ((difffirstxmax - difffirstxmin)/0.0001);
896
897       difffirstymin = -1;
898       difffirstymax = +1;
899       difffirstybins = (Int_t) ((difffirstymax - difffirstymin)/0.0001);
900
901       difffirstzmin = -1;
902       difffirstzmax = +1;
903       difffirstzbins = (Int_t) ((difffirstzmax - difffirstzmin)/0.0001);
904
905       difflastxmin = -1;
906       difflastxmax = +1;
907       difflastxbins = (Int_t) ((difflastxmax - difflastxmin)/0.0001);
908
909       difflastymin = -1;
910       difflastymax = +1;
911       difflastybins = (Int_t) ((difflastymax - difflastymin)/0.0001);
912
913       difflastzmin = -1;
914       difflastzmax = +1;
915       difflastzbins = (Int_t) ((difflastzmax - difflastzmin)/0.0001);
916
917       diffptmin = -1;
918       diffptmax = +1;
919       diffptbins = (Int_t) ((diffptmax - diffptmin)/0.0001);
920       
921       diffpsimin = -1;
922       diffpsimax = +1;
923       diffpsibins = (Int_t) ((diffpsimax - diffpsimin)/0.0001);
924
925       difftglmin = -1;
926       difftglmax = +1;
927       difftglbins = (Int_t) ((difftglmax - difftglmin)/0.0001);
928
929       diffclustermin = -1;
930       diffclustermax = +1;
931       diffclusterbins = (Int_t) ((diffclustermax - diffclustermin)/0.0001);
932
933       //#if 0
934       diffpterrmin = -1;
935       diffpterrmax = 1;
936       diffpterrbins = (Int_t) ((diffpterrmax - diffpterrmin)/0.0001);
937
938       diffpsierrmin = -1;
939       diffpsierrmax = 1;
940       diffpsierrbins = (Int_t) ((diffpsierrmax - diffpsierrmin)/0.0001);
941
942       difftglerrmin = -1;
943       difftglerrmax = 1;
944       difftglerrbins = (Int_t) ((difftglerrmax - difftglerrmin)/0.0001);
945       //#endif
946     }
947
948   // intialise histogramms
949   TH1F* firstxhisto = new TH1F("Differences of first x (original - secondary) track", "Differences of first x (original - secondary) track", difffirstxbins, difffirstxmin, difffirstxmax);
950   TH1F* firstyhisto = new TH1F("Differences of first y (original - secondary) track", "Differences of first y (original - secondary) track", difffirstybins, difffirstymin, difffirstymax);
951  TH1F* firstzhisto = new TH1F("Differences of first z (original - secondary) track", "Differences of first z (original - secondary) track", difffirstzbins, difffirstzmin, difffirstzmax);
952  TH1F* lastxhisto = new TH1F("Differences of last x (original - secondary) track", "Differences of last x (original - secondary) track", difflastxbins, difflastxmin, difflastxmax);
953  TH1F* lastyhisto = new TH1F("Differences of last y (original - secondary) track", "Differences of last y (original - secondary) track", difflastybins, difflastymin, difflastymax);
954  TH1F* lastzhisto = new TH1F("Differences of last z (original - secondary) track", "Differences of last z (original - secondary) track", difflastzbins, difflastzmin, difflastzmax);
955  TH1F* pthisto = new TH1F("Differences of pt (original - secondary) track", "Differences of pt (original - secondary) track", diffptbins, diffptmin, diffptmax);
956  TH1F* psihisto = new TH1F("Differences of psi (original - secondary) track", "Differences of psi (original - secondary) track", diffpsibins, diffpsimin, diffpsimax);
957  TH1F* tglhisto = new TH1F("Differences of tgl (original - secondary) track", "Differences of tgl (original - secondary) track", difftglbins, difftglmin, difftglmax);
958  TH1F* clusterhisto = new TH1F("Differences of asserted clusters (original - secondary) track", "Differences of asserted clusters (original - secondary) track", diffclusterbins, diffclustermin, diffclustermax);
959
960  //#if 0
961  // commented out since pterr is zero for every track!
962  TH1F* pterrhisto = new TH1F("Differences of pt error (original - secondary) track", "Differences of pt error (original - secondary) track", diffpterrbins, diffpterrmin, diffpterrmax);
963  TH1F* psierrhisto = new TH1F("Differences of psi error (original - secondary) track", "Differences of psi error (original - secondary) track", diffpsierrbins, diffpsierrmin, diffpsierrmax);
964  TH1F* tglerrhisto = new TH1F("Differences of tgl error (original - secondary) track", "Differences of tgl error (original - secondary) track", difftglerrbins, difftglerrmin, difftglerrmax);
965
966  //#endif
967
968  // initialise histograms for tracking efficiency against pt
969  // relative p_t resolution: 1.2% -> take 0.1GeV * 1.2 % --> binsize 0.001 sufficient to grant correct resolution for low pt
970  TH1F* firsttracks = new TH1F("pt occurrence original", "Occurrence of pt in original tracks", 10000, 0, 10);
971  TH1F* matchedtrackeff = new TH1F("matchedtreffeff", "Occurrence of 2ndary tracks with good pt", 10000, 0, 10);
972  TH1F* trackeff = new TH1F("tracking efficiency vs. pt", "Tracking efficiency vs. pt", 10000, 0, 10);
973
974  // evaluate quality of fit:
975  TH1I* matchinghisto = new TH1I("Matching indicator (5 - 10)", "Matching indicator (5 - 10)", 11, 0, 11);
976
977   while(tracklistpointer != NULL)
978     {
979       // if currently processed track is matched, store differences of their parameters in histogram
980       if(tracklistpointer->fMatchingindicator > 0)
981         {
982           // matching track
983           trackmatchingpointer = tracklistpointer->fMatchingtrack;
984
985           // fill histograms for trackingefficiency vs. pt
986           firsttracks->Fill(tracklistpointer->fTrack.GetPt(),1);
987           
988           if(abs(tracklistpointer->fTrack.GetPt()-trackmatchingpointer->fTrack.GetPt()) < 0.012*tracklistpointer->fTrack.GetPt())
989             {
990               matchedtrackeff->Fill(tracklistpointer->fTrack.GetPt(),1);
991             }
992
993           //tracklistpointer = tracklistpointer->fNext; // only efficiency is considered...
994           //continue; // only efficiency is considered...
995
996           if(relativedifferences == 1) // fill histogram with relative differences
997             {
998
999               // check if first track parameters are not zero!
1000               if (tracklistpointer->fTrack.GetFirstPointX()==0)
1001                 {
1002                   HLTWarning("Warning! First x of original track is zero, relative differences cannot be calculated!");
1003                 };
1004               if (tracklistpointer->fTrack.GetFirstPointY()==0)
1005                 {
1006                   HLTWarning("Warning! First y of original track is zero, relative differences cannot be calculated!");
1007                 };
1008               if (tracklistpointer->fTrack.GetFirstPointZ()==0)
1009                 {
1010                   HLTWarning("Warning! First z of original track is zero, relative differences cannot be calculated!");
1011                 };
1012               if (tracklistpointer->fTrack.GetLastPointX()==0)
1013                 {
1014                   HLTWarning("Warning! Last x of original track is zero, relative differences cannot be calculated!");
1015                 };
1016               if (tracklistpointer->fTrack.GetLastPointY()==0)
1017                 {
1018                   HLTWarning("Warning! Last y of original track is zero, relative differences cannot be calculated!");
1019                 };
1020               if (tracklistpointer->fTrack.GetLastPointZ()==0)
1021                 {
1022                   HLTWarning("Warning! Last z of original track is zero, relative differences cannot be calculated!");
1023                 };
1024               if (tracklistpointer->fTrack.GetPt()==0)
1025                 {
1026                   HLTWarning("Warning! Pt of original track is zero, relative differences cannot be calculated!");
1027                 };
1028               if (tracklistpointer->fTrack.GetPsi()==0)
1029                 {
1030                   HLTWarning("Warning! Psi of original track is zero, relative differences cannot be calculated!");
1031                 };
1032               if (tracklistpointer->fTrack.GetTgl()==0)
1033                 {
1034                   HLTWarning("Warning! Tgl of original track is zero, relative differences cannot be calculated!");
1035                 };
1036               firstxhisto->Fill((tracklistpointer->fTrack.GetFirstPointX()-trackmatchingpointer->fTrack.GetFirstPointX())/tracklistpointer->fTrack.GetFirstPointX(),1);
1037               firstyhisto->Fill((tracklistpointer->fTrack.GetFirstPointY()-trackmatchingpointer->fTrack.GetFirstPointY())/tracklistpointer->fTrack.GetFirstPointY(),1);
1038               firstzhisto->Fill((tracklistpointer->fTrack.GetFirstPointZ()-trackmatchingpointer->fTrack.GetFirstPointZ())/tracklistpointer->fTrack.GetFirstPointZ(),1);
1039               lastxhisto->Fill((tracklistpointer->fTrack.GetLastPointX()-trackmatchingpointer->fTrack.GetLastPointX())/tracklistpointer->fTrack.GetLastPointX(),1);
1040               lastyhisto->Fill((tracklistpointer->fTrack.GetLastPointY()-trackmatchingpointer->fTrack.GetLastPointY())/tracklistpointer->fTrack.GetLastPointY(),1);
1041               lastzhisto->Fill((tracklistpointer->fTrack.GetLastPointZ()-trackmatchingpointer->fTrack.GetLastPointZ())/tracklistpointer->fTrack.GetLastPointZ(),1);
1042               pthisto->Fill((tracklistpointer->fTrack.GetPt()-trackmatchingpointer->fTrack.GetPt())/tracklistpointer->fTrack.GetPt(),1);
1043               psihisto->Fill((tracklistpointer->fTrack.GetPsi()-trackmatchingpointer->fTrack.GetPsi())/tracklistpointer->fTrack.GetPsi(),1);
1044               tglhisto->Fill((tracklistpointer->fTrack.GetTgl()-trackmatchingpointer->fTrack.GetTgl())/tracklistpointer->fTrack.GetTgl(),1);
1045               clusterhisto->Fill((tracklistpointer->fTrack.GetNHits()-trackmatchingpointer->fTrack.GetNHits())/tracklistpointer->fTrack.GetNHits(),1);
1046
1047               //#if 0
1048               pterrhisto->Fill((tracklistpointer->fTrack.GetPterr()-trackmatchingpointer->fTrack.GetPterr())/tracklistpointer->fTrack.GetPterr(),1);
1049               psierrhisto->Fill((tracklistpointer->fTrack.GetPsierr()-trackmatchingpointer->fTrack.GetPsierr())/tracklistpointer->fTrack.GetPsierr(),1);
1050               tglerrhisto->Fill((tracklistpointer->fTrack.GetTglerr()-trackmatchingpointer->fTrack.GetTglerr())/tracklistpointer->fTrack.GetTglerr(),1);
1051
1052               //HLTInfo("Pterr: 1st:%f  2nd:%f   value:%f",tracklistpointer->fTrack.GetPterr(),trackmatchingpointer->fTrack.GetPterr(),(tracklistpointer->fTrack.GetPterr()-trackmatchingpointer->fTrack.GetPterr())/tracklistpointer->fTrack.GetPterr());
1053               //HLTInfo("Psierr: 1st:%f  2nd:%f   value:%f",tracklistpointer->fTrack.GetPsierr(),trackmatchingpointer->fTrack.GetPsierr(),(tracklistpointer->fTrack.GetPsierr()-trackmatchingpointer->fTrack.GetPsierr())/tracklistpointer->fTrack.GetPsierr());
1054               //HLTInfo("Tglerr: 1st:%f  2nd:%f   value:%f",tracklistpointer->fTrack.GetTglerr(),trackmatchingpointer->fTrack.GetTglerr(),(tracklistpointer->fTrack.GetTglerr()-trackmatchingpointer->fTrack.GetTglerr())/tracklistpointer->fTrack.GetTglerr());
1055               //#endif
1056             }
1057           else // otherwise fill histogram with absolute differences
1058             {
1059               firstxhisto->Fill(tracklistpointer->fTrack.GetFirstPointX()-trackmatchingpointer->fTrack.GetFirstPointX(),1);
1060               firstyhisto->Fill(tracklistpointer->fTrack.GetFirstPointY()-trackmatchingpointer->fTrack.GetFirstPointY(),1);
1061               firstzhisto->Fill(tracklistpointer->fTrack.GetFirstPointZ()-trackmatchingpointer->fTrack.GetFirstPointZ(),1);
1062               lastxhisto->Fill(tracklistpointer->fTrack.GetLastPointX()-trackmatchingpointer->fTrack.GetLastPointX(),1);
1063               lastyhisto->Fill(tracklistpointer->fTrack.GetLastPointY()-trackmatchingpointer->fTrack.GetLastPointY(),1);
1064               lastzhisto->Fill(tracklistpointer->fTrack.GetLastPointZ()-trackmatchingpointer->fTrack.GetLastPointZ(),1);
1065               pthisto->Fill(tracklistpointer->fTrack.GetPt()-trackmatchingpointer->fTrack.GetPt(),1);
1066               psihisto->Fill(tracklistpointer->fTrack.GetPsi()-trackmatchingpointer->fTrack.GetPsi(),1);
1067               tglhisto->Fill(tracklistpointer->fTrack.GetTgl()-trackmatchingpointer->fTrack.GetTgl(),1);
1068               clusterhisto->Fill(tracklistpointer->fTrack.GetNHits()-trackmatchingpointer->fTrack.GetNHits(),1);
1069               //#if 0
1070               // commented out since pterr is always zero for every track!
1071               pterrhisto->Fill(tracklistpointer->fTrack.GetPterr()-trackmatchingpointer->fTrack.GetPterr(),1);
1072               psierrhisto->Fill(tracklistpointer->fTrack.GetPsierr()-trackmatchingpointer->fTrack.GetPsierr(),1);
1073               tglerrhisto->Fill(tracklistpointer->fTrack.GetTglerr()-trackmatchingpointer->fTrack.GetTglerr(),1);
1074               //#endif
1075             }
1076
1077           // fill histogram that determines the quality of the fit
1078           matchinghisto->Fill(tracklistpointer->fMatchingindicator,1);
1079         }
1080
1081       tracklistpointer = tracklistpointer->fNext;
1082     }
1083
1084   trackeff->Divide(matchedtrackeff, firsttracks,1,1,"");
1085
1086   // write histograms to root file specified in command line argument -graphs <filename>.ROOT
1087   if(!fGraphFileName.IsNull())
1088     {
1089       TFile* graphrootfile = new TFile(fGraphFileName, "update");
1090        graphrootfile->WriteObject(firstxhisto,"firstxhistogram");
1091       //firstxhisto->Write();
1092       graphrootfile->WriteObject(firstyhisto,"firstyhistogram");
1093       //firstyhisto->Write();
1094       graphrootfile->WriteObject(firstzhisto,"firstzhistogram");
1095       //firstzhisto->Write();
1096       graphrootfile->WriteObject(lastxhisto,"lastxhistogram");
1097       //lastxhisto->Write();
1098       graphrootfile->WriteObject(lastyhisto,"lastyhistogram");
1099       //lastyhisto->Write();
1100       graphrootfile->WriteObject(lastzhisto,"lastzhistogram");
1101       //lastzhisto->Write();
1102       graphrootfile->WriteObject(pthisto,"pthistogram");
1103       //pthisto->Write();
1104       graphrootfile->WriteObject(psihisto,"psihistogram");
1105       //psihisto->Write();
1106       graphrootfile->WriteObject(tglhisto,"tglhistogram");
1107       //tglhisto->Write();
1108       graphrootfile->WriteObject(clusterhisto,"clusterhistogram");
1109       //clusterhisto->Write();
1110       graphrootfile->WriteObject(matchinghisto,"matchinghistogram");
1111       //matchinghisto->Write();
1112       graphrootfile->WriteObject(firsttracks, "firsttrackeff");
1113       graphrootfile->WriteObject(matchedtrackeff, "secondtrackeff");
1114       graphrootfile->WriteObject(trackeff, "trackingefficiency");
1115
1116       // errors in tracking (commented out since pterr is always 0 for every track!)
1117       graphrootfile->WriteObject(pterrhisto, "pterrhistogram");
1118       //pterrhisto->Write();
1119       graphrootfile->WriteObject(psierrhisto, "psierrhistogram");
1120       //psierrhisto->Write();
1121       graphrootfile->WriteObject(tglerrhisto, "tglerrhistogram");
1122       //tglerrhisto->Write();
1123       graphrootfile->Close();
1124     }
1125   else
1126     {
1127       HLTError("Error! No file for graphical output specified.");
1128       return EINVAL;
1129     }
1130
1131   return 0;
1132 }
1133
1134 Int_t AliHLTTPCCompModelAnalysis::DisplayModelResults()
1135 {
1136   // see header file for class documentation
1137
1138   //print out parameters of discarded track:
1139   AliHLTTPCTrackList* trackprintpointer = fTrackListPointer;
1140   AliHLTTPCTrackList* tracklistdeleter= trackprintpointer;
1141
1142  FILE* dumpfile = NULL;
1143
1144  if(!fDumpFileName.IsNull())
1145    {
1146      // open new file specified by command line argument
1147      dumpfile = fopen(fDumpFileName.Data(),"a");
1148      fprintf(dumpfile,"---------------MODEL ANALYSIS--------------- \n");
1149      fprintf(dumpfile,"---------------DISCARDED TRACKS: %d --------------- \n", fTrashTracks);
1150    };
1151  
1152  if(fTrackListPointer == NULL)
1153    {
1154      HLTInfo("No tracks discarded");
1155      HLTInfo("--------------");
1156    }
1157  else
1158    {
1159      
1160      HLTInfo("---------------DISCARDED TRACKS: %d ---------------", fTrashTracks);
1161      
1162      Int_t trashtrackcounter = 1;
1163      
1164      while(trackprintpointer != NULL)
1165        {
1166          
1167          // infos about found track
1168          HLTInfo("%d : Discarding track with %d clusters.", trashtrackcounter, trackprintpointer->fTrack.GetNHits());
1169          //PYTHIA INFORMATION ABOUT PARTICLE ID
1170          HLTInfo("Track parameters of discarded track:");       
1171          HLTInfo("First x: %9.6f \t first y: %9.6f \t first z: %9.6f",trackprintpointer->fTrack.GetFirstPointX(),trackprintpointer->fTrack.GetFirstPointY(),trackprintpointer->fTrack.GetFirstPointZ()); 
1172          HLTInfo(" Last x: %9.6f \t  last y: %9.6f \t  last z: %9.6f", trackprintpointer->fTrack.GetLastPointX(),trackprintpointer->fTrack.GetLastPointY(),trackprintpointer->fTrack.GetLastPointZ());
1173          HLTInfo("     Pt: %9.6f \t     Psi: %9.6f \t     Tgl: %9.6f", trackprintpointer->fTrack.GetPt(), trackprintpointer->fTrack.GetPsi(), trackprintpointer->fTrack.GetTgl());
1174          
1175          // write results to file if specified by command line argument
1176          if(!fDumpFileName.IsNull())
1177            {
1178              fprintf(dumpfile,"%d : Discarding track with %d clusters. \n", trashtrackcounter, trackprintpointer->fTrack.GetNHits());
1179              fprintf(dumpfile,"Track parameters of discarded track: \n");       
1180              fprintf(dumpfile,"First x: %9.6f \t first y: %9.6f \t first z: %9.6f \n",trackprintpointer->fTrack.GetFirstPointX(),trackprintpointer->fTrack.GetFirstPointY(),trackprintpointer->fTrack.GetFirstPointZ()); 
1181              fprintf(dumpfile," Last x: %9.6f \t  last y: %9.6f \t  last z: %9.6f \n", trackprintpointer->fTrack.GetLastPointX(),trackprintpointer->fTrack.GetLastPointY(),trackprintpointer->fTrack.GetLastPointZ());
1182              fprintf(dumpfile,"     Pt: %9.6f \t     Psi: %9.6f \t     Tgl: %9.6f \n", trackprintpointer->fTrack.GetPt(), trackprintpointer->fTrack.GetPsi(), trackprintpointer->fTrack.GetTgl());
1183            };
1184          
1185          // comparison with pythia information
1186          if(trackprintpointer->fWronglydiscarded)
1187            {
1188              HLTInfo("Found track has been wrongly discarded according to Pythia information.");
1189              
1190              // write results to file if specified by command line argument
1191              if(!fDumpFileName.IsNull())
1192                {
1193                   fprintf(dumpfile,"Found track has been wrongly discarded accoring to Pythia information. \n");
1194                };
1195              
1196             };
1197          
1198          // found pt must be in range pt_pythia \pm 10% to be accepted
1199          Double_t ptmin = trackprintpointer->fPythiatrack.GetPt() - 0.1*trackprintpointer->fPythiatrack.GetPt();
1200          Double_t ptmax = trackprintpointer->fPythiatrack.GetPt() + 0.1*trackprintpointer->fPythiatrack.GetPt();
1201          
1202          if( (trackprintpointer->fTrack.GetPt() < ptmin) ||(trackprintpointer->fTrack.GetPt() > ptmax) )
1203            {
1204              HLTInfo("Pt of found track %f differs more than 10 %% from pt of pythia track %f.",trackprintpointer->fTrack.GetPt(), trackprintpointer->fPythiatrack.GetPt());
1205              
1206              if(!fDumpFileName.IsNull())
1207                {
1208                  // write result to file if specified by command line argument
1209                  fprintf(dumpfile,"Pt of found track %f differs more than 10 %% from pt of pythia track %f. \n",trackprintpointer->fTrack.GetPt(), trackprintpointer->fPythiatrack.GetPt());
1210                };
1211              
1212            };
1213          
1214          HLTInfo("--------------");
1215          
1216          if(!fDumpFileName.IsNull())
1217            {
1218              fprintf(dumpfile,"-------------- \n");
1219            }; 
1220          
1221          
1222           // go to next element in trash track list
1223          tracklistdeleter = trackprintpointer;
1224          trackprintpointer = trackprintpointer->fNext;
1225          
1226          ++trashtrackcounter;
1227          
1228          // free space 
1229          delete tracklistdeleter;
1230          tracklistdeleter = NULL;
1231          
1232        } // end of while(trackpointer != NULL)
1233
1234    } // end of else
1235  
1236  // print out number of noise clusters (not assigned to any track and not valid)
1237  HLTInfo("Number of discarded clusters not assigned to any track: %d", fTotalDiscardedClusters);
1238  HLTInfo("--------------");
1239
1240  // write results to file if specified by command line argument
1241  if(!fDumpFileName.IsNull())
1242    {
1243      fprintf(dumpfile,"Number of discarded clusters not assigned to any track: %d \n", fTotalDiscardedClusters);
1244      fprintf(dumpfile,"-------------- \n");
1245    };
1246  
1247  // print out paramters of discarded valuable clusters
1248  HLTInfo("Number of discarded VALUABLE clusters: %d", fValuableDiscardedClusters);
1249
1250  HLTInfo("--------------");
1251  
1252  if(!fDumpFileName.IsNull())
1253    {
1254
1255      fprintf(dumpfile,"Number of discarded VALUABLE clusters: %d \n", fValuableDiscardedClusters);
1256      fclose(dumpfile);
1257    };
1258  
1259   return 0;
1260 }
1261
1262 Int_t AliHLTTPCCompModelAnalysis::DisplayTrackResults()
1263 {
1264   // see header file for class documentation
1265   HLTInfo("---------------CLUSTER ANALYSIS---------------");
1266   if(CompareClusters(1) != 0)
1267     { 
1268       return EINVAL;
1269     }
1270
1271   // start with comparison
1272   if(CompareTracks() != 0)
1273     {
1274       return EINVAL;
1275     };
1276
1277   // if dumptofile is activated, append results to output analysis file
1278   FILE* dumpfile = NULL;
1279   
1280   if(!fDumpFileName.IsNull())
1281     {
1282       // open new file specified by command line argument
1283       dumpfile = fopen(fDumpFileName.Data(),"a");
1284       
1285       fprintf(dumpfile,"---------------TRACK ANALYSIS--------------- \n");
1286       
1287     }
1288   
1289   // print out number of compared tracks
1290   HLTInfo("---------------TRACK ANALYSIS---------------");
1291   HLTInfo("---------------ORIGINAL TRACKS: %d ---------------", fFirstTrackArray.GetNTracks());
1292   HLTInfo("Number of tracks with pt < 0.1 GeV: %d", fFirstTrashTracks);
1293   HLTInfo("Number of matched tracks with pt < 0.1 GeV: %d", fMatchedFirstTrashTracks);
1294   //PYTHIA INFORMATION ABOUT PARTICLE IDs
1295   HLTInfo("---------------2NDARY TRACKS: %d ---------------", fSecondTrackArray.GetNTracks());
1296   HLTInfo("Number of tracks with pt < 0.1 GeV: %d", fSecondTrashTracks);
1297   HLTInfo("Number of matched tracks with pt < 0.1 GeV: %d", fMatchedSecondTrashTracks);
1298   //PYTHIA INFORMATION ABOUT PARTICLE IDs
1299   HLTInfo("--------------");
1300   HLTInfo("Comparison of tracks within parameter precision: %f", fToleranceDeviation);
1301   HLTInfo("Number of compared tracks: %d", fTotalComparedTracks);
1302   HLTInfo("Number of unmatched original tracks: %d", fFirstUnmatchedTracks);
1303   HLTInfo("Number of unmatched secondary tracks: %d", fSecondUnmatchedTracks);
1304   HLTInfo("--------------");
1305   
1306   // print results to file
1307   if(!fDumpFileName.IsNull())
1308     {
1309       fprintf(dumpfile, "---------------%d ORIGINAL TRACKS---------------\n", fFirstTrackArray.GetNTracks());
1310       fprintf(dumpfile,"Number of tracks with pt < 0.1 GeV: %d \n", fFirstTrashTracks);
1311       fprintf(dumpfile,"Number of matched tracks with pt < 0.1 GeV: %d \n", fMatchedFirstTrashTracks);
1312       fprintf(dumpfile,"---------------%d 2NDARY TRACKS---------------\n", fSecondTrackArray.GetNTracks());
1313       fprintf(dumpfile,"Number of tracks with pt < 0.1 GeV: %d \n", fSecondTrashTracks);
1314       fprintf(dumpfile,"Number of matched tracks with pt < 0.1 GeV: %d \n", fMatchedSecondTrashTracks);
1315       fprintf(dumpfile,"--------------\n");
1316       fprintf(dumpfile,"Comparison of tracks within parameter precision: %f \n", fToleranceDeviation);
1317       fprintf(dumpfile,"Number of compared tracks: %d \n", fTotalComparedTracks); 
1318       fprintf(dumpfile,"Number of unmatched original tracks: %d \n", fFirstUnmatchedTracks); 
1319       fprintf(dumpfile,"Number of unmatched secondary tracks: %d \n", fSecondUnmatchedTracks); 
1320       fprintf(dumpfile,"--------------\n");
1321     }
1322
1323   //////////////////////////////////////////////////////////////////////
1324   // additional files temporarily necessary for output information
1325   FILE* infofile = fopen("/afsuser/jwagner/TrackerTest_25092007/pp-ca/fullanalysis08012008/trackingefficiency.out","a");
1326   FILE* info2file = fopen("/afsuser/jwagner/TrackerTest_25092007/pp-ca/fullanalysis08012008/parameters.out", "a");
1327
1328   // consistent = 0 if tracks and second tracks match perfectly, i.e. no doubly assigned tracks,etc.
1329   Int_t consistentfirst = fFirstTrackArray.GetNTracks() - fTotalComparedTracks - fFirstUnmatchedTracks;
1330   Int_t consistentsecond = fSecondTrackArray.GetNTracks() - fTotalComparedTracks - fSecondUnmatchedTracks;
1331
1332   //fprintf(infofile, "1st tracks, 2nd tracks, compared, 1st uncompared, 2nd uncompared, cons.1st, cons.2nd \n");
1333   fprintf(info2file, " %d \t %d \t %d \t %d \t %d \t %d \t %d \n", fFirstTrackArray.GetNTracks(), fSecondTrackArray.GetNTracks(), fTotalComparedTracks, fFirstUnmatchedTracks, fSecondUnmatchedTracks, consistentfirst, consistentsecond);
1334
1335   //fprintf(infofile, "1st trash tracks, 2nd trash tracks, 1st matched trash tracks, 2nd matched trash tracks \n");
1336   fprintf(infofile, "%d \t %d \t %d \t %d \n", fFirstTrashTracks, fSecondTrashTracks, fMatchedFirstTrashTracks, fMatchedSecondTrashTracks);
1337
1338   fclose(infofile);
1339   fclose(info2file);
1340   ////////////////////////////////////////////////////////////////////
1341
1342   // print out deviations
1343   Int_t tracknumber = 1;
1344   AliHLTTPCTrackList* listprintpointer = fFirstTrackList;
1345   while(listprintpointer != NULL)
1346     {
1347       // print out parameters of original track in comparison to secondary track:
1348       
1349       if(listprintpointer->fMatchingindicator != 0)
1350         {
1351           
1352 #if 0
1353           HLTInfo("Track %d:", tracknumber);
1354           HLTInfo("Original track matched to secondary track with matchingindicator %d", listprintpointer->fMatchingindicator);
1355           
1356           HLTInfo("Parameter comparison: Original vs. Secondary");
1357
1358           HLTInfo("Clusters: %d \t %d ",listprintpointer->fTrack.GetNHits(), listprintpointer->fMatchingtrack->fTrack.GetNHits());
1359           HLTInfo("First x: %9.6f \t %9.6f", listprintpointer->fTrack.GetFirstPointX(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointX());
1360           HLTInfo("First y: %9.6f \t %9.6f", listprintpointer->fTrack.GetFirstPointY(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointY());
1361           HLTInfo("First z: %9.6f \t %9.6f", listprintpointer->fTrack.GetFirstPointZ(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointZ()); 
1362           HLTInfo("Last x: %9.6f \t %9.6f", listprintpointer->fTrack.GetLastPointX(), listprintpointer->fMatchingtrack->fTrack.GetLastPointX());
1363           HLTInfo("Last y: %9.6f \t %9.6f", listprintpointer->fTrack.GetLastPointY(), listprintpointer->fMatchingtrack->fTrack.GetLastPointY());
1364           HLTInfo("Last z: %9.6f \t %9.6f", listprintpointer->fTrack.GetLastPointZ(), listprintpointer->fMatchingtrack->fTrack.GetLastPointZ());
1365           HLTInfo("    Pt: %9.6f \t %9.6f", listprintpointer->fTrack.GetPt(), listprintpointer->fMatchingtrack->fTrack.GetPt());  
1366           HLTInfo("   Psi: %9.6f \t %9.6f", listprintpointer->fTrack.GetPsi(), listprintpointer->fMatchingtrack->fTrack.GetPsi());
1367           HLTInfo("   Tgl: %9.6f \t %9.6f", listprintpointer->fTrack.GetTgl(), listprintpointer->fMatchingtrack->fTrack.GetTgl());
1368
1369           HLTInfo(" Pterr: %9.6f \t %9.6f", listprintpointer->fTrack.GetPterr(), listprintpointer->fMatchingtrack->fTrack.GetPterr());  
1370           HLTInfo("Psierr: %9.6f \t %9.6f", listprintpointer->fTrack.GetPsierr(), listprintpointer->fMatchingtrack->fTrack.GetPsierr());
1371           HLTInfo("Tglerr: %9.6f \t %9.6f", listprintpointer->fTrack.GetTglerr(), listprintpointer->fMatchingtrack->fTrack.GetTglerr());
1372          
1373           HLTInfo("--------------");
1374 #endif
1375           // print these results to file
1376           if(!fDumpFileName.IsNull())
1377             {
1378               fprintf(dumpfile, "Track %d: \n", tracknumber);
1379               fprintf(dumpfile, "Original track matched to secondary track with matchingindicator %d \n", listprintpointer->fMatchingindicator); 
1380               fprintf(dumpfile, "Parameter comparison: Original vs. Secondary \n"); 
1381               fprintf(dumpfile, "Clusters: %d \t %d \n ",listprintpointer->fTrack.GetNHits(), listprintpointer->fMatchingtrack->fTrack.GetNHits());
1382               fprintf(dumpfile, "First x: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetFirstPointX(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointX());
1383               fprintf(dumpfile, "First y: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetFirstPointY(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointY()); 
1384               fprintf(dumpfile, "First z: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetFirstPointZ(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointZ()); 
1385               fprintf(dumpfile, "Last x: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetLastPointX(), listprintpointer->fMatchingtrack->fTrack.GetLastPointX());  
1386               fprintf(dumpfile, "Last y: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetLastPointY(), listprintpointer->fMatchingtrack->fTrack.GetLastPointY());
1387               fprintf(dumpfile, "Last z: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetLastPointZ(), listprintpointer->fMatchingtrack->fTrack.GetLastPointZ());
1388               fprintf(dumpfile, "    Pt: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetPt(), listprintpointer->fMatchingtrack->fTrack.GetPt());  
1389               fprintf(dumpfile, "   Psi: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetPsi(), listprintpointer->fMatchingtrack->fTrack.GetPsi()); 
1390               fprintf(dumpfile, "   Tgl: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetTgl(), listprintpointer->fMatchingtrack->fTrack.GetTgl());
1391               fprintf(dumpfile, "--------------\n");
1392             }
1393           ++tracknumber;
1394         }
1395      
1396       listprintpointer = listprintpointer->fNext;
1397     }
1398
1399
1400   // print out not matched tracks from first track array:
1401   listprintpointer = fFirstTrackList;
1402   Int_t notmatchedtracknumber = 1;
1403
1404   while(listprintpointer != NULL)
1405     {
1406       if(listprintpointer->fMatchingindicator == 0)
1407         {
1408 #if 0     
1409           HLTInfo("Original Track, not matched with secondary track %d:", notmatchedtracknumber);
1410           HLTInfo("Clusters: %d",listprintpointer->fTrack.GetNHits());
1411           //PYTHIA INFORMATION ABOUT PARTICLE ID
1412           HLTInfo("First x: %9.6f \t first y: %9.6f \t first z: %9.6f",listprintpointer->fTrack.GetFirstPointX(),listprintpointer->fTrack.GetFirstPointY(),listprintpointer->fTrack.GetFirstPointZ()); 
1413           HLTInfo(" Last x: %9.6f \t  last y: %9.6f \t  last z: %9.6f", listprintpointer->fTrack.GetLastPointX(),listprintpointer->fTrack.GetLastPointY(),listprintpointer->fTrack.GetLastPointZ());
1414           HLTInfo("     Pt: %9.6f \t     Psi: %9.6f \t     Tgl: %9.6f", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl());
1415           HLTInfo("--------------");
1416 #endif
1417  
1418           // print these results to file
1419           if(!fDumpFileName.IsNull())
1420             {
1421               fprintf(dumpfile, "Original Track, not matched with secondary track %d: \n", notmatchedtracknumber); 
1422               fprintf(dumpfile, "Clusters: %d \n",listprintpointer->fTrack.GetNHits());
1423               fprintf(dumpfile, "First x: %9.6f \t first y: %9.6f \t first z: %9.6f \n",listprintpointer->fTrack.GetFirstPointX(),listprintpointer->fTrack.GetFirstPointY(),listprintpointer->fTrack.GetFirstPointZ());
1424               fprintf(dumpfile, " Last x: %9.6f \t  last y: %9.6f \t  last z: %9.6f \n", listprintpointer->fTrack.GetLastPointX(),listprintpointer->fTrack.GetLastPointY(),listprintpointer->fTrack.GetLastPointZ());
1425               fprintf(dumpfile, "     Pt: %9.6f \t     Psi: %9.6f \t     Tgl: %9.6f \n", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl()); 
1426               fprintf(dumpfile, "--------------\n");
1427             }
1428
1429           ++notmatchedtracknumber;
1430         }
1431
1432       listprintpointer = listprintpointer->fNext;
1433     }
1434
1435   // print out not matched tracks from second track array:
1436   listprintpointer = fSecondTrackList;
1437   notmatchedtracknumber = 1;
1438
1439   while(listprintpointer != NULL)
1440     {
1441       if(listprintpointer->fMatchingindicator == 0)
1442         {
1443 #if 0     
1444           HLTInfo("Secondary Track, not matched with original track %d:", notmatchedtracknumber);
1445           HLTInfo("Clusters: %d",listprintpointer->fTrack.GetNHits());
1446           //PYTHIA INFORMATION ABOUT PARTICLE ID        
1447           HLTInfo("First x: %9.6f \t first y: %9.6f \t first z: %9.6f",listprintpointer->fTrack.GetFirstPointX(),listprintpointer->fTrack.GetFirstPointY(),listprintpointer->fTrack.GetFirstPointZ()); 
1448           HLTInfo(" Last x: %9.6f \t  last y: %9.6f \t  last z: %9.6f", listprintpointer->fTrack.GetLastPointX(),listprintpointer->fTrack.GetLastPointY(),listprintpointer->fTrack.GetLastPointZ());
1449           HLTInfo("     Pt: %9.6f \t     Psi: %9.6f \t     Tgl: %9.6f", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl());
1450           HLTInfo("--------------");
1451 #endif
1452           // print these results to file
1453           if(!fDumpFileName.IsNull())
1454             {
1455               fprintf(dumpfile, "Secondary Track, not matched with original track %d: \n", notmatchedtracknumber);
1456               fprintf(dumpfile, "Clusters: %d \n",listprintpointer->fTrack.GetNHits());
1457               fprintf(dumpfile, "First x: %9.6f \t first y: %9.6f \t first z: %9.6f \n",listprintpointer->fTrack.GetFirstPointX(),listprintpointer->fTrack.GetFirstPointY(),listprintpointer->fTrack.GetFirstPointZ()); 
1458               fprintf(dumpfile, " Last x: %9.6f \t  last y: %9.6f \t  last z: %9.6f \n", listprintpointer->fTrack.GetLastPointX(),listprintpointer->fTrack.GetLastPointY(),listprintpointer->fTrack.GetLastPointZ());
1459               fprintf(dumpfile, "     Pt: %9.6f \t     Psi: %9.6f \t     Tgl: %9.6f \n", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl());
1460               fprintf(dumpfile, "--------------\n");
1461             }
1462
1463           ++notmatchedtracknumber;
1464         }
1465
1466       listprintpointer = listprintpointer->fNext;
1467     }
1468
1469   // close output analysis file
1470  if(!fDumpFileName.IsNull())
1471    {
1472      fclose(dumpfile);
1473    };
1474
1475  // if results should be written to graphical output:
1476  if(!fGraphFileName.IsNull())
1477    {
1478      CreateGraphs(1); // specifiy if absolute or rel. differences should be saved (CreateGraphs(0) or CreateGraphs()/ CreateGraphs(1)
1479    };
1480  
1481   // free reserved space
1482   
1483   return 0;
1484 }