3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project *
5 //* ALICE Experiment at CERN, All rights reserved. *
7 //* Primary Authors: J. Wagner <jwagner@cern.ch> *
8 //* for The ALICE HLT Project. *
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 //**************************************************************************
19 /** @file AliHLTTPCCompModelAnalysis.cxx
20 @author J. Wagner jwagner@cern.ch
22 @brief A processing analysis component for the HLT */
28 #include "AliHLTTPCCompModelAnalysis.h"
29 #include "AliHLTTPCTransform.h"
30 #include "AliHLTTPCModelTrack.h"
31 #include "AliHLTTPCCompDataCompressorHelper.h"
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),
44 fFirstTrackList(NULL),
45 fSecondTrackList(NULL),
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),
59 // see header file for class documentation
61 // refer to README to build package
63 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
67 AliHLTTPCCompModelAnalysis::~AliHLTTPCCompModelAnalysis()
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++ )
73 if ( fDiscardedClusters[slice][patch]!=NULL )
75 delete [] fDiscardedClusters[slice][patch];
76 fDiscardedClusters[slice][patch]=NULL;
81 /** initialise arrays for tracks/ discarded clusters depending on model-flags **/
82 Int_t AliHLTTPCCompModelAnalysis::Init()
84 // see header file for class documentation
86 if(fTrackAnalysis) // track quantities are to be initialised
88 fFirstTrackArray.Reset();
89 fSecondTrackArray.Reset();
91 fFirstTrackList = NULL;
92 fSecondTrackList = NULL;
94 fFirstTrashTracks = 0;
95 fSecondTrashTracks = 0;
97 fTotalComparedTracks = 0;
98 fMatchedFirstTrashTracks = 0;
99 fMatchedSecondTrashTracks = 0;
101 fFirstUnmatchedTracks = 0;
102 fSecondUnmatchedTracks = 0;
104 fToleranceDeviation = 0.001;
108 if(fModelAnalysis) // cluster array to be initialised
110 for ( UInt_t slice=0; slice<36; slice++ )
112 for ( UInt_t patch=0; patch<6; patch++ )
114 fDiscardedClusters[slice][patch] = NULL;
118 // initialise trash track list to store discarded tracks
119 fTrackListPointer = NULL;
121 // set all counters to zero:
123 fTotalDiscardedClusters = 0;
124 fValuableDiscardedClusters = 0;
130 Int_t AliHLTTPCCompModelAnalysis::DisplayResults()
132 // see header file for class documentation
133 HLTInfo("--------------------DISPLAYING RESULTS---------------------");
134 // if model loss analysis, then display these results, else display track results
137 DisplayModelResults();
142 DisplayTrackResults();
145 // error message: if no analysis flag is switched on
146 if( (!fModelAnalysis) && (!fTrackAnalysis) )
148 HLTError("Error! Display Results called without any analysis flag switched on.");
155 Int_t AliHLTTPCCompModelAnalysis::SetTracks(AliHLTTPCTrackletData* tracklets, Bool_t fillingfirsttracks)
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
162 if(fillingfirsttracks)
164 HLTDebug( "Reading %u tracks in first array", (unsigned) tracklets->fTrackletCnt );
166 if(tracklets->fTrackletCnt == 0)
168 HLTError("Error! No tracklets to fill into first track array!");
173 fFirstTrackArray.FillTracks(tracklets->fTrackletCnt, tracklets->fTracklets );
179 if(tracklets->fTrackletCnt == 0)
181 HLTError("Error! No tracklets to fill into second track array!");
186 fSecondTrackArray.FillTracks(tracklets->fTrackletCnt, tracklets->fTracklets );
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 );
196 Int_t AliHLTTPCCompModelAnalysis::SetClusters(AliHLTTPCClusterData* clusters, UInt_t slice, UInt_t patch, Bool_t fillingfirstclusters)
198 // see header file for class documentation
199 if(fillingfirstclusters == 1)
201 if ( slice>=36 || patch>=6 )
203 if ( fOriginalClusters[slice][patch] )
205 fOriginalClusters[slice][patch] = clusters;
207 //HLTDebug( "Filling %u clusters in first array", (unsigned)clusters->fSpacePointCnt);
211 if ( slice>=36 || patch>=6 )
213 if ( fSecondaryClusters[slice][patch] )
215 fSecondaryClusters[slice][patch] = clusters;
217 //HLTDebug( "Filling %u clusters in second array", (unsigned)clusters->fSpacePointCnt);
222 Int_t AliHLTTPCCompModelAnalysis::CompareTracks()
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();
229 // error checking: if second array has been filled or not:
232 HLTError("No tracks in first track array!");
237 if(secondtracks == 0)
239 HLTError("No tracks in second track array!");
243 // take track from first tracking,
244 for(Int_t ii=0; ii < firsttracks; ii++)
246 if (fFirstTrackArray.GetCheckedTrack(ii)==NULL) continue;
248 // build track list for all tracks in first array
249 // FIXME: I can't find the cleanup of the linked list fFirstTrackList
250 AliHLTTPCTrackList* currenttrackentry = new AliHLTTPCTrackList;
251 if (!currenttrackentry) return -ENOMEM;
252 currenttrackentry->fTrack = *(fFirstTrackArray.GetCheckedTrack(ii));
254 // get its pythia information,
255 currenttrackentry->fPythiatrack = GetComparableTrackPythiaInfo(currenttrackentry->fTrack);
257 currenttrackentry->fMatchingindicator = 0;
259 // put this element as first in list
260 currenttrackentry->fNext = fFirstTrackList;
261 fFirstTrackList = currenttrackentry;
263 // count tracks below 0.1GeV
264 if(currenttrackentry->fTrack.GetPt()<0.1)
271 // take track from second tracking,
272 for(Int_t ii=0; ii < secondtracks; ii++)
274 if (fSecondTrackArray.GetCheckedTrack(ii)==NULL) continue;
276 // build track list for all tracks in second array
277 // FIXME: I can't find the cleanup of the linked list fSecondTrackArray
278 AliHLTTPCTrackList* currenttrackentry = new AliHLTTPCTrackList;
279 if (!currenttrackentry) return -ENOMEM;
280 currenttrackentry->fTrack = *(fSecondTrackArray.GetCheckedTrack(ii));
282 // get its pythia information,
283 currenttrackentry->fPythiatrack = GetComparableTrackPythiaInfo(currenttrackentry->fTrack);
285 // put this element as first in list
286 currenttrackentry->fNext = fSecondTrackList;
287 fSecondTrackList = currenttrackentry;
289 // count tracks below 0.1GeV
290 if(currenttrackentry->fTrack.GetPt()<0.1)
292 ++fSecondTrashTracks;
297 // search for matching track from secondary tracking
298 AliHLTTPCTrackList* firstmatchpointer = fFirstTrackList;
300 while(firstmatchpointer != NULL)
302 AliHLTTPCTrackList* secondmatchpointer = fSecondTrackList;
304 while(secondmatchpointer != NULL)
307 // compare paramters of the two tracks,
308 // match only when coincidence >= 50% in fToleranceDeviation range!
309 if ((CompareTrackInfo(firstmatchpointer, secondmatchpointer) > 4))
312 if((CompareTrackInfo(firstmatchpointer, secondmatchpointer) > firstmatchpointer->fMatchingindicator))
314 // look if current better matching track has already been matched before
315 if((CompareTrackInfo(firstmatchpointer, secondmatchpointer) > secondmatchpointer->fMatchingindicator))
318 // set previously assigned matchingindicator (if there was one) of secondary track back to zero
319 if(firstmatchpointer->fMatchingindicator > 0)
321 firstmatchpointer->fMatchingtrack->fMatchingindicator = 0;
322 firstmatchpointer->fMatchingtrack->fMatchingtrack = NULL;
325 if(secondmatchpointer->fMatchingindicator > 0)
327 secondmatchpointer->fMatchingtrack->fMatchingindicator = 0;
328 secondmatchpointer->fMatchingtrack->fMatchingtrack = NULL;
331 // compare according to tracks themselves (other possibility: compare pythiatracks - better!)
332 secondmatchpointer->fMatchingindicator = CompareTrackInfo(firstmatchpointer, secondmatchpointer) ;
333 firstmatchpointer->fMatchingindicator = CompareTrackInfo(firstmatchpointer, secondmatchpointer);
335 // remember which track matches which
336 secondmatchpointer->fMatchingtrack = firstmatchpointer;
337 firstmatchpointer->fMatchingtrack = secondmatchpointer;
339 } // end if compare > second matching indicator
341 } // end if compare > first matching indicator
343 }// end if compare > 4
345 secondmatchpointer = secondmatchpointer->fNext;
348 // go on with next original track
349 firstmatchpointer = firstmatchpointer->fNext;
352 // count not matched tracks in first and second track list
353 AliHLTTPCTrackList* nomatchcounter = fFirstTrackList;
355 while(nomatchcounter != NULL)
357 if(nomatchcounter->fMatchingindicator == 0)
359 ++fFirstUnmatchedTracks;
363 ++fTotalComparedTracks;
365 // count matched trash tracks
366 if(nomatchcounter->fTrack.GetPt() < 0.1)
368 ++fMatchedFirstTrashTracks;
371 nomatchcounter = nomatchcounter->fNext;
374 nomatchcounter = fSecondTrackList;
375 while(nomatchcounter != NULL)
377 if(nomatchcounter->fMatchingindicator == 0)
379 ++fSecondUnmatchedTracks;
383 // count matched trash tracks
384 if(nomatchcounter->fTrack.GetPt() < 0.1)
386 ++fMatchedSecondTrashTracks;
390 nomatchcounter = nomatchcounter->fNext;
393 // consistency check: fFirstUnmatchedTracks + fTotalComparedTracks = # of tracks in first array
394 // ...and analogously for second array:
395 if(fFirstUnmatchedTracks + fTotalComparedTracks != firsttracks)
397 HLTWarning("Warning! Possible inconsistency in original track array: Number of compared and unmatched tracks not equal to total number of tracks!");
400 if(fSecondUnmatchedTracks + fTotalComparedTracks != secondtracks)
402 HLTWarning("Warning! Possible inconsistency in second track array: Number of compared and unmatched tracks not equal to total number of tracks!");
408 Int_t AliHLTTPCCompModelAnalysis::CompareClusters(Bool_t relativedifferences)
411 Int_t totaloriginal = 0;
412 Int_t totalsecondary = 0;
413 Int_t usedclusters = 0;
414 Int_t comparedclusters = 0;
415 Int_t notcomparedclusters = 0;
417 // create graphs out of differences and leave loop
418 TFile* clustergraphrootfile = NULL;
419 if(!fGraphFileName.IsNull())
421 clustergraphrootfile = new TFile(fGraphFileName, "recreate");
424 // specifications of histograms
425 Double_t clusterdifffxmin, clusterdifffxmax;
426 Double_t clusterdifffymin, clusterdifffymax;
427 Double_t clusterdifffzmin, clusterdifffzmax;
428 Int_t clusterdifffxbins, clusterdifffybins, clusterdifffzbins;
430 Double_t clusterdifffsigmay2min, clusterdifffsigmay2max;
431 Double_t clusterdifffsigmaz2min, clusterdifffsigmaz2max;
432 Int_t clusterdifffsigmay2bins, clusterdifffsigmaz2bins;
434 if (!relativedifferences) // not tested yet!
436 clusterdifffxmin = -1;
437 clusterdifffxmax = +1;
438 clusterdifffxbins = (Int_t) ((clusterdifffxmax - clusterdifffxmin)/0.0001);
440 clusterdifffymin = -1;
441 clusterdifffymax = +1;
442 clusterdifffybins = (Int_t) ((clusterdifffymax - clusterdifffymin)/0.0001);
444 clusterdifffzmin = -1;
445 clusterdifffzmax = +1;
446 clusterdifffzbins = (Int_t) ((clusterdifffzmax - clusterdifffzmin)/0.0001);
448 clusterdifffsigmay2min = -1;
449 clusterdifffsigmay2max = +1;
450 clusterdifffsigmay2bins = (Int_t) ((clusterdifffsigmay2max - clusterdifffsigmay2min)/0.0001);
452 clusterdifffsigmaz2min = -1;
453 clusterdifffsigmaz2max = +1;
454 clusterdifffsigmaz2bins = (Int_t) ((clusterdifffsigmaz2max - clusterdifffsigmaz2min)/0.0001);
458 clusterdifffxmin = -1;
459 clusterdifffxmax = +1;
460 clusterdifffxbins = (Int_t) ((clusterdifffxmax - clusterdifffxmin)/0.0001);
462 clusterdifffymin = -1;
463 clusterdifffymax = +1;
464 clusterdifffybins = (Int_t) ((clusterdifffymax - clusterdifffymin)/0.0001);
466 clusterdifffzmin = -1;
467 clusterdifffzmax = +1;
468 clusterdifffzbins = (Int_t) ((clusterdifffzmax - clusterdifffzmin)/0.0001);
470 clusterdifffsigmay2min = -1;
471 clusterdifffsigmay2max = +1;
472 clusterdifffsigmay2bins = (Int_t) ((clusterdifffsigmay2max - clusterdifffsigmay2min)/0.0001);
474 clusterdifffsigmaz2min = -1;
475 clusterdifffsigmaz2max = +1;
476 clusterdifffsigmaz2bins = (Int_t) ((clusterdifffsigmaz2max - clusterdifffsigmaz2min)/0.0001);
479 // intialise histogramms
480 TH1F* clusterfxhisto = new TH1F("Differences of x (original - secondary) clusters", "Differences of x (original - secondary) clusters", clusterdifffxbins, clusterdifffxmin, clusterdifffxmax);
481 TH1F* clusterfyhisto = new TH1F("Differences of y (original - secondary) clusters", "Differences of y (original - secondary) clusters", clusterdifffybins, clusterdifffymin, clusterdifffymax);
482 TH1F* clusterfzhisto = new TH1F("Differences of z (original - secondary) clusters", "Differences of z (original - secondary) clusters", clusterdifffzbins, clusterdifffzmin, clusterdifffzmax);
483 TH1F* clusterfsigmay2histo = new TH1F("Differences of sigmay2 (original - secondary) clusters", "Differences of sigmay2 (original - secondary) clusters", clusterdifffsigmay2bins, clusterdifffsigmay2min, clusterdifffsigmay2max);
484 TH1F* clusterfsigmaz2histo = new TH1F("Differences of sigmaz2 (original - secondary) clusters", "Differences of sigmaz2 (original - secondary) clusters", clusterdifffsigmaz2bins, clusterdifffsigmaz2min, clusterdifffsigmaz2max);
487 // see headerfile for class documentation
488 // compare for each slice and patch the clusters in the original cluster array
489 // to the ones of the secondary cluster array
490 for(Int_t slicecntr = 0; slicecntr < 36; slicecntr++)
492 for (Int_t patchcntr = 0; patchcntr < 6; patchcntr++)
494 if(!fOriginalClusters[slicecntr][patchcntr])
496 // HLTDebug("No original clusters for slice %d patch %d", slicecntr, patchcntr);
500 if(!fSecondaryClusters[slicecntr][patchcntr])
502 //HLTDebug("No secondary clusters for slice %d patch %d", slicecntr, patchcntr);
506 for ( unsigned long ii=0; ii<fOriginalClusters[slicecntr][patchcntr]->fSpacePointCnt; ii++ )
510 // search matching secondary cluster by charge and padrow,
511 // fill histograms if cluster has been used in tracking process
514 for(unsigned long jj=0; jj<fSecondaryClusters[slicecntr][patchcntr]->fSpacePointCnt; jj++)
516 // if fTrackN != -1 -> draw histograms out of used clusters
517 // if fTrackN == -1 -> draw histograms out of unused clusters
518 if (fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].GetTrackNumber() == -1)
520 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) )
523 if (relativedifferences == 1)
525 // check whether first entries in cluster array are zero
526 if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX == 0)
528 HLTWarning("Warning! x value of original cluster is zero, relative differences cannot be calculated!");
531 if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY == 0)
533 HLTWarning("Warning! y value of original cluster is zero, relative differences cannot be calculated!");
536 if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ == 0)
538 HLTWarning("Warning! z value of original cluster is zero, relative differences cannot be calculated!");
541 if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2 == 0)
543 HLTWarning("Warning! sigmay2 value of original cluster is zero, relative differences cannot be calculated!");
546 if(fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2 == 0)
548 HLTWarning("Warning! sigmaz2 value of original cluster is zero, relative differences cannot be calculated!");
551 // fill relative differences in histograms
552 clusterfxhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fX)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX,1);
553 clusterfyhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fY)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY,1);
554 clusterfzhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fZ)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ,1);
555 clusterfsigmay2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaY2)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2,1);
556 clusterfsigmaz2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaZ2)/fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2,1);
560 // fill absolute differences histograms
561 clusterfxhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fX - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fX),1);
562 clusterfyhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fY - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fY),1);
563 clusterfzhisto->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fZ - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fZ),1);
564 clusterfsigmay2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaY2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaY2),1);
565 clusterfsigmaz2histo->Fill((fOriginalClusters[slicecntr][patchcntr]->fSpacePoints[ii].fSigmaZ2 - fSecondaryClusters[slicecntr][patchcntr]->fSpacePoints[jj].fSigmaZ2),1);
577 ++notcomparedclusters;
584 // write graphs to rootfile
585 if(!fGraphFileName.IsNull())
587 clustergraphrootfile->WriteObject(clusterfxhisto, "clusterfxhistogram");
588 clustergraphrootfile->WriteObject(clusterfyhisto, "clusterfyhistogram");
589 clustergraphrootfile->WriteObject(clusterfzhisto, "clusterfzhistogram");
590 clustergraphrootfile->WriteObject(clusterfsigmay2histo, "clusterfsigmay2histogram");
591 clustergraphrootfile->WriteObject(clusterfsigmaz2histo, "clusterfsigmaz2histogram");
592 clustergraphrootfile->Close();
595 // count clusters used for tracking
596 for (Int_t slicecount=0; slicecount<36; slicecount++)
598 for(Int_t patchcount=0; patchcount<6; patchcount++)
601 if(!fSecondaryClusters[slicecount][patchcount])
603 //HLTDebug("No secondary clusters for slice %d patch %d", slicecntr, patchcntr);
607 for(Int_t count=0; count < (Int_t) fSecondaryClusters[slicecount][patchcount]->fSpacePointCnt; count++)
612 if(fSecondaryClusters[slicecount][patchcount]->fSpacePoints[count].GetTrackNumber() != -1)
620 // Display results of cluster analysis
621 HLTInfo("Number of original clusters: %d", totaloriginal);
622 HLTInfo("Number of 2ndary clusters: %d", totalsecondary);
623 HLTInfo("Number of 2ndary clusters used for tracking: %d", usedclusters);
624 HLTInfo("Number of compared (tracked) original clusters: %d", comparedclusters);
625 HLTInfo("Number of uncompared (tracked) original clusters: %d", notcomparedclusters);
627 //////////////////////////////////////////////////////
628 FILE* clusterfile = fopen("/afsuser/jwagner/TrackerTest_25092007/cosmics/fullanalysis/clusteranalysis.out", "a");
630 fprintf(clusterfile, "%d \t %d \t %d \t %d \t %d \t %d \n", totaloriginal, totalsecondary, usedclusters, comparedclusters, notcomparedclusters, totaloriginal-comparedclusters-notcomparedclusters);
633 ////////////////////////////////////////////////////////
636 if(comparedclusters + notcomparedclusters != totaloriginal)
638 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);
644 AliHLTTPCTrack AliHLTTPCCompModelAnalysis::GetComparableTrackPythiaInfo(const AliHLTTPCTrack& comparabletrack) const
646 // see headerfile for class documentation
648 AliHLTTPCTrack pythiatrack = comparabletrack;
653 Int_t AliHLTTPCCompModelAnalysis::MarkTrashTrack(AliHLTTPCTrack *lowpttrack)
655 // see header file for class documentation
657 // save track first in lowpttrack list (all lowpttracks are displayed in display function altogether)
658 AliHLTTPCTrackList* tracklistentry = new AliHLTTPCTrackList;
659 tracklistentry->fTrack = *lowpttrack;
660 tracklistentry->fWronglydiscarded = GetTrashTrackPythiaInfo(lowpttrack);
661 tracklistentry->fMatchingindicator = 0; // not needed here, therefore initialised to zero
662 tracklistentry->fNext = fTrackListPointer;
663 tracklistentry->fMatchingtrack = NULL; // not needed here, therefore initialised to NULL
664 fTrackListPointer = tracklistentry;
671 Int_t AliHLTTPCCompModelAnalysis::MarkTrashCluster(AliHLTTPCClusterData *discardedcluster, UInt_t slice, UInt_t patch)
673 // see header file for class documentation
675 // get Pythia information of discarded cluster
676 Bool_t wronglydiscarded = GetClusterPythiaInfo(discardedcluster);
678 // if cluster has been discarded wrongly, save information
682 fDiscardedClusters[slice][patch] = discardedcluster;
683 // increase number of valuable discarded clusters
684 ++fValuableDiscardedClusters;
687 ++fTotalDiscardedClusters;
692 Bool_t AliHLTTPCCompModelAnalysis::GetTrashTrackPythiaInfo(const AliHLTTPCTrack* /*discardedtrack*/ ) const
694 // see header file for class documentation
695 // store information from pythia in current track list entry
696 // fTrackListPointer.pythiatrack = FillFromPythia...
701 Bool_t AliHLTTPCCompModelAnalysis::GetClusterPythiaInfo(const AliHLTTPCClusterData* /*discardedcluster*/) const
703 // see header file for class documentation
704 // Pythia information can be
705 // either: cluster belongs to discarded track with pt < 0.1 Gev (--> cluster correctly discarded)
706 // or: cluster is discarded and does not belong to any pythia track (--> correctly discarded)
707 // or: cluster is discarded but belongs to pythia track (--> cluster WRONGLY discarded!!!)
712 Int_t AliHLTTPCCompModelAnalysis::CompareTrackInfo(AliHLTTPCTrackList* firsttracklistelement, AliHLTTPCTrackList* secondtracklistelement)
714 // see header file for class documentation
715 // calculate matching indicator accoring to the track information
716 // ++matchingindicator for every paramter that matches
718 Int_t currentmatchingindicator = 0;
720 // tolerance range of 1 percent deviation for each quantity
722 // compare start point (x,y,z)
723 if(abs((firsttracklistelement->fTrack.GetFirstPointX() - secondtracklistelement->fTrack.GetFirstPointX()))/firsttracklistelement->fTrack.GetFirstPointX() <= fToleranceDeviation)
724 ++currentmatchingindicator;
726 if(abs((firsttracklistelement->fTrack.GetFirstPointY() - secondtracklistelement->fTrack.GetFirstPointY()))/firsttracklistelement->fTrack.GetFirstPointY() <= fToleranceDeviation)
727 ++currentmatchingindicator;
729 if(abs((firsttracklistelement->fTrack.GetFirstPointZ() - secondtracklistelement->fTrack.GetFirstPointZ()))/firsttracklistelement->fTrack.GetFirstPointZ() <= fToleranceDeviation)
730 ++currentmatchingindicator;
733 if(abs((firsttracklistelement->fTrack.GetLastPointX() - secondtracklistelement->fTrack.GetLastPointX()))/firsttracklistelement->fTrack.GetLastPointX() <= fToleranceDeviation)
734 ++currentmatchingindicator;
736 if(abs((firsttracklistelement->fTrack.GetLastPointY() - secondtracklistelement->fTrack.GetLastPointY()))/firsttracklistelement->fTrack.GetLastPointY() <= fToleranceDeviation)
737 ++currentmatchingindicator;
739 if(abs((firsttracklistelement->fTrack.GetLastPointZ() - secondtracklistelement->fTrack.GetLastPointZ()))/firsttracklistelement->fTrack.GetLastPointZ() <= fToleranceDeviation)
740 ++currentmatchingindicator;
742 // compare pt, psi, tgl
743 if(abs((firsttracklistelement->fTrack.GetPt() - secondtracklistelement->fTrack.GetPt()))/firsttracklistelement->fTrack.GetPt() <= fToleranceDeviation)
744 ++currentmatchingindicator;
746 if(abs((firsttracklistelement->fTrack.GetPsi() - secondtracklistelement->fTrack.GetPsi()))/firsttracklistelement->fTrack.GetPsi() <= fToleranceDeviation)
747 ++currentmatchingindicator;
749 if(abs((firsttracklistelement->fTrack.GetTgl() - secondtracklistelement->fTrack.GetTgl()))/firsttracklistelement->fTrack.GetTgl() <= fToleranceDeviation)
750 ++currentmatchingindicator;
752 // compare number of assigned cluster hits
753 if(abs((firsttracklistelement->fTrack.GetNHits() - secondtracklistelement->fTrack.GetNHits()))/firsttracklistelement->fTrack.GetNHits() <= fToleranceDeviation)
754 ++currentmatchingindicator;
756 return currentmatchingindicator;
759 Int_t AliHLTTPCCompModelAnalysis::ComparePythiaTrackInfo(AliHLTTPCTrackList* firsttracklistelement, AliHLTTPCTrackList* secondtracklistelement)
761 // see header file for class documentation
762 // calculate matching indicator accoring to the track information
763 // ++matchingindicator for every paramter that matches
765 Int_t currentmatchingindicator = 0;
767 // tolerance range of 1 percent deviation for each quantity
769 // compare start point (x,y,z)
770 if(firsttracklistelement->fPythiatrack.GetFirstPointX() == secondtracklistelement->fPythiatrack.GetFirstPointX())
771 ++currentmatchingindicator;
773 if(firsttracklistelement->fPythiatrack.GetFirstPointY() == secondtracklistelement->fPythiatrack.GetFirstPointY())
774 ++currentmatchingindicator;
776 if(firsttracklistelement->fPythiatrack.GetFirstPointZ() == secondtracklistelement->fPythiatrack.GetFirstPointZ())
777 ++currentmatchingindicator;
780 if(firsttracklistelement->fPythiatrack.GetLastPointX() == secondtracklistelement->fPythiatrack.GetLastPointX())
781 ++currentmatchingindicator;
783 if(firsttracklistelement->fPythiatrack.GetLastPointY() == secondtracklistelement->fPythiatrack.GetLastPointY())
784 ++currentmatchingindicator;
786 if(firsttracklistelement->fPythiatrack.GetLastPointZ() == secondtracklistelement->fPythiatrack.GetLastPointZ())
787 ++currentmatchingindicator;
789 // compare pt, psi, tgl
790 if(firsttracklistelement->fPythiatrack.GetPt() == secondtracklistelement->fPythiatrack.GetPt())
791 ++currentmatchingindicator;
793 if(firsttracklistelement->fPythiatrack.GetPsi() == secondtracklistelement->fPythiatrack.GetPsi())
794 ++currentmatchingindicator;
796 if(firsttracklistelement->fPythiatrack.GetTgl() == secondtracklistelement->fPythiatrack.GetTgl())
797 ++currentmatchingindicator;
799 // compare number of assigned cluster hits
800 if(firsttracklistelement->fPythiatrack.GetNHits() == secondtracklistelement->fPythiatrack.GetNHits())
801 ++currentmatchingindicator;
803 return currentmatchingindicator;
806 Int_t AliHLTTPCCompModelAnalysis::CreateGraphs(Bool_t relativedifferences)
808 // see header file for class documentation
809 AliHLTTPCTrackList* tracklistpointer = fFirstTrackList;
811 AliHLTTPCTrackList* trackmatchingpointer;
813 // set up histogram ranges
814 Double_t difffirstxmin, difffirstxmax;
815 Double_t difffirstymin, difffirstymax;
816 Double_t difffirstzmin, difffirstzmax;
817 Int_t difffirstxbins, difffirstybins, difffirstzbins;
819 Double_t difflastxmin, difflastxmax;
820 Double_t difflastymin, difflastymax;
821 Double_t difflastzmin, difflastzmax;
822 Int_t difflastxbins, difflastybins, difflastzbins;
824 Double_t diffptmin, diffptmax;
825 Double_t diffpsimin, diffpsimax;
826 Double_t difftglmin, difftglmax;
827 Int_t diffptbins, diffpsibins, difftglbins;
829 Double_t diffclustermin, diffclustermax;
830 Int_t diffclusterbins;
832 // resolution histograms (currently not working since pterr = 0 for every track!)
833 Double_t diffpterrmin, diffpterrmax;
834 Double_t diffpsierrmin, diffpsierrmax;
835 Double_t difftglerrmin, difftglerrmax;
836 Int_t diffpterrbins, diffpsierrbins, difftglerrbins;
838 if(!relativedifferences)
842 difffirstxbins = (Int_t) ((difffirstxmax - difffirstxmin)/0.0001);
846 difffirstybins = (Int_t) ((difffirstymax - difffirstymin)/0.0001);
850 difffirstzbins = (Int_t) ((difffirstzmax - difffirstzmin)/0.0001);
854 difflastxbins = (Int_t) ((difflastxmax - difflastxmin)/0.0001);
858 difflastybins = (Int_t) ((difflastymax - difflastymin)/0.0001);
862 difflastzbins = (Int_t) ((difflastzmax - difflastzmin)/0.0001);
866 diffptbins = (Int_t) ((diffptmax - diffptmin)/0.0001);
870 diffpsibins = (Int_t) ((diffpsimax - diffpsimin)/0.0001);
874 difftglbins = (Int_t) ((difftglmax - difftglmin)/0.0001);
876 diffclustermin = -50;
877 diffclustermax = +50;
878 diffclusterbins = (Int_t) ((diffclustermax - diffclustermin)/1);
883 diffpterrbins = (Int_t) ((diffpterrmax - diffpterrmin)/1);
887 diffpsierrbins = (Int_t) ((diffpsierrmax - diffpsierrmin)/1);
891 difftglerrbins = (Int_t) ((difftglerrmax - difftglerrmin)/1);
899 difffirstxbins = (Int_t) ((difffirstxmax - difffirstxmin)/0.0001);
903 difffirstybins = (Int_t) ((difffirstymax - difffirstymin)/0.0001);
907 difffirstzbins = (Int_t) ((difffirstzmax - difffirstzmin)/0.0001);
911 difflastxbins = (Int_t) ((difflastxmax - difflastxmin)/0.0001);
915 difflastybins = (Int_t) ((difflastymax - difflastymin)/0.0001);
919 difflastzbins = (Int_t) ((difflastzmax - difflastzmin)/0.0001);
923 diffptbins = (Int_t) ((diffptmax - diffptmin)/0.0001);
927 diffpsibins = (Int_t) ((diffpsimax - diffpsimin)/0.0001);
931 difftglbins = (Int_t) ((difftglmax - difftglmin)/0.0001);
935 diffclusterbins = (Int_t) ((diffclustermax - diffclustermin)/0.0001);
940 diffpterrbins = (Int_t) ((diffpterrmax - diffpterrmin)/0.0001);
944 diffpsierrbins = (Int_t) ((diffpsierrmax - diffpsierrmin)/0.0001);
948 difftglerrbins = (Int_t) ((difftglerrmax - difftglerrmin)/0.0001);
952 // intialise histogramms
953 TH1F* firstxhisto = new TH1F("Differences of first x (original - secondary) track", "Differences of first x (original - secondary) track", difffirstxbins, difffirstxmin, difffirstxmax);
954 TH1F* firstyhisto = new TH1F("Differences of first y (original - secondary) track", "Differences of first y (original - secondary) track", difffirstybins, difffirstymin, difffirstymax);
955 TH1F* firstzhisto = new TH1F("Differences of first z (original - secondary) track", "Differences of first z (original - secondary) track", difffirstzbins, difffirstzmin, difffirstzmax);
956 TH1F* lastxhisto = new TH1F("Differences of last x (original - secondary) track", "Differences of last x (original - secondary) track", difflastxbins, difflastxmin, difflastxmax);
957 TH1F* lastyhisto = new TH1F("Differences of last y (original - secondary) track", "Differences of last y (original - secondary) track", difflastybins, difflastymin, difflastymax);
958 TH1F* lastzhisto = new TH1F("Differences of last z (original - secondary) track", "Differences of last z (original - secondary) track", difflastzbins, difflastzmin, difflastzmax);
959 TH1F* pthisto = new TH1F("Differences of pt (original - secondary) track", "Differences of pt (original - secondary) track", diffptbins, diffptmin, diffptmax);
960 TH1F* psihisto = new TH1F("Differences of psi (original - secondary) track", "Differences of psi (original - secondary) track", diffpsibins, diffpsimin, diffpsimax);
961 TH1F* tglhisto = new TH1F("Differences of tgl (original - secondary) track", "Differences of tgl (original - secondary) track", difftglbins, difftglmin, difftglmax);
962 TH1F* clusterhisto = new TH1F("Differences of asserted clusters (original - secondary) track", "Differences of asserted clusters (original - secondary) track", diffclusterbins, diffclustermin, diffclustermax);
965 // commented out since pterr is zero for every track!
966 TH1F* pterrhisto = new TH1F("Differences of pt error (original - secondary) track", "Differences of pt error (original - secondary) track", diffpterrbins, diffpterrmin, diffpterrmax);
967 TH1F* psierrhisto = new TH1F("Differences of psi error (original - secondary) track", "Differences of psi error (original - secondary) track", diffpsierrbins, diffpsierrmin, diffpsierrmax);
968 TH1F* tglerrhisto = new TH1F("Differences of tgl error (original - secondary) track", "Differences of tgl error (original - secondary) track", difftglerrbins, difftglerrmin, difftglerrmax);
972 // initialise histograms for tracking efficiency against pt
973 // relative p_t resolution: 1.2% -> take 0.1GeV * 1.2 % --> binsize 0.001 sufficient to grant correct resolution for low pt
974 TH1F* firsttracks = new TH1F("pt occurrence original", "Occurrence of pt in original tracks", 10000, 0, 10);
975 TH1F* matchedtrackeff = new TH1F("matchedtreffeff", "Occurrence of 2ndary tracks with good pt", 10000, 0, 10);
976 TH1F* trackeff = new TH1F("tracking efficiency vs. pt", "Tracking efficiency vs. pt", 10000, 0, 10);
978 // evaluate quality of fit:
979 TH1I* matchinghisto = new TH1I("Matching indicator (5 - 10)", "Matching indicator (5 - 10)", 11, 0, 11);
981 while(tracklistpointer != NULL)
983 // if currently processed track is matched, store differences of their parameters in histogram
984 if(tracklistpointer->fMatchingindicator > 0)
987 trackmatchingpointer = tracklistpointer->fMatchingtrack;
989 // fill histograms for trackingefficiency vs. pt
990 firsttracks->Fill(tracklistpointer->fTrack.GetPt(),1);
992 if(abs(tracklistpointer->fTrack.GetPt()-trackmatchingpointer->fTrack.GetPt()) < 0.012*tracklistpointer->fTrack.GetPt())
994 matchedtrackeff->Fill(tracklistpointer->fTrack.GetPt(),1);
997 //tracklistpointer = tracklistpointer->fNext; // only efficiency is considered...
998 //continue; // only efficiency is considered...
1000 if(relativedifferences == 1) // fill histogram with relative differences
1003 // check if first track parameters are not zero!
1004 if (tracklistpointer->fTrack.GetFirstPointX()==0)
1006 HLTWarning("Warning! First x of original track is zero, relative differences cannot be calculated!");
1008 if (tracklistpointer->fTrack.GetFirstPointY()==0)
1010 HLTWarning("Warning! First y of original track is zero, relative differences cannot be calculated!");
1012 if (tracklistpointer->fTrack.GetFirstPointZ()==0)
1014 HLTWarning("Warning! First z of original track is zero, relative differences cannot be calculated!");
1016 if (tracklistpointer->fTrack.GetLastPointX()==0)
1018 HLTWarning("Warning! Last x of original track is zero, relative differences cannot be calculated!");
1020 if (tracklistpointer->fTrack.GetLastPointY()==0)
1022 HLTWarning("Warning! Last y of original track is zero, relative differences cannot be calculated!");
1024 if (tracklistpointer->fTrack.GetLastPointZ()==0)
1026 HLTWarning("Warning! Last z of original track is zero, relative differences cannot be calculated!");
1028 if (tracklistpointer->fTrack.GetPt()==0)
1030 HLTWarning("Warning! Pt of original track is zero, relative differences cannot be calculated!");
1032 if (tracklistpointer->fTrack.GetPsi()==0)
1034 HLTWarning("Warning! Psi of original track is zero, relative differences cannot be calculated!");
1036 if (tracklistpointer->fTrack.GetTgl()==0)
1038 HLTWarning("Warning! Tgl of original track is zero, relative differences cannot be calculated!");
1040 firstxhisto->Fill((tracklistpointer->fTrack.GetFirstPointX()-trackmatchingpointer->fTrack.GetFirstPointX())/tracklistpointer->fTrack.GetFirstPointX(),1);
1041 firstyhisto->Fill((tracklistpointer->fTrack.GetFirstPointY()-trackmatchingpointer->fTrack.GetFirstPointY())/tracklistpointer->fTrack.GetFirstPointY(),1);
1042 firstzhisto->Fill((tracklistpointer->fTrack.GetFirstPointZ()-trackmatchingpointer->fTrack.GetFirstPointZ())/tracklistpointer->fTrack.GetFirstPointZ(),1);
1043 lastxhisto->Fill((tracklistpointer->fTrack.GetLastPointX()-trackmatchingpointer->fTrack.GetLastPointX())/tracklistpointer->fTrack.GetLastPointX(),1);
1044 lastyhisto->Fill((tracklistpointer->fTrack.GetLastPointY()-trackmatchingpointer->fTrack.GetLastPointY())/tracklistpointer->fTrack.GetLastPointY(),1);
1045 lastzhisto->Fill((tracklistpointer->fTrack.GetLastPointZ()-trackmatchingpointer->fTrack.GetLastPointZ())/tracklistpointer->fTrack.GetLastPointZ(),1);
1046 pthisto->Fill((tracklistpointer->fTrack.GetPt()-trackmatchingpointer->fTrack.GetPt())/tracklistpointer->fTrack.GetPt(),1);
1047 psihisto->Fill((tracklistpointer->fTrack.GetPsi()-trackmatchingpointer->fTrack.GetPsi())/tracklistpointer->fTrack.GetPsi(),1);
1048 tglhisto->Fill((tracklistpointer->fTrack.GetTgl()-trackmatchingpointer->fTrack.GetTgl())/tracklistpointer->fTrack.GetTgl(),1);
1049 clusterhisto->Fill((tracklistpointer->fTrack.GetNHits()-trackmatchingpointer->fTrack.GetNHits())/tracklistpointer->fTrack.GetNHits(),1);
1052 pterrhisto->Fill((tracklistpointer->fTrack.GetPterr()-trackmatchingpointer->fTrack.GetPterr())/tracklistpointer->fTrack.GetPterr(),1);
1053 psierrhisto->Fill((tracklistpointer->fTrack.GetPsierr()-trackmatchingpointer->fTrack.GetPsierr())/tracklistpointer->fTrack.GetPsierr(),1);
1054 tglerrhisto->Fill((tracklistpointer->fTrack.GetTglerr()-trackmatchingpointer->fTrack.GetTglerr())/tracklistpointer->fTrack.GetTglerr(),1);
1056 //HLTInfo("Pterr: 1st:%f 2nd:%f value:%f",tracklistpointer->fTrack.GetPterr(),trackmatchingpointer->fTrack.GetPterr(),(tracklistpointer->fTrack.GetPterr()-trackmatchingpointer->fTrack.GetPterr())/tracklistpointer->fTrack.GetPterr());
1057 //HLTInfo("Psierr: 1st:%f 2nd:%f value:%f",tracklistpointer->fTrack.GetPsierr(),trackmatchingpointer->fTrack.GetPsierr(),(tracklistpointer->fTrack.GetPsierr()-trackmatchingpointer->fTrack.GetPsierr())/tracklistpointer->fTrack.GetPsierr());
1058 //HLTInfo("Tglerr: 1st:%f 2nd:%f value:%f",tracklistpointer->fTrack.GetTglerr(),trackmatchingpointer->fTrack.GetTglerr(),(tracklistpointer->fTrack.GetTglerr()-trackmatchingpointer->fTrack.GetTglerr())/tracklistpointer->fTrack.GetTglerr());
1061 else // otherwise fill histogram with absolute differences
1063 firstxhisto->Fill(tracklistpointer->fTrack.GetFirstPointX()-trackmatchingpointer->fTrack.GetFirstPointX(),1);
1064 firstyhisto->Fill(tracklistpointer->fTrack.GetFirstPointY()-trackmatchingpointer->fTrack.GetFirstPointY(),1);
1065 firstzhisto->Fill(tracklistpointer->fTrack.GetFirstPointZ()-trackmatchingpointer->fTrack.GetFirstPointZ(),1);
1066 lastxhisto->Fill(tracklistpointer->fTrack.GetLastPointX()-trackmatchingpointer->fTrack.GetLastPointX(),1);
1067 lastyhisto->Fill(tracklistpointer->fTrack.GetLastPointY()-trackmatchingpointer->fTrack.GetLastPointY(),1);
1068 lastzhisto->Fill(tracklistpointer->fTrack.GetLastPointZ()-trackmatchingpointer->fTrack.GetLastPointZ(),1);
1069 pthisto->Fill(tracklistpointer->fTrack.GetPt()-trackmatchingpointer->fTrack.GetPt(),1);
1070 psihisto->Fill(tracklistpointer->fTrack.GetPsi()-trackmatchingpointer->fTrack.GetPsi(),1);
1071 tglhisto->Fill(tracklistpointer->fTrack.GetTgl()-trackmatchingpointer->fTrack.GetTgl(),1);
1072 clusterhisto->Fill(tracklistpointer->fTrack.GetNHits()-trackmatchingpointer->fTrack.GetNHits(),1);
1074 // commented out since pterr is always zero for every track!
1075 pterrhisto->Fill(tracklistpointer->fTrack.GetPterr()-trackmatchingpointer->fTrack.GetPterr(),1);
1076 psierrhisto->Fill(tracklistpointer->fTrack.GetPsierr()-trackmatchingpointer->fTrack.GetPsierr(),1);
1077 tglerrhisto->Fill(tracklistpointer->fTrack.GetTglerr()-trackmatchingpointer->fTrack.GetTglerr(),1);
1081 // fill histogram that determines the quality of the fit
1082 matchinghisto->Fill(tracklistpointer->fMatchingindicator,1);
1085 tracklistpointer = tracklistpointer->fNext;
1088 trackeff->Divide(matchedtrackeff, firsttracks,1,1,"");
1090 // write histograms to root file specified in command line argument -graphs <filename>.ROOT
1091 if(!fGraphFileName.IsNull())
1093 TFile* graphrootfile = new TFile(fGraphFileName, "update");
1094 graphrootfile->WriteObject(firstxhisto,"firstxhistogram");
1095 //firstxhisto->Write();
1096 graphrootfile->WriteObject(firstyhisto,"firstyhistogram");
1097 //firstyhisto->Write();
1098 graphrootfile->WriteObject(firstzhisto,"firstzhistogram");
1099 //firstzhisto->Write();
1100 graphrootfile->WriteObject(lastxhisto,"lastxhistogram");
1101 //lastxhisto->Write();
1102 graphrootfile->WriteObject(lastyhisto,"lastyhistogram");
1103 //lastyhisto->Write();
1104 graphrootfile->WriteObject(lastzhisto,"lastzhistogram");
1105 //lastzhisto->Write();
1106 graphrootfile->WriteObject(pthisto,"pthistogram");
1108 graphrootfile->WriteObject(psihisto,"psihistogram");
1109 //psihisto->Write();
1110 graphrootfile->WriteObject(tglhisto,"tglhistogram");
1111 //tglhisto->Write();
1112 graphrootfile->WriteObject(clusterhisto,"clusterhistogram");
1113 //clusterhisto->Write();
1114 graphrootfile->WriteObject(matchinghisto,"matchinghistogram");
1115 //matchinghisto->Write();
1116 graphrootfile->WriteObject(firsttracks, "firsttrackeff");
1117 graphrootfile->WriteObject(matchedtrackeff, "secondtrackeff");
1118 graphrootfile->WriteObject(trackeff, "trackingefficiency");
1120 // errors in tracking (commented out since pterr is always 0 for every track!)
1121 graphrootfile->WriteObject(pterrhisto, "pterrhistogram");
1122 //pterrhisto->Write();
1123 graphrootfile->WriteObject(psierrhisto, "psierrhistogram");
1124 //psierrhisto->Write();
1125 graphrootfile->WriteObject(tglerrhisto, "tglerrhistogram");
1126 //tglerrhisto->Write();
1127 graphrootfile->Close();
1131 HLTError("Error! No file for graphical output specified.");
1138 Int_t AliHLTTPCCompModelAnalysis::DisplayModelResults()
1140 // see header file for class documentation
1142 //print out parameters of discarded track:
1143 AliHLTTPCTrackList* trackprintpointer = fTrackListPointer;
1144 AliHLTTPCTrackList* tracklistdeleter= trackprintpointer;
1146 FILE* dumpfile = NULL;
1148 if(!fDumpFileName.IsNull())
1150 // open new file specified by command line argument
1151 dumpfile = fopen(fDumpFileName.Data(),"a");
1152 fprintf(dumpfile,"---------------MODEL ANALYSIS--------------- \n");
1153 fprintf(dumpfile,"---------------DISCARDED TRACKS: %d --------------- \n", fTrashTracks);
1156 if(fTrackListPointer == NULL)
1158 HLTInfo("No tracks discarded");
1159 HLTInfo("--------------");
1164 HLTInfo("---------------DISCARDED TRACKS: %d ---------------", fTrashTracks);
1166 Int_t trashtrackcounter = 1;
1168 while(trackprintpointer != NULL)
1171 // infos about found track
1172 HLTInfo("%d : Discarding track with %d clusters.", trashtrackcounter, trackprintpointer->fTrack.GetNHits());
1173 //PYTHIA INFORMATION ABOUT PARTICLE ID
1174 HLTInfo("Track parameters of discarded track:");
1175 HLTInfo("First x: %9.6f \t first y: %9.6f \t first z: %9.6f",trackprintpointer->fTrack.GetFirstPointX(),trackprintpointer->fTrack.GetFirstPointY(),trackprintpointer->fTrack.GetFirstPointZ());
1176 HLTInfo(" Last x: %9.6f \t last y: %9.6f \t last z: %9.6f", trackprintpointer->fTrack.GetLastPointX(),trackprintpointer->fTrack.GetLastPointY(),trackprintpointer->fTrack.GetLastPointZ());
1177 HLTInfo(" Pt: %9.6f \t Psi: %9.6f \t Tgl: %9.6f", trackprintpointer->fTrack.GetPt(), trackprintpointer->fTrack.GetPsi(), trackprintpointer->fTrack.GetTgl());
1179 // write results to file if specified by command line argument
1180 if(!fDumpFileName.IsNull())
1182 fprintf(dumpfile,"%d : Discarding track with %d clusters. \n", trashtrackcounter, trackprintpointer->fTrack.GetNHits());
1183 fprintf(dumpfile,"Track parameters of discarded track: \n");
1184 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());
1185 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());
1186 fprintf(dumpfile," Pt: %9.6f \t Psi: %9.6f \t Tgl: %9.6f \n", trackprintpointer->fTrack.GetPt(), trackprintpointer->fTrack.GetPsi(), trackprintpointer->fTrack.GetTgl());
1189 // comparison with pythia information
1190 if(trackprintpointer->fWronglydiscarded)
1192 HLTInfo("Found track has been wrongly discarded according to Pythia information.");
1194 // write results to file if specified by command line argument
1195 if(!fDumpFileName.IsNull())
1197 fprintf(dumpfile,"Found track has been wrongly discarded accoring to Pythia information. \n");
1202 // found pt must be in range pt_pythia \pm 10% to be accepted
1203 Double_t ptmin = trackprintpointer->fPythiatrack.GetPt() - 0.1*trackprintpointer->fPythiatrack.GetPt();
1204 Double_t ptmax = trackprintpointer->fPythiatrack.GetPt() + 0.1*trackprintpointer->fPythiatrack.GetPt();
1206 if( (trackprintpointer->fTrack.GetPt() < ptmin) ||(trackprintpointer->fTrack.GetPt() > ptmax) )
1208 HLTInfo("Pt of found track %f differs more than 10 %% from pt of pythia track %f.",trackprintpointer->fTrack.GetPt(), trackprintpointer->fPythiatrack.GetPt());
1210 if(!fDumpFileName.IsNull())
1212 // write result to file if specified by command line argument
1213 fprintf(dumpfile,"Pt of found track %f differs more than 10 %% from pt of pythia track %f. \n",trackprintpointer->fTrack.GetPt(), trackprintpointer->fPythiatrack.GetPt());
1218 HLTInfo("--------------");
1220 if(!fDumpFileName.IsNull())
1222 fprintf(dumpfile,"-------------- \n");
1226 // go to next element in trash track list
1227 tracklistdeleter = trackprintpointer;
1228 trackprintpointer = trackprintpointer->fNext;
1230 ++trashtrackcounter;
1233 delete tracklistdeleter;
1234 tracklistdeleter = NULL;
1236 } // end of while(trackpointer != NULL)
1240 // print out number of noise clusters (not assigned to any track and not valid)
1241 HLTInfo("Number of discarded clusters not assigned to any track: %d", fTotalDiscardedClusters);
1242 HLTInfo("--------------");
1244 // write results to file if specified by command line argument
1245 if(!fDumpFileName.IsNull())
1247 fprintf(dumpfile,"Number of discarded clusters not assigned to any track: %d \n", fTotalDiscardedClusters);
1248 fprintf(dumpfile,"-------------- \n");
1251 // print out paramters of discarded valuable clusters
1252 HLTInfo("Number of discarded VALUABLE clusters: %d", fValuableDiscardedClusters);
1254 HLTInfo("--------------");
1256 if(!fDumpFileName.IsNull())
1259 fprintf(dumpfile,"Number of discarded VALUABLE clusters: %d \n", fValuableDiscardedClusters);
1266 Int_t AliHLTTPCCompModelAnalysis::DisplayTrackResults()
1268 // see header file for class documentation
1269 HLTInfo("---------------CLUSTER ANALYSIS---------------");
1270 if(CompareClusters(1) != 0)
1275 // start with comparison
1276 if(CompareTracks() != 0)
1281 // if dumptofile is activated, append results to output analysis file
1282 FILE* dumpfile = NULL;
1284 if(!fDumpFileName.IsNull())
1286 // open new file specified by command line argument
1287 dumpfile = fopen(fDumpFileName.Data(),"a");
1289 fprintf(dumpfile,"---------------TRACK ANALYSIS--------------- \n");
1293 // print out number of compared tracks
1294 HLTInfo("---------------TRACK ANALYSIS---------------");
1295 HLTInfo("---------------ORIGINAL TRACKS: %d ---------------", fFirstTrackArray.GetNTracks());
1296 HLTInfo("Number of tracks with pt < 0.1 GeV: %d", fFirstTrashTracks);
1297 HLTInfo("Number of matched tracks with pt < 0.1 GeV: %d", fMatchedFirstTrashTracks);
1298 //PYTHIA INFORMATION ABOUT PARTICLE IDs
1299 HLTInfo("---------------2NDARY TRACKS: %d ---------------", fSecondTrackArray.GetNTracks());
1300 HLTInfo("Number of tracks with pt < 0.1 GeV: %d", fSecondTrashTracks);
1301 HLTInfo("Number of matched tracks with pt < 0.1 GeV: %d", fMatchedSecondTrashTracks);
1302 //PYTHIA INFORMATION ABOUT PARTICLE IDs
1303 HLTInfo("--------------");
1304 HLTInfo("Comparison of tracks within parameter precision: %f", fToleranceDeviation);
1305 HLTInfo("Number of compared tracks: %d", fTotalComparedTracks);
1306 HLTInfo("Number of unmatched original tracks: %d", fFirstUnmatchedTracks);
1307 HLTInfo("Number of unmatched secondary tracks: %d", fSecondUnmatchedTracks);
1308 HLTInfo("--------------");
1310 // print results to file
1311 if(!fDumpFileName.IsNull())
1313 fprintf(dumpfile, "---------------%d ORIGINAL TRACKS---------------\n", fFirstTrackArray.GetNTracks());
1314 fprintf(dumpfile,"Number of tracks with pt < 0.1 GeV: %d \n", fFirstTrashTracks);
1315 fprintf(dumpfile,"Number of matched tracks with pt < 0.1 GeV: %d \n", fMatchedFirstTrashTracks);
1316 fprintf(dumpfile,"---------------%d 2NDARY TRACKS---------------\n", fSecondTrackArray.GetNTracks());
1317 fprintf(dumpfile,"Number of tracks with pt < 0.1 GeV: %d \n", fSecondTrashTracks);
1318 fprintf(dumpfile,"Number of matched tracks with pt < 0.1 GeV: %d \n", fMatchedSecondTrashTracks);
1319 fprintf(dumpfile,"--------------\n");
1320 fprintf(dumpfile,"Comparison of tracks within parameter precision: %f \n", fToleranceDeviation);
1321 fprintf(dumpfile,"Number of compared tracks: %d \n", fTotalComparedTracks);
1322 fprintf(dumpfile,"Number of unmatched original tracks: %d \n", fFirstUnmatchedTracks);
1323 fprintf(dumpfile,"Number of unmatched secondary tracks: %d \n", fSecondUnmatchedTracks);
1324 fprintf(dumpfile,"--------------\n");
1327 //////////////////////////////////////////////////////////////////////
1328 // additional files temporarily necessary for output information
1329 FILE* infofile = fopen("/afsuser/jwagner/TrackerTest_25092007/pp-ca/fullanalysis08012008/trackingefficiency.out","a");
1330 FILE* info2file = fopen("/afsuser/jwagner/TrackerTest_25092007/pp-ca/fullanalysis08012008/parameters.out", "a");
1332 // consistent = 0 if tracks and second tracks match perfectly, i.e. no doubly assigned tracks,etc.
1333 Int_t consistentfirst = fFirstTrackArray.GetNTracks() - fTotalComparedTracks - fFirstUnmatchedTracks;
1334 Int_t consistentsecond = fSecondTrackArray.GetNTracks() - fTotalComparedTracks - fSecondUnmatchedTracks;
1336 //fprintf(infofile, "1st tracks, 2nd tracks, compared, 1st uncompared, 2nd uncompared, cons.1st, cons.2nd \n");
1337 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);
1339 //fprintf(infofile, "1st trash tracks, 2nd trash tracks, 1st matched trash tracks, 2nd matched trash tracks \n");
1340 fprintf(infofile, "%d \t %d \t %d \t %d \n", fFirstTrashTracks, fSecondTrashTracks, fMatchedFirstTrashTracks, fMatchedSecondTrashTracks);
1344 ////////////////////////////////////////////////////////////////////
1346 // print out deviations
1347 Int_t tracknumber = 1;
1348 AliHLTTPCTrackList* listprintpointer = fFirstTrackList;
1349 while(listprintpointer != NULL)
1351 // print out parameters of original track in comparison to secondary track:
1353 if(listprintpointer->fMatchingindicator != 0)
1357 HLTInfo("Track %d:", tracknumber);
1358 HLTInfo("Original track matched to secondary track with matchingindicator %d", listprintpointer->fMatchingindicator);
1360 HLTInfo("Parameter comparison: Original vs. Secondary");
1362 HLTInfo("Clusters: %d \t %d ",listprintpointer->fTrack.GetNHits(), listprintpointer->fMatchingtrack->fTrack.GetNHits());
1363 HLTInfo("First x: %9.6f \t %9.6f", listprintpointer->fTrack.GetFirstPointX(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointX());
1364 HLTInfo("First y: %9.6f \t %9.6f", listprintpointer->fTrack.GetFirstPointY(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointY());
1365 HLTInfo("First z: %9.6f \t %9.6f", listprintpointer->fTrack.GetFirstPointZ(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointZ());
1366 HLTInfo("Last x: %9.6f \t %9.6f", listprintpointer->fTrack.GetLastPointX(), listprintpointer->fMatchingtrack->fTrack.GetLastPointX());
1367 HLTInfo("Last y: %9.6f \t %9.6f", listprintpointer->fTrack.GetLastPointY(), listprintpointer->fMatchingtrack->fTrack.GetLastPointY());
1368 HLTInfo("Last z: %9.6f \t %9.6f", listprintpointer->fTrack.GetLastPointZ(), listprintpointer->fMatchingtrack->fTrack.GetLastPointZ());
1369 HLTInfo(" Pt: %9.6f \t %9.6f", listprintpointer->fTrack.GetPt(), listprintpointer->fMatchingtrack->fTrack.GetPt());
1370 HLTInfo(" Psi: %9.6f \t %9.6f", listprintpointer->fTrack.GetPsi(), listprintpointer->fMatchingtrack->fTrack.GetPsi());
1371 HLTInfo(" Tgl: %9.6f \t %9.6f", listprintpointer->fTrack.GetTgl(), listprintpointer->fMatchingtrack->fTrack.GetTgl());
1373 HLTInfo(" Pterr: %9.6f \t %9.6f", listprintpointer->fTrack.GetPterr(), listprintpointer->fMatchingtrack->fTrack.GetPterr());
1374 HLTInfo("Psierr: %9.6f \t %9.6f", listprintpointer->fTrack.GetPsierr(), listprintpointer->fMatchingtrack->fTrack.GetPsierr());
1375 HLTInfo("Tglerr: %9.6f \t %9.6f", listprintpointer->fTrack.GetTglerr(), listprintpointer->fMatchingtrack->fTrack.GetTglerr());
1377 HLTInfo("--------------");
1379 // print these results to file
1380 if(!fDumpFileName.IsNull())
1382 fprintf(dumpfile, "Track %d: \n", tracknumber);
1383 fprintf(dumpfile, "Original track matched to secondary track with matchingindicator %d \n", listprintpointer->fMatchingindicator);
1384 fprintf(dumpfile, "Parameter comparison: Original vs. Secondary \n");
1385 fprintf(dumpfile, "Clusters: %d \t %d \n ",listprintpointer->fTrack.GetNHits(), listprintpointer->fMatchingtrack->fTrack.GetNHits());
1386 fprintf(dumpfile, "First x: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetFirstPointX(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointX());
1387 fprintf(dumpfile, "First y: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetFirstPointY(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointY());
1388 fprintf(dumpfile, "First z: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetFirstPointZ(), listprintpointer->fMatchingtrack->fTrack.GetFirstPointZ());
1389 fprintf(dumpfile, "Last x: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetLastPointX(), listprintpointer->fMatchingtrack->fTrack.GetLastPointX());
1390 fprintf(dumpfile, "Last y: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetLastPointY(), listprintpointer->fMatchingtrack->fTrack.GetLastPointY());
1391 fprintf(dumpfile, "Last z: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetLastPointZ(), listprintpointer->fMatchingtrack->fTrack.GetLastPointZ());
1392 fprintf(dumpfile, " Pt: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetPt(), listprintpointer->fMatchingtrack->fTrack.GetPt());
1393 fprintf(dumpfile, " Psi: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetPsi(), listprintpointer->fMatchingtrack->fTrack.GetPsi());
1394 fprintf(dumpfile, " Tgl: %9.6f \t %9.6f \n", listprintpointer->fTrack.GetTgl(), listprintpointer->fMatchingtrack->fTrack.GetTgl());
1395 fprintf(dumpfile, "--------------\n");
1400 listprintpointer = listprintpointer->fNext;
1404 // print out not matched tracks from first track array:
1405 listprintpointer = fFirstTrackList;
1406 Int_t notmatchedtracknumber = 1;
1408 while(listprintpointer != NULL)
1410 if(listprintpointer->fMatchingindicator == 0)
1413 HLTInfo("Original Track, not matched with secondary track %d:", notmatchedtracknumber);
1414 HLTInfo("Clusters: %d",listprintpointer->fTrack.GetNHits());
1415 //PYTHIA INFORMATION ABOUT PARTICLE ID
1416 HLTInfo("First x: %9.6f \t first y: %9.6f \t first z: %9.6f",listprintpointer->fTrack.GetFirstPointX(),listprintpointer->fTrack.GetFirstPointY(),listprintpointer->fTrack.GetFirstPointZ());
1417 HLTInfo(" Last x: %9.6f \t last y: %9.6f \t last z: %9.6f", listprintpointer->fTrack.GetLastPointX(),listprintpointer->fTrack.GetLastPointY(),listprintpointer->fTrack.GetLastPointZ());
1418 HLTInfo(" Pt: %9.6f \t Psi: %9.6f \t Tgl: %9.6f", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl());
1419 HLTInfo("--------------");
1422 // print these results to file
1423 if(!fDumpFileName.IsNull())
1425 fprintf(dumpfile, "Original Track, not matched with secondary track %d: \n", notmatchedtracknumber);
1426 fprintf(dumpfile, "Clusters: %d \n",listprintpointer->fTrack.GetNHits());
1427 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());
1428 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());
1429 fprintf(dumpfile, " Pt: %9.6f \t Psi: %9.6f \t Tgl: %9.6f \n", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl());
1430 fprintf(dumpfile, "--------------\n");
1433 ++notmatchedtracknumber;
1436 listprintpointer = listprintpointer->fNext;
1439 // print out not matched tracks from second track array:
1440 listprintpointer = fSecondTrackList;
1441 notmatchedtracknumber = 1;
1443 while(listprintpointer != NULL)
1445 if(listprintpointer->fMatchingindicator == 0)
1448 HLTInfo("Secondary Track, not matched with original track %d:", notmatchedtracknumber);
1449 HLTInfo("Clusters: %d",listprintpointer->fTrack.GetNHits());
1450 //PYTHIA INFORMATION ABOUT PARTICLE ID
1451 HLTInfo("First x: %9.6f \t first y: %9.6f \t first z: %9.6f",listprintpointer->fTrack.GetFirstPointX(),listprintpointer->fTrack.GetFirstPointY(),listprintpointer->fTrack.GetFirstPointZ());
1452 HLTInfo(" Last x: %9.6f \t last y: %9.6f \t last z: %9.6f", listprintpointer->fTrack.GetLastPointX(),listprintpointer->fTrack.GetLastPointY(),listprintpointer->fTrack.GetLastPointZ());
1453 HLTInfo(" Pt: %9.6f \t Psi: %9.6f \t Tgl: %9.6f", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl());
1454 HLTInfo("--------------");
1456 // print these results to file
1457 if(!fDumpFileName.IsNull())
1459 fprintf(dumpfile, "Secondary Track, not matched with original track %d: \n", notmatchedtracknumber);
1460 fprintf(dumpfile, "Clusters: %d \n",listprintpointer->fTrack.GetNHits());
1461 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());
1462 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());
1463 fprintf(dumpfile, " Pt: %9.6f \t Psi: %9.6f \t Tgl: %9.6f \n", listprintpointer->fTrack.GetPt(), listprintpointer->fTrack.GetPsi(), listprintpointer->fTrack.GetTgl());
1464 fprintf(dumpfile, "--------------\n");
1467 ++notmatchedtracknumber;
1470 listprintpointer = listprintpointer->fNext;
1473 // close output analysis file
1474 if(!fDumpFileName.IsNull())
1479 // if results should be written to graphical output:
1480 if(!fGraphFileName.IsNull())
1482 CreateGraphs(1); // specifiy if absolute or rel. differences should be saved (CreateGraphs(0) or CreateGraphs()/ CreateGraphs(1)
1485 // free reserved space