]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/benchmark/AliHLTBenchExternalTrackComponent.cxx
removing the HLT autoconf build system, however keep on using that for the
[u/mrichter/AliRoot.git] / HLT / benchmark / AliHLTBenchExternalTrackComponent.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: Matthias Richter <Matthias.Richter@ift.uib.no>        *
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   AliHLTBenchExternalTrackComponent.cxx
20     @author Matthias Richter
21     @date   2008-10-30
22     @brief  Benchmark component for AliExternalTrackParam transportation.
23 */
24
25 // see header file for class documentation
26 // or
27 // refer to README to build package
28 // or
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
30
31 #include <cassert>
32 #include "AliHLTBenchExternalTrackComponent.h"
33 #include "AliExternalTrackParam.h"
34 #include "AliHLTExternalTrackParam.h"
35 #include "TString.h"
36 #include "TDatime.h"
37 #include "TRandom.h"
38 #include "TMath.h"
39 #include "TObjArray.h"
40 #include "TClonesArray.h"
41 #include "TClass.h"
42 #include "TList.h"
43
44 /** global object for component registration */
45 AliHLTBenchExternalTrackComponent gAliHLTBenchExternalTrackComponent;
46
47 /** ROOT macro for the implementation of ROOT specific class methods */
48 ClassImp(AliHLTBenchExternalTrackComponent)
49
50 AliHLTBenchExternalTrackComponent::AliHLTBenchExternalTrackComponent()
51   : AliHLTProcessor()
52   , fVerbosity(0)
53   , fDisableRegistry(false)
54   , fMode(kPublishingOff)
55   , fMaxSize(10000)
56   , fMinSize(100)
57   , fEventModulo(-1)
58   , fRangeOffset(0)
59   , fRangeMultiplicator(1.0)
60   , fpTcArray(NULL)
61   , fpTObjArray(NULL)
62   , fpDice(NULL)
63 {
64   // see header file for class documentation
65   // or
66   // refer to README to build package
67   // or
68   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
69 }
70
71 AliHLTBenchExternalTrackComponent::~AliHLTBenchExternalTrackComponent()
72 {
73   // see header file for class documentation
74 }
75
76 const char* AliHLTBenchExternalTrackComponent::GetComponentID()
77 {
78   // see header file for class documentation
79   return "BenchmarkAliExternalTrackParam";
80 }
81
82 void AliHLTBenchExternalTrackComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
83 {
84   // see header file for class documentation
85   list.clear();
86   list.push_back(kAliHLTAnyDataType);
87 }
88
89 AliHLTComponentDataType AliHLTBenchExternalTrackComponent::GetOutputDataType()
90 {
91   // see header file for class documentation
92   return kAliHLTMultipleDataType;
93 }
94
95 int AliHLTBenchExternalTrackComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
96 {
97   // see header file for class documentation
98   tgtList.clear();
99   tgtList.push_back(kAliHLTDataTypeTObject);
100   return tgtList.size();
101 }
102
103 void AliHLTBenchExternalTrackComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
104 {
105   // see header file for class documentation
106   constBase=sizeof(AliHLTExternalTrackParam)*fMaxSize;
107   inputMultiplier=0.0;
108 }
109
110 AliHLTComponent* AliHLTBenchExternalTrackComponent::Spawn()
111 {
112   // see header file for class documentation
113   return new AliHLTBenchExternalTrackComponent;
114 }
115
116 int AliHLTBenchExternalTrackComponent::DoInit(int argc, const char** argv)
117 {
118   // see header file for class documentation
119   int iResult=0;
120   TString argument="";
121   bool bMissingParam=0;
122   char* cpErr=NULL;
123   int i=0;
124   for (; i<argc && iResult>=0; i++) {
125     cpErr=NULL;
126     argument=argv[i];
127     if (argument.IsNull()) continue;
128
129     // -tobjarray
130     if (argument.CompareTo("-tobjarray")==0) {
131       fMode=ktobjarray;
132
133     // -tclonesarray
134     } else if (argument.CompareTo("-tclonesarray")==0) {
135       fMode=ktclonesarray;
136
137     // -carray
138     } else if (argument.CompareTo("-carray")==0) {
139       fMode=kcarray;
140
141     // -maxsize
142     } else if (argument.CompareTo("-maxsize")==0) {
143       if ((bMissingParam=(++i>=argc))) break;
144       fMaxSize=strtoul( argv[i], &cpErr ,0);
145       if ( *cpErr ) break;
146
147     // -minsize
148     } else if (argument.CompareTo("-minsize")==0) {
149       if ((bMissingParam=(++i>=argc))) break;
150       fMinSize=strtoul( argv[i], &cpErr ,0);
151       if ( *cpErr ) break;
152
153     // -rangemodulo
154     } else if (argument.CompareTo("-rangemodulo")==0) {
155       if ((bMissingParam=(++i>=argc))) break;
156       fEventModulo=strtoul( argv[i], &cpErr ,0);
157       if ( *cpErr ) break;
158
159     // -rangeoffset
160     } else if (argument.CompareTo("-rangeoffset")==0) {
161       if ((bMissingParam=(++i>=argc))) break;
162       fRangeOffset=strtoul( argv[i], &cpErr ,0);
163       if ( *cpErr ) break;
164
165     // -rangefactor
166     } else if (argument.CompareTo("-rangefactor")==0) {
167       if ((bMissingParam=(++i>=argc))) break;
168       fRangeMultiplicator=strtof( argv[i], &cpErr);
169       if ( *cpErr ) break;
170
171     // -verbosity
172     } else if (argument.CompareTo("-verbosity")==0) {
173       if ((bMissingParam=(++i>=argc))) break;
174       fVerbosity=strtoul( argv[i], &cpErr ,0);
175       if ( *cpErr ) break;
176
177     // -nocheck
178     } else if (argument.CompareTo("-nocheck")==0) {
179       fDisableRegistry=true;
180
181     } else {
182       HLTError("unknown argument %s", argument.Data());
183       iResult=-EINVAL;
184     }
185   }
186
187   if (cpErr && *cpErr) {
188     HLTError("Cannot convert specifier '%s' for argument '%s'", argv[i], argument.Data());
189     iResult=-EINVAL;
190   } else if (bMissingParam) {
191     HLTError("missing parameter for argument %s", argument.Data());
192     iResult=-EINVAL;
193   }
194
195   // use a TClonesArray internally and for publishing in case of
196   // ktclonesarray. The maximum size and all members are allocated
197   fpTcArray=new TClonesArray("AliExternalTrackParam", fMaxSize);
198   if (fpTcArray) {
199     fpTcArray->ExpandCreate(fMaxSize);
200     switch (fMode) {
201     case kPublishingOff:
202       // nothing to do
203       break;
204     case ktclonesarray:
205     case kcarray:
206       // nothing to do
207       break;
208     case ktobjarray:
209       // use a TObjArray which owns the objects
210       fpTObjArray=new TObjArray;
211       if (!fpTObjArray) {
212         iResult=-ENOMEM;
213       }
214       break;
215     default:
216       HLTError("unknown publishing mode %d", fMode);
217       iResult=-EINVAL;
218     }
219   } else {
220     iResult=-ENOMEM;
221   }
222
223   if (iResult>=0) {
224     fpDice=new TRandom;
225     if (fpDice) {
226       TDatime dt;
227       fpDice->SetSeed(dt.Get());
228     } else {
229        iResult=-ENOMEM;
230     }
231   }
232
233   return iResult;
234 }
235
236 int AliHLTBenchExternalTrackComponent::DoDeinit()
237 {
238   // see header file for class documentation
239   if (fpTObjArray) delete fpTObjArray;
240   fpTObjArray=NULL;
241
242   if (fpTcArray) {
243     fpTcArray->Delete();
244     delete fpTcArray;
245     fpTcArray=NULL;
246   }
247
248   if (fpDice) delete fpDice;
249   fpDice=NULL;
250
251   return 0;
252 }
253
254 int AliHLTBenchExternalTrackComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
255                                                         const AliHLTComponentBlockData* /*blocks*/,
256                                                         AliHLTComponentTriggerData& /*trigData*/,
257                                                         AliHLTUInt8_t* outputPtr, 
258                                                         AliHLTUInt32_t& size,
259                                                         AliHLTComponentBlockDataList& outputBlocks)
260 {
261   // see header file for class documentation
262
263   // scan the input data blocks for types kAliHLTDataTypeTrack and
264   // kAliHLTDataTypeTObjArray with AliExternalTrackParam objects
265   // 
266   // If data is available, the array is extracted into the internal
267   // TClonesArray fpTcArray. If no input data is available, the array
268   // is filled randomly, both in size and content.
269   //
270   // The array is published according to the fMode member either as
271   // TClonesArray, TObjArray or C-structure (AliHLTExternalTrackParam)
272
273   int iResult=0;
274   AliHLTUInt32_t capacity=size;
275   size=0;
276   if (!IsDataEvent()) return 0;
277   if (!fpTcArray || !fpDice) return -ENODEV;
278
279   unsigned int arraySize=0;
280   const AliHLTComponentBlockData* pBlock=NULL;
281   const TObject* pObject=NULL;
282   const TObjArray* pObjArray=NULL;
283   AliHLTUInt32_t spec=0;
284   if ((pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack))!=NULL) {
285     if (pBlock->fSize%sizeof(AliHLTExternalTrackParam)==0) {
286       arraySize=pBlock->fSize/sizeof(AliHLTExternalTrackParam);
287       if (fpTcArray) {
288         spec=GetSpecification(pBlock);
289         fpTcArray->ExpandCreate(arraySize);
290         if ((iResult=ReadFromStruct(fpTcArray, reinterpret_cast<AliHLTExternalTrackParam*>(pBlock->fPtr), arraySize))<0) {
291           HLTError("can not convert struct to TObjArray, error %d", iResult);
292         } else {
293           pObjArray=fpTcArray;
294           HLTDebug("created TObjArray from C struct, specification 0x%08x", spec);
295         }
296       } else {
297         iResult=-ENOMEM;
298       }
299     } else {
300       HLTWarning("data block does not match the AliHLTExternalTrackParam struct size");
301     }
302   } else if ((pObject=GetFirstInputObject(kAliHLTDataTypeTObjArray))!=NULL) {
303     pObjArray=dynamic_cast<const TObjArray*>(pObject);
304     if (pObjArray) {
305       spec=GetSpecification(pObject);
306       HLTDebug("extracted TObjArray, specification 0x%08x", spec);
307     } else {
308       HLTWarning("object of data type %s is not of type TObjArray, ignoring it",
309                  DataType2Text(GetDataType()).c_str());
310     }
311   }
312   
313   if (pObjArray) {
314     // we already have the object array, if publishing mode is >0 the
315     // block is forwarded according to the mode
316     // 
317     // if specification is >0 the source object is searched in the list
318     // and compared with the extracted one.
319     if (spec>0) {
320       TObject* pSrcObject=FindObject(spec);
321       TObjArray* pSrcArray=NULL;
322       if (!pSrcObject) {
323         HLTError("can not find source object of id 0x%08x", spec);
324       } else if ((pSrcArray=dynamic_cast<TObjArray*>(pSrcObject))==NULL) {
325         HLTError("type cast failed for object");
326       } else if (!Compare(pObjArray, pSrcArray)) {
327         HLTError("extracted array differs from source array");
328         if (fVerbosity>1) {
329           for (int i=0; i<pObjArray->GetEntries(); i++) {
330             pObjArray->At(i)->Print();
331           }
332           HLTError("dump of original array");
333           for (int i=0; i<pSrcArray->GetEntries(); i++) {
334             pSrcArray->At(i)->Print();
335           }
336         }
337       } else {
338         if (fVerbosity>0) {
339           HLTInfo("extracted array of size %d", pObjArray->GetEntries());
340         }
341       }
342     }
343   } else {
344     // no array to forward, create a new one
345     arraySize=fpDice->Integer(fMaxSize-fMinSize)+fMinSize;
346     fpTcArray->ExpandCreate(arraySize);
347     for (unsigned int track=0; track<arraySize; track++) {
348       FillRandom(reinterpret_cast<AliExternalTrackParam*>(fpTcArray->At(track)));
349     }
350     pObjArray=fpTcArray;
351   }
352
353   switch (fMode) {
354   case kPublishingOff:
355     break;
356   case ktclonesarray:
357       // register the object if it is not a forwarded array
358     if (spec==0 && pObjArray==fpTcArray) spec=Register(fpTcArray);
359     if ((iResult=PushBack(const_cast<TObjArray*>(pObjArray), kAliHLTDataTypeTObjArray, spec))>=0) {
360       if (fVerbosity>0) {
361         HLTInfo("publishing TClonesarray, %d elements, specification 0x%08x", pObjArray->GetEntries(), spec);
362       }
363     }
364     break;
365   case kcarray:
366     iResult=SerializeToStruct(pObjArray, outputPtr, capacity);
367     if (iResult>=0) {
368       size=iResult;
369       AliHLTComponentBlockData bd;
370       FillBlockData(bd);
371       bd.fPtr=NULL;
372       bd.fOffset=0;
373       bd.fSize=size;
374       bd.fDataType=kAliHLTDataTypeTrack;
375       // register the object if it is not a forwarded array
376       if (spec==0 && pObjArray==fpTcArray) spec=Register(fpTcArray);
377       bd.fSpecification=spec;
378       outputBlocks.push_back(bd);
379       if (fVerbosity>0) {
380         HLTInfo("publishing C array, %d elements, specification 0x%08x", bd.fSize/sizeof(AliHLTExternalTrackParam), bd.fSpecification);
381       }
382     }
383     break;
384   case ktobjarray:
385     // use a TObjArray which owns the objects
386     if (fpTObjArray) {
387       fpTObjArray->Clear();
388       int entries=pObjArray->GetEntries();
389       for (int i=0; i<entries; i++) {
390         fpTObjArray->Add(pObjArray->At(i));
391       }
392
393       // register the object if it is not a forwarded array
394       if (spec==0 && pObjArray==fpTcArray) spec=Register(fpTObjArray);
395       if ((iResult=PushBack(fpTObjArray, kAliHLTDataTypeTObjArray, spec))>=0) {
396         if (fVerbosity>0) {
397           HLTInfo("publishing TObjArray, %d elements, specification 0x%08x", pObjArray->GetEntries(), spec);
398         }
399       }
400     } else {
401       HLTError("object array not initialized");
402       iResult=-EFAULT;
403     }
404     break;
405   default:
406     HLTError("unknown publishing mode %d", fMode);
407     iResult=-EINVAL;
408   }
409
410   if (fEventModulo>0 &&
411       ((GetEventCount()+1)%fEventModulo)==0 &&
412       (fRangeMultiplicator!=1.0 || fRangeOffset!=0)) {
413     fMaxSize+=fRangeOffset;
414     fMaxSize=(int)(fMaxSize*fRangeMultiplicator);
415     fMinSize+=fRangeOffset;
416     fMinSize=(int)(fMinSize*fRangeMultiplicator);
417     if (fMaxSize<1) fMaxSize=1;
418     if (fMinSize<1) fMinSize=1;
419   }
420   return iResult;
421 }
422
423 int AliHLTBenchExternalTrackComponent::FillRandom(AliExternalTrackParam* track, int fillCov)
424 {
425   // see header file for class documentation
426   int iResult=0;
427   if (!track) return -EINVAL;
428   Double_t param[5];
429   Double_t covar[15];
430
431   static TRandom* rand=NULL;
432   if (!rand) {
433     rand=new TRandom;
434     if (!rand) return -ENOMEM;
435     TDatime dt;
436     rand->SetSeed(dt.Get());
437   }
438
439   param[0]=(Double_t)rand->Integer(100);
440   param[1]=(Double_t)rand->Integer(100);
441   param[2]=(Double_t)rand->Integer(100)/100;
442   param[3]=(Double_t)rand->Integer(100)/100;
443   param[4]=(Double_t)rand->Integer(100);
444   for (int i=0; i<15 && i<fillCov; i++) {
445     covar[i]=(Double_t)rand->Integer(1000)/1000;
446   }
447   track->Set((Double_t)rand->Integer(1000),(Double_t)rand->Integer(100)/100, param, covar );
448   return iResult;
449 }
450
451 int AliHLTBenchExternalTrackComponent::SerializeToStruct(const TObjArray* pArray, AliHLTUInt8_t* buffer, unsigned int size)
452 {
453   // see header file for class documentation
454   int iResult=0;
455   if (!pArray || !buffer) return -EINVAL;
456   int entries=pArray->GetEntries();
457   if (entries*sizeof(AliHLTExternalTrackParam)>size) return -ENOSPC;
458
459   AliHLTLogging log;
460   AliHLTExternalTrackParam* tgtParam=reinterpret_cast<AliHLTExternalTrackParam*>(buffer);
461   for (int track=0; track<entries; track++, tgtParam++) {
462     TObject* pObj=pArray->At(track);
463     if (!pObj) {
464       log.LoggingVarargs(kHLTLogError, "AliHLTBenchExternalTrackComponent", "SerializeToStruct" , __FILE__ , __LINE__ ,
465                           "internal mismatch: can not get object %d from array. aborting", track);
466       iResult=-ENOENT;
467       break;
468     }
469     AliExternalTrackParam* srcParam=dynamic_cast<AliExternalTrackParam*>(pObj);
470     if (!srcParam) {
471       log.LoggingVarargs(kHLTLogError, "AliHLTBenchExternalTrackComponent", "SerializeToStruct" , __FILE__ , __LINE__ ,
472                           "object %d has wrong type %s, expecting %s", track, pObj->Class()->GetName(), "AliExternalTrackParam");
473       iResult=-ENOENT;
474       break;
475     }
476     
477     //= srcParam->GetAlpha();
478     tgtParam->fX = srcParam->GetX();
479     tgtParam->fY = srcParam->GetY();
480     tgtParam->fZ = srcParam->GetZ();
481     tgtParam->fLastX=0;
482     tgtParam->fLastY=0;
483     tgtParam->fLastZ=0;
484     tgtParam->fSinPsi = srcParam->GetSnp();
485     tgtParam->fTgl = srcParam->GetTgl();
486     tgtParam->fq1Pt = srcParam->GetSigned1Pt();
487     const Double_t *cov=srcParam->GetCovariance();
488
489     for (int i=0; i<15; i++) {
490       tgtParam->fC[i]=cov[i];
491     }
492     tgtParam->fNPoints=0;
493   }
494
495   if (iResult>=0) {
496     iResult=entries*sizeof(AliHLTExternalTrackParam);
497   }
498
499   return iResult;
500 }
501
502 int AliHLTBenchExternalTrackComponent::ReadFromStruct(TObjArray* pTgtArray, AliHLTExternalTrackParam* pArray, unsigned int arraySize)
503 {
504   // see header file for class documentation
505   AliHLTLogging log;
506   if ((unsigned)pTgtArray->GetEntries()<arraySize) {
507     log.LoggingVarargs(kHLTLogError, "AliHLTBenchExternalTrackComponent", "ReadFromStruct" , __FILE__ , __LINE__ ,
508                        "not enough space in target array: %d, required %d", pTgtArray->GetEntries(), arraySize);
509     return -ENOSPC;
510   }
511
512   Double_t param[5];
513   Double_t covar[15];
514
515   for (unsigned int track=0; track<arraySize; track++) {
516     AliExternalTrackParam* tp=dynamic_cast<AliExternalTrackParam*>(pTgtArray->At(track));
517     if (tp) {
518       // enable if AliExternalTrackParam allows templates
519       //tp->Set(pArray[track].fX, (Float_t)0.0, &pArray[track].fY, pArray[track].fC);
520       param[0]=pArray[track].fY;
521       param[1]=pArray[track].fZ;
522       param[2]=pArray[track].fSinPsi;
523       param[3]=pArray[track].fTgl;
524       param[4]=pArray[track].fq1Pt;
525       for (int i=0; i<15; i++) {
526         covar[i]=pArray[track].fC[i];
527       }
528       tp->Set((Double_t)pArray[track].fX, (Double_t)0.0, param, covar);
529     } else {
530       log.LoggingVarargs(kHLTLogError, "AliHLTBenchExternalTrackComponent", "ReadFromStruct" , __FILE__ , __LINE__ ,
531                          "invalid object type %s", pTgtArray->At(track)->Class()->GetName());
532       return -EFAULT;
533     }
534   }
535   return arraySize;
536 }
537
538 AliHLTUInt32_t AliHLTBenchExternalTrackComponent::CalcChecksum(const TObjArray* pArray)
539 {
540   // see header file for class documentation
541   AliHLTUInt32_t crc=0;
542   int entries=0;
543   if (pArray && (entries=pArray->GetEntries())>0) {
544     AliHLTLogging log;
545     unsigned int bufferSize=entries*sizeof(AliHLTExternalTrackParam);
546     AliHLTUInt8_t* buffer=new AliHLTUInt8_t[bufferSize];
547     if (buffer && SerializeToStruct(pArray, buffer, bufferSize)>0) {
548       crc=CalculateChecksum(buffer, bufferSize);
549     } else if (buffer) {
550       log.LoggingVarargs(kHLTLogError, "AliHLTBenchExternalTrackComponent", "SerializeToStruct" , __FILE__ , __LINE__ ,
551                          "failed to serialize TObjArray");
552     }
553   }
554   return crc;
555 }
556
557 bool AliHLTBenchExternalTrackComponent::Compare(const TObjArray* array1, const TObjArray* array2) 
558 {
559   // see header file for class documentation
560   if (!array1 || !array2) return false;
561   int entries=array1->GetEntries();
562   if (entries!=array2->GetEntries()) {
563     return false;
564   }
565
566   for (int i=0; i<entries; i++) {
567     TObject* object1=array1->At(i);
568     TObject* object2=array2->At(i);
569     if (!object1 || !object2) return false;
570
571     AliExternalTrackParam* param1=dynamic_cast<AliExternalTrackParam*>(object1);
572     AliExternalTrackParam* param2=dynamic_cast<AliExternalTrackParam*>(object2);
573     if (!param1 || !param2) return false;
574
575     if (TMath::Abs(param1->GetX()-param2->GetX())>0.0001) return false;
576     if (TMath::Abs(param1->GetY()-param2->GetY())>0.0001) return false;
577     if (TMath::Abs(param1->GetZ()-param2->GetZ())>0.0001) return false;
578     if (TMath::Abs(param1->GetSnp()-param2->GetSnp())>0.0001) return false;
579     if (TMath::Abs(param1->GetTgl()-param2->GetTgl())>0.0001) return false;
580     if (TMath::Abs(param1->GetSigned1Pt()-param2->GetSigned1Pt())>0.0001) return false;
581
582     if (TMath::Abs(param1->GetSigmaY2()-param2->GetSigmaY2())>0.0001) return false;
583     if (TMath::Abs(param1->GetSigmaZY()-param2->GetSigmaZY())>0.0001) return false;
584     if (TMath::Abs(param1->GetSigmaZ2()-param2->GetSigmaZ2())>0.0001) return false;
585     if (TMath::Abs(param1->GetSigmaSnpY()-param2->GetSigmaSnpY())>0.0001) return false;
586     if (TMath::Abs(param1->GetSigmaSnpZ()-param2->GetSigmaSnpZ())>0.0001) return false;
587     if (TMath::Abs(param1->GetSigmaSnp2()-param2->GetSigmaSnp2())>0.0001) return false;
588     if (TMath::Abs(param1->GetSigmaTglY()-param2->GetSigmaTglY())>0.0001) return false;
589     if (TMath::Abs(param1->GetSigmaTglZ()-param2->GetSigmaTglZ())>0.0001) return false;
590     if (TMath::Abs(param1->GetSigmaTglSnp()-param2->GetSigmaTglSnp())>0.0001) return false;
591     if (TMath::Abs(param1->GetSigmaTgl2()-param2->GetSigmaTgl2())>0.0001) return false;
592     if (TMath::Abs(param1->GetSigma1PtY()-param2->GetSigma1PtY())>0.0001) return false;
593     if (TMath::Abs(param1->GetSigma1PtZ()-param2->GetSigma1PtZ())>0.0001) return false;
594     if (TMath::Abs(param1->GetSigma1PtSnp()-param2->GetSigma1PtSnp())>0.0001) return false;
595     if (TMath::Abs(param1->GetSigma1PtTgl()-param2->GetSigma1PtTgl())>0.0001) return false;
596     if (TMath::Abs(param1->GetSigma1Pt2()-param2->GetSigma1Pt2())>0.0001) return false;
597   }
598
599   return true;
600 }
601
602 TList* AliHLTBenchExternalTrackComponent::fgpRegistry=NULL;
603
604 TObject* AliHLTBenchExternalTrackComponent::FindObject(AliHLTUInt32_t id)
605 {
606   // see header file for class documentation
607   if (!fgpRegistry) return NULL;
608
609   TIter iter(fgpRegistry);
610   TObject* pObj=NULL;
611   while ((pObj=iter.Next())) {
612     if (pObj->GetUniqueID()==id) break;
613   }
614   return pObj;
615 }
616
617 AliHLTUInt32_t AliHLTBenchExternalTrackComponent::Register(TObject* pObject)
618 {
619   // see header file for class documentation
620   if (fDisableRegistry) return 0;
621
622   if (!fgpRegistry) {
623     fgpRegistry=new TList;
624   }
625   if (!fgpRegistry) return 0;
626
627   TObject* pExist=fgpRegistry->FindObject(pObject);
628   if (pExist) return pExist->GetUniqueID();
629
630   AliHLTUInt32_t id=AliHLTComponent::CalculateChecksum(reinterpret_cast<AliHLTUInt8_t*>(&pObject), sizeof(TObject*));
631   pExist=FindObject(id);
632   assert(pExist==NULL);
633   pObject->SetUniqueID(id);
634   fgpRegistry->Add(pObject);
635   HLTInfo("adding object %p with id 0x%08x to registry", pObject, id);
636   return id;
637 }
638
639 int AliHLTBenchExternalTrackComponent::Unregister(TObject* pObject)
640 {
641   // see header file for class documentation
642   if (!fgpRegistry) return 0;
643
644   if (fgpRegistry->Remove(pObject)==pObject) return 0;
645
646   return -ENOENT;
647 }