]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/HWCFemulator/AliHLTTPCHWCFConsistencyControlComponent.cxx
removing unused classes from repository
[u/mrichter/AliRoot.git] / HLT / TPCLib / HWCFemulator / AliHLTTPCHWCFConsistencyControlComponent.cxx
1 //****************************************************************************
2 //* This file is property of and copyright by the ALICE HLT Project          * 
3 //* ALICE Experiment at CERN, All rights reserved.                           *
4 //*                                                                          *
5 //* Primary Authors: Sergey Gorbunov, Torsten Alt                            *
6 //* Developers:      Sergey Gorbunov <sergey.gorbunov@fias.uni-frankfurt.de> *
7 //*                  Torsten Alt <talt@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   AliHLTTPCHWCFConsistencyControlComponent.cxx
20 //  @author Sergey Gorbunov <sergey.gorbunov@fias.uni-frankfurt.de>
21 //  @author Torsten Alt <talt@cern.ch> 
22 //  @brief  Comparison of TPC clusters produced by FPGA clusterfinder and by FPGA Emulator 
23 //  @brief  ( see AliHLTTPCHWCFEmulator class )
24 //  @note
25
26
27 #if __GNUC__>= 3
28 using namespace std;
29 #endif
30 #include "AliHLTTPCHWCFConsistencyControlComponent.h"
31
32 #include "AliHLTTPCDefinitions.h"
33 #include "AliHLTTPCHWCFDataTypes.h"
34 #include "AliHLTTPCClusterMCData.h"
35
36 #include "AliGRPObject.h"
37 #include "AliCDBEntry.h"
38 #include "AliCDBManager.h"
39 #include "AliRawDataHeader.h"
40 #include <cstdlib>
41 #include <cerrno>
42 #include "TString.h"
43 #include "TObjString.h"
44 #include "TObjArray.h"
45 #include "TH1F.h"
46
47 #include <sys/time.h>
48 #include "TFile.h"
49
50 AliHLTTPCHWCFConsistencyControlComponent::AliHLTTPCHWCFConsistencyControlComponent()
51   :
52   AliHLTProcessor(),
53   fNDismatch(0),
54   fNBlocks(0),
55   fBenchmark("TPCHWConsistencyControl"), 
56   fHistHeaderAll(0),
57   fHistHeaderGood(0),
58   fHistClusterAll(0),
59   fHistClusterGood(0),
60   fProfHeader(0),
61   fProfCluster(0)
62 {
63   // see header file for class documentation
64   // or
65   // refer to README to build package
66   // or
67   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
68 }
69
70
71 AliHLTTPCHWCFConsistencyControlComponent::AliHLTTPCHWCFConsistencyControlComponent(const AliHLTTPCHWCFConsistencyControlComponent&)
72   :
73   AliHLTProcessor(),
74   fNDismatch(0),
75   fNBlocks(0),
76   fBenchmark("TPCHWConsistencyControl"),
77   fHistHeaderAll(0),
78   fHistHeaderGood(0),
79   fHistClusterAll(0),
80   fHistClusterGood(0),
81   fProfHeader(0),
82   fProfCluster(0)
83 {
84   // dummy
85 }
86
87 AliHLTTPCHWCFConsistencyControlComponent& AliHLTTPCHWCFConsistencyControlComponent::operator=(const AliHLTTPCHWCFConsistencyControlComponent&)
88 {
89   // dummy
90   return *this;
91 }
92
93 AliHLTTPCHWCFConsistencyControlComponent::~AliHLTTPCHWCFConsistencyControlComponent()
94 {
95   // see header file for class documentation
96   DoDeinit();
97 }
98
99 // Public functions to implement AliHLTComponent's interface.
100 // These functions are required for the registration process
101
102 const char* AliHLTTPCHWCFConsistencyControlComponent::GetComponentID()
103 {
104   // see header file for class documentation
105   return "TPCHWCFConsistenyControl";
106 }
107
108 void AliHLTTPCHWCFConsistencyControlComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
109 {
110   // see header file for class documentation
111   list.clear();
112   list.push_back( AliHLTTPCDefinitions::fgkHWClustersDataType| kAliHLTDataOriginTPC );
113 }
114
115 AliHLTComponentDataType AliHLTTPCHWCFConsistencyControlComponent::GetOutputDataType()
116 {
117   // see header file for class documentation
118   return kAliHLTDataTypeHistogram  | kAliHLTDataOriginTPC;
119 }
120
121
122 void AliHLTTPCHWCFConsistencyControlComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
123 {
124   // see header file for class documentation
125   // XXX TODO: Find more realistic values.  
126   constBase = 1000;
127   inputMultiplier = 0;
128 }
129
130
131 AliHLTComponent* AliHLTTPCHWCFConsistencyControlComponent::Spawn()
132 {
133   // see header file for class documentation
134   return new AliHLTTPCHWCFConsistencyControlComponent();
135 }
136
137 void AliHLTTPCHWCFConsistencyControlComponent::GetOCDBObjectDescription( TMap* const /*targetMap*/){
138 // Get a list of OCDB object description needed for the particular component
139   /*
140   if (!targetMap) return;
141   
142   // OCDB entries for component arguments
143   targetMap->Add(new TObjString("HLT/ConfigTPC/TPCHWClusterFinder"), new TObjString("component arguments, empty at the moment"));
144   */
145 }
146
147
148 int AliHLTTPCHWCFConsistencyControlComponent::DoInit( int argc, const char** argv )
149 {
150   // see header file for class documentation
151
152   fHistHeaderAll = new TH1F("hHWCFHeaderAll", "fHistHeaderAll",6,0.,6.);
153   fHistHeaderGood = new TH1F("hHWCFHeaderGood", "fHistHeaderGood",6,0.,6.);
154   fHistClusterAll = new TH1F("hHWCFClusterAll", "fHistClusterAll",8,0.,8.);
155   fHistClusterGood = new TH1F("hHWCFClusterGood", "fHistClusterGood",8,0.,8.);
156   fProfHeader = new TH1F("pHWCFHeader", "HWCF: Consistency of header data", 6, 0., 6.);
157   fProfCluster = new TH1F("pHWCFClusters", "HWCF: Consisteny of cluster data", 8, 0., 8.);
158
159   TString arguments = "";
160   for ( int i = 0; i < argc; i++ ) {
161     if ( !arguments.IsNull() ) arguments += " ";
162     arguments += argv[i];
163   }
164
165   return Configure( NULL, NULL, arguments.Data()  );
166 }
167
168 int AliHLTTPCHWCFConsistencyControlComponent::Reconfigure( const char* cdbEntry, const char* chainId )
169 {
170   // Reconfigure the component from OCDB
171
172   return Configure( cdbEntry, chainId, NULL );
173 }
174
175 int AliHLTTPCHWCFConsistencyControlComponent::ScanConfigurationArgument(int argc, const char** argv)
176 {
177   // see header file for class documentation
178   TString arguments = "";
179   for ( int i = 0; i < argc; i++ ) {
180     if ( !arguments.IsNull() ) arguments += " ";
181     arguments += argv[i];
182   }
183   return ReadConfigurationString(arguments);
184 }
185
186
187
188 void AliHLTTPCHWCFConsistencyControlComponent::SetDefaultConfiguration()
189 {
190   // Set default configuration for the FPGA ClusterFinder Emulator component
191   // Some parameters can be later overwritten from the OCDB
192
193   fNDismatch = 0;
194   fNBlocks = 0;
195   fBenchmark.Reset();
196   fBenchmark.SetTimer(0,"total");
197   fBenchmark.SetTimer(1,"reco");    
198 }
199
200 int AliHLTTPCHWCFConsistencyControlComponent::ReadConfigurationString(  const char* arguments )
201 {
202   // Set configuration parameters for the FPGA ClusterFinder Emulator component
203   // from the string
204
205   int iResult = 0;
206   if ( !arguments ) return iResult;
207
208   TString allArgs = arguments;
209   TString argument;
210   int bMissingParam = 0;
211
212   TObjArray* pTokens = allArgs.Tokenize( " " );
213
214   int nArgs =  pTokens ? pTokens->GetEntries() : 0;
215
216   for ( int i = 0; i < nArgs; i++ ) {
217     argument = ( ( TObjString* )pTokens->At( i ) )->GetString();
218     if ( argument.IsNull() ) continue;
219
220     //if ( argument.CompareTo( "-deconvolute-time" ) == 0 ) {
221     //if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
222     //fDoDeconvTime  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
223     //HLTInfo( "Time deconvolution is set to: %d", fDoDeconvTime );
224     //continue;
225     //}
226
227     HLTError( "Unknown option \"%s\"", argument.Data() );
228     iResult = -EINVAL;
229   }
230   delete pTokens;
231
232   if ( bMissingParam ) {
233     HLTError( "Specifier missed for parameter \"%s\"", argument.Data() );
234     iResult = -EINVAL;
235   }
236
237   return iResult;
238 }
239
240
241 int AliHLTTPCHWCFConsistencyControlComponent::ReadCDBEntry( const char* /*cdbEntry*/, const char* /*chainId*/ )
242 {
243   // Read configuration from OCDB
244   return 0;
245   /*
246   const char* defaultNotify = "";
247
248   if ( !cdbEntry ){
249     cdbEntry = "HLT/ConfigTPC/TPCHWClusterFinder";
250     defaultNotify = " (default)";
251     chainId = 0;
252   }
253
254   HLTInfo( "configure from entry \"%s\"%s, chain id %s", cdbEntry, defaultNotify, ( chainId != NULL && chainId[0] != 0 ) ? chainId : "<none>" );
255   AliCDBEntry *pEntry = AliCDBManager::Instance()->Get( cdbEntry );//,GetRunNo());
256   
257   if ( !pEntry ) {
258     HLTError( "cannot fetch object \"%s\" from CDB", cdbEntry );
259     return -EINVAL;
260   }
261
262   TObjString* pString = dynamic_cast<TObjString*>( pEntry->GetObject() );
263
264   if ( !pString ) {
265     HLTError( "configuration object \"%s\" has wrong type, required TObjString", cdbEntry );
266     return -EINVAL;
267   }
268
269   HLTInfo( "received configuration object string: \"%s\"", pString->GetString().Data() );
270
271   return  ReadConfigurationString( pString->GetString().Data() );
272   */
273 }
274
275
276 int AliHLTTPCHWCFConsistencyControlComponent::Configure( const char* cdbEntry, const char* chainId, const char *commandLine )
277 {
278   // Configure the component
279   // There are few levels of configuration,
280   // parameters which are set on one step can be overwritten on the next step
281
282   //* read hard-coded values
283
284   SetDefaultConfiguration();
285
286   //* read the default CDB entry
287
288   int iResult1 = ReadCDBEntry( NULL, chainId );
289
290   //* read the actual CDB entry if required
291
292   int iResult2 = ( cdbEntry ) ? ReadCDBEntry( cdbEntry, chainId ) : 0;
293
294   //* read extra parameters from input (if they are)
295
296   int iResult3 = 0;
297
298   if ( commandLine && commandLine[0] != '\0' ) {
299     HLTInfo( "received configuration string from HLT framework: \"%s\"", commandLine );
300     iResult3 = ReadConfigurationString( commandLine );
301   }
302
303   return iResult1 ? iResult1 : ( iResult2 ? iResult2 : iResult3 );
304 }
305
306
307 int AliHLTTPCHWCFConsistencyControlComponent::DoDeinit()
308 {
309   // see header file for class documentation 
310   if( fHistHeaderAll ) delete fHistHeaderAll;
311   if( fHistHeaderGood ) delete fHistHeaderGood;
312   if( fHistClusterAll ) delete fHistClusterAll;
313   if( fHistClusterGood ) delete fHistClusterGood;
314   if( fProfHeader ) delete fProfHeader;
315   if( fProfCluster ) delete fProfCluster;
316
317   fHistHeaderAll = 0;
318   fHistHeaderGood = 0;
319   fHistClusterAll = 0;
320   fHistClusterGood = 0;
321   fProfHeader = 0;
322   fProfCluster = 0;
323   return 0;
324 }
325
326
327 int AliHLTTPCHWCFConsistencyControlComponent::DoEvent( const AliHLTComponentEventData& evtData, 
328                                                        const AliHLTComponentBlockData* blocks, 
329                                                        AliHLTComponentTriggerData& /*trigData*/, AliHLTUInt8_t* /*outputPtr*/, 
330                                                        AliHLTUInt32_t& size, 
331                                                        vector<AliHLTComponentBlockData>& /*outputBlocks*/ )
332 {
333   // see header file for class documentation
334
335   int iResult=0;
336   //AliHLTUInt32_t maxSize = size;
337   size = 0;
338   
339   if(!IsDataEvent()){
340     return 0;
341   }
342
343   fBenchmark.StartNewEvent();
344   fBenchmark.Start(0);
345   bool *checkedBlocks = new bool[evtData.fBlockCnt];
346   for( unsigned long i=0; i<evtData.fBlockCnt; i++ ) checkedBlocks[i] = 0;
347   
348   for ( unsigned long ndx1 = 0; ndx1 < evtData.fBlockCnt; ndx1++ ){
349     const AliHLTComponentBlockData* iter1 = blocks+ndx1;
350     fBenchmark.AddInput(iter1->fSize);
351     
352     if (  iter1->fDataType != (AliHLTTPCDefinitions::fgkHWClustersDataType | kAliHLTDataOriginTPC) 
353           ) continue;
354  
355     if( checkedBlocks[ndx1] ) continue;// block already checked
356
357     checkedBlocks[ndx1] = 1;
358
359     int slice1 = AliHLTTPCDefinitions::GetMinSliceNr( *iter1 );
360     int patch1 = AliHLTTPCDefinitions::GetMinPatchNr( *iter1 );
361    
362     bool header1OK = 1;
363     if( iter1->fSize>0 ){
364       header1OK = ( iter1->fSize % sizeof(AliHLTUInt32_t) == 0 ) 
365         && ( iter1->fSize >= sizeof(AliRawDataHeader) ) 
366         && (iter1->fPtr != NULL );      
367     }
368     
369     int  nMatchedBlocks=0;
370         
371     bool sameHLTHd = 1;     
372     bool sameCDH = 1;     
373     bool sameRCU = 1;      
374     bool sameSize = 1;
375    
376     bool sameFlag = 1;
377     bool sameQMax = 1;
378     bool sameCharge = 1;
379     bool sameRow = 1;
380     bool sameFloat[4] = {1,1,1,1};
381     
382     for ( unsigned long ndx2 = ndx1+1; ndx2 < evtData.fBlockCnt; ndx2++ ){      
383       const AliHLTComponentBlockData* iter2 = blocks+ndx2;
384       if (  iter2->fDataType != (AliHLTTPCDefinitions::fgkHWClustersDataType | kAliHLTDataOriginTPC)
385             ) continue;
386             
387       int slice2 = AliHLTTPCDefinitions::GetMinSliceNr( *iter2 );
388       int patch2 = AliHLTTPCDefinitions::GetMinPatchNr( *iter2 );
389       if( slice1!=slice2 || patch1!=patch2 ) continue;
390       
391       if( checkedBlocks[ndx2] ) continue;
392       checkedBlocks[ndx2] = 1;
393
394       nMatchedBlocks++;
395
396       bool header2OK = 1;
397       if( iter2->fSize>0 ){
398         header2OK = ( iter2->fSize % sizeof(AliHLTUInt32_t) == 0 ) 
399           && ( iter2->fSize >= sizeof(AliRawDataHeader) ) 
400           && (iter2->fPtr != NULL );
401       }
402  
403       fHistHeaderAll->Fill(2);
404       sameHLTHd = ( header1OK == header2OK );
405       if( sameHLTHd ) fHistHeaderGood->Fill(2);
406       if( !header1OK || !header2OK ) continue;
407       
408       int nWordsHeader =  sizeof(AliRawDataHeader)/sizeof(AliHLTUInt32_t);
409       int nWords1 = iter1->fSize/sizeof(AliHLTUInt32_t);
410       int nWords2 = iter2->fSize/sizeof(AliHLTUInt32_t);
411       const AliHLTUInt32_t *p1 = (const AliHLTUInt32_t *) iter1->fPtr;
412       const AliHLTUInt32_t *p2 = (const AliHLTUInt32_t *) iter2->fPtr;
413
414       // compare CDH headers
415       
416        for(AliHLTInt32_t i=0; i<nWordsHeader; i++ ){
417         if( p1[i] != p2[i] ) sameCDH = 0;
418       }
419
420       fHistHeaderAll->Fill(3);
421       if( sameCDH ) fHistHeaderGood->Fill(3);
422
423       // find rcu headers
424       
425       int startRCU1 = nWordsHeader;
426       int startRCU2 = nWordsHeader;
427
428       for( ; startRCU1<nWords1; startRCU1+=6 ){
429         if( p1[startRCU1]>>30 == 0x2 ) break;     
430       }
431
432       for( ; startRCU2<nWords2; startRCU2+=6 ){
433         if( p2[startRCU2]>>30 == 0x2 ) break;     
434       } 
435       
436
437       // compare RCU headers
438
439       if( startRCU1 < nWords1 || startRCU2 < nWords2 ){
440         if( nWords1 - startRCU1 != nWords2 - startRCU2 ) sameRCU = 0;
441         else{
442           for( AliHLTInt32_t i1=startRCU1, i2=startRCU2; (i1<nWords1) && (i2<nWords2); i1++,i2++ ){
443             if( p1[i1]!=p2[i2] ) sameRCU = 0;
444           }
445         }
446         fHistHeaderAll->Fill(4);
447         if( sameRCU ) fHistHeaderGood->Fill(4);
448       }
449       
450       sameSize = ( startRCU1 == startRCU2 );
451
452       fHistHeaderAll->Fill(5);
453       if( sameSize ) fHistHeaderGood->Fill(5);
454
455       // compare clusters
456
457       if( startRCU1 == startRCU2 ){
458         for( AliHLTInt32_t i=nWordsHeader; i<startRCU1; i+=6){
459           AliHLTUInt32_t header1 = p1[i];
460           AliHLTUInt32_t charge1 = p1[i+1];
461           AliHLTUInt32_t header2 = p2[i];
462           AliHLTUInt32_t charge2 = p2[i+1];
463
464
465           AliHLTUInt32_t flag1 = header1 >> 30;
466           AliHLTUInt32_t qmax1 = ( header1 & 0xFFFFFF ) >> 6;
467           AliHLTUInt32_t row1 = (header1>>24) & 0x3f;
468           AliHLTUInt32_t flag2 = header2 >> 30;
469           AliHLTUInt32_t qmax2 = ( header2 & 0xFFFFFF ) >> 6;
470           AliHLTUInt32_t row2 = (header2>>24) & 0x3f;
471           
472           fHistClusterAll->Fill(0);
473           if( flag1 == flag2 ) fHistClusterGood->Fill(0);
474           else sameFlag = 0;
475           
476           if( flag1!=0x3 || flag2!=0x3 ) continue;
477           
478           // compare cluster qMax
479           
480           fHistClusterAll->Fill(1);
481           if( qmax1 == qmax2 ) fHistClusterGood->Fill(1);
482           else sameQMax = 0;
483           
484           // compare cluster row index
485           
486           fHistClusterAll->Fill(2);
487           if( row1 == row2 ) fHistClusterGood->Fill(2);
488           else sameRow = 0;
489           
490           // compare cluster charge
491           
492           fHistClusterAll->Fill(3);
493           if( charge1 == charge2 ) fHistClusterGood->Fill(3);
494           else sameCharge = 0;
495
496           // compare floating point data
497
498           for( int j=0; j<4; j++ ){
499             AliHLTFloat64_t f1     = *((AliHLTFloat32_t*)&p1[i+2+j]);
500             AliHLTFloat64_t f2     = *((AliHLTFloat32_t*)&p2[i+2+j]);
501             double w = (fabs(f1) + fabs(f2))/2;
502             if( w>1.e-20 ){
503               fHistClusterAll->Fill(4+j);
504               if( fabs(f1 - f2) < 1.e-6*w ) fHistClusterGood->Fill(4+j);
505               else sameFloat[j] = 0;
506             }
507           }
508         }
509       }
510     }
511         
512     fNBlocks++;
513       
514     bool err = 0;
515     TString warn;
516     warn.Form("HWCF consistency check for slice %d patch %d. Values  not matched:", slice1, patch1 );
517     
518     if( nMatchedBlocks!=1 ){ 
519       err=1; 
520       TString x;
521       x.Form(" NMatchedBlocks(%d)",nMatchedBlocks);
522       warn+=x;
523     }
524     if( !sameHLTHd ){ err=1; warn+=", HLT header";}
525     if( !sameCDH ){ err=1; warn+=", CDH header";}
526     if( !sameRCU ){ err=1; warn+=", RCU header";}
527     if( !sameSize ){ err=1; warn+=", N Clusters";}
528     if( !sameFlag ){ err=1; warn+=", Cluster Header";}
529     if( !sameQMax ){ err=1; warn+=", Cluster Qmax";}
530     if( !sameCharge ){ err=1; warn+=", Cluster Charge";}
531     if( !sameRow ){ err=1; warn+=", Cluster Row";}
532     if( !sameFloat[0] ){ err=1; warn+=", Cluster Pad value";}
533     if( !sameFloat[1] ){ err=1; warn+=", Cluster Time value";}
534     if( !sameFloat[2] ){ err=1; warn+=", Cluster PadErr value";}
535     if( !sameFloat[3] ){ err=1; warn+=", Cluster TimeErr value";}
536       
537     if( err ){
538       fNDismatch++;             
539       HLTWarning(warn.Data());
540     } else {
541       //warn+=" NO ";
542       //HLTWarning(warn.Data());     
543     }
544        
545     fHistHeaderAll->Fill(0);
546     if( nMatchedBlocks >=1 ) fHistHeaderGood->Fill(0);
547     fHistHeaderAll->Fill(1);
548     if( nMatchedBlocks <=1 ) fHistHeaderGood->Fill(1);
549
550   } // first block
551     
552   
553   delete[] checkedBlocks;
554    
555   HLTInfo("HWCF consistency check: %.10f %s of %ld data blocks are OK",
556              (double)(fNBlocks-fNDismatch)/(double)fNBlocks*100.,"%",fNBlocks);
557   
558   //if( fNDismatch>0 ){
559   //HLTWarning("HWCF inconsistency: %ld of %ld data blocks are not OK",
560   //fNDismatch,fNBlocks);      
561   //}  
562
563   fProfHeader->Divide(fHistHeaderGood, fHistHeaderAll,1,1.,"b");
564   fProfCluster->Divide(fHistClusterGood, fHistClusterAll, 1, 1, "b");
565
566   PushBack( (TObject*) fProfHeader, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC,0);
567   fBenchmark.AddOutput(GetLastObjectSize());
568   PushBack( (TObject*) fProfCluster, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC,0);
569   fBenchmark.AddOutput(GetLastObjectSize());
570   
571
572   fBenchmark.Stop(0);  
573   HLTInfo(fBenchmark.GetStatistics());
574   return iResult;
575 }
576
577