]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/HWCFemulator/AliHLTTPCHWCFEmulatorComponent.cxx
4cb5a66f35a30bd4efd453c9339fcadd2c24f579
[u/mrichter/AliRoot.git] / HLT / TPCLib / HWCFemulator / AliHLTTPCHWCFEmulatorComponent.cxx
1 // $Id$
2 //****************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project          * 
4 //* ALICE Experiment at CERN, All rights reserved.                           *
5 //*                                                                          *
6 //* Primary Authors: Sergey Gorbunov, Torsten Alt                            *
7 //* Developers:      Sergey Gorbunov <sergey.gorbunov@fias.uni-frankfurt.de> *
8 //*                  Torsten Alt <talt@cern.ch>                              *
9 //*                  for The ALICE HLT Project.                              *
10 //*                                                                          *
11 //* Permission to use, copy, modify and distribute this software and its     *
12 //* documentation strictly for non-commercial purposes is hereby granted     *
13 //* without fee, provided that the above copyright notice appears in all     *
14 //* copies and that both the copyright notice and this permission notice     *
15 //* appear in the supporting documentation. The authors make no claims       *
16 //* about the suitability of this software for any purpose. It is            *
17 //* provided "as is" without express or implied warranty.                    *
18 //****************************************************************************
19
20 //  @file   AliHLTTPCHWCFEmulatorComponent.cxx
21 //  @author Sergey Gorbunov <sergey.gorbunov@fias.uni-frankfurt.de>
22 //  @author Torsten Alt <talt@cern.ch> 
23 //  @brief  HLT Component interface for for FPGA ClusterFinder Emulator for TPC
24 //  @brief  ( see AliHLTTPCHWCFEmulator class )
25 //  @note
26
27
28 #if __GNUC__>= 3
29 using namespace std;
30 #endif
31 #include "AliHLTTPCHWCFEmulatorComponent.h"
32
33 #include "AliHLTTPCDefinitions.h"
34 #include "AliHLTTPCHWCFDataTypes.h"
35 #include "AliHLTTPCClusterMCData.h"
36 #include "AliHLTTPCHWCFData.h"
37
38 #include "AliGRPObject.h"
39 #include "AliCDBEntry.h"
40 #include "AliCDBManager.h"
41 #include "AliRawDataHeader.h"
42 #include <cstdlib>
43 #include <cerrno>
44 #include <memory>
45 #include "TString.h"
46 #include "TObjString.h"
47 #include "TObjArray.h"
48
49
50 #include <sys/time.h>
51 #include "TFile.h"
52
53 AliHLTTPCHWCFEmulatorComponent::AliHLTTPCHWCFEmulatorComponent()
54   :
55   AliHLTProcessor(),
56   fDoDeconvTime(0),
57   fDoDeconvPad(0),
58   fDoMC(0),
59   fDoFlowControl(0),
60   fDoSinglePadSuppression(0),
61   fBypassMerger(0),
62   fClusterLowerLimit(0),
63   fSingleSeqLimit(0),
64   fMergerDistance(0),
65   fUseTimeBinWindow(0),
66   fUseTimeFollow(0),
67   fChargeFluctuation(0),
68   fDebug(0),
69   fCFSupport(),
70   fCFEmulator(),
71   fBenchmark("TPCHWClusterFinderEmulator")
72 {
73   // see header file for class documentation
74   // or
75   // refer to README to build package
76   // or
77   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
78 }
79
80
81 AliHLTTPCHWCFEmulatorComponent::AliHLTTPCHWCFEmulatorComponent(const AliHLTTPCHWCFEmulatorComponent&)
82   :
83   AliHLTProcessor(),
84   fDoDeconvTime(0),
85   fDoDeconvPad(0),
86   fDoMC(0),
87   fDoFlowControl(0),
88   fDoSinglePadSuppression(0),
89   fBypassMerger(0),
90   fClusterLowerLimit(0),
91   fSingleSeqLimit(0),
92   fMergerDistance(0),
93   fUseTimeBinWindow(0),
94   fUseTimeFollow(0),
95   fChargeFluctuation(0),
96   fDebug(0),
97   fCFSupport(),
98   fCFEmulator(),
99   fBenchmark("TPCHWClusterFinderEmulator")
100 {
101   // dummy
102 }
103
104 AliHLTTPCHWCFEmulatorComponent& AliHLTTPCHWCFEmulatorComponent::operator=(const AliHLTTPCHWCFEmulatorComponent&)
105 {
106   // dummy
107   return *this;
108 }
109
110 AliHLTTPCHWCFEmulatorComponent::~AliHLTTPCHWCFEmulatorComponent()
111 {
112   // see header file for class documentation
113 }
114
115 // Public functions to implement AliHLTComponent's interface.
116 // These functions are required for the registration process
117
118 const char* AliHLTTPCHWCFEmulatorComponent::GetComponentID()
119 {
120   // see header file for class documentation
121   return "TPCHWClusterFinderEmulator";
122 }
123
124 void AliHLTTPCHWCFEmulatorComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
125 {
126   // see header file for class documentation
127   list.clear();
128   list.push_back( AliHLTTPCDefinitions::fgkUnpackedRawDataType );        
129   list.push_back( kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC );
130 }
131
132 AliHLTComponentDataType AliHLTTPCHWCFEmulatorComponent::GetOutputDataType()
133 {
134   // see header file for class documentation
135   return kAliHLTMultipleDataType;
136 }
137
138 int AliHLTTPCHWCFEmulatorComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
139 {
140   // see header file for class documentation
141   tgtList.clear();
142   tgtList.push_back(AliHLTTPCDefinitions::fgkHWClustersDataType | kAliHLTDataOriginTPC );
143   tgtList.push_back(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo | kAliHLTDataOriginTPC );
144   return tgtList.size();
145 }
146
147
148 void AliHLTTPCHWCFEmulatorComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
149 {
150   // see header file for class documentation
151   // XXX TODO: Find more realistic values.  
152   constBase = 0;
153   inputMultiplier = (6 * 0.7);
154 }
155
156
157 AliHLTComponent* AliHLTTPCHWCFEmulatorComponent::Spawn()
158 {
159   // see header file for class documentation
160   return new AliHLTTPCHWCFEmulatorComponent();
161 }
162
163 void AliHLTTPCHWCFEmulatorComponent::GetOCDBObjectDescription( TMap* const targetMap){
164 // Get a list of OCDB object description needed for the particular component
165   
166   if (!targetMap) return;
167   
168   // OCDB entries for component arguments
169   targetMap->Add(new TObjString("HLT/ConfigTPC/TPCHWClusterFinder"), new TObjString("component arguments, empty at the moment"));
170 }
171
172
173 int AliHLTTPCHWCFEmulatorComponent::DoInit( int argc, const char** argv )
174 {
175   // see header file for class documentation
176
177   TString arguments = "";
178   for ( int i = 0; i < argc; i++ ) {
179     if ( !arguments.IsNull() ) arguments += " ";
180     arguments += argv[i];
181   }
182
183   return Configure( NULL, NULL, arguments.Data()  );
184 }
185
186 int AliHLTTPCHWCFEmulatorComponent::Reconfigure( const char* cdbEntry, const char* chainId )
187 {
188   // Reconfigure the component from OCDB
189
190   return Configure( cdbEntry, chainId, NULL );
191 }
192
193 int AliHLTTPCHWCFEmulatorComponent::ScanConfigurationArgument(int argc, const char** argv)
194 {
195   // see header file for class documentation
196   TString arguments = "";
197   for ( int i = 0; i < argc; i++ ) {
198     if ( !arguments.IsNull() ) arguments += " ";
199     arguments += argv[i];
200   }
201   return ReadConfigurationString(arguments);
202 }
203
204
205
206 void AliHLTTPCHWCFEmulatorComponent::SetDefaultConfiguration()
207 {
208   // Set default configuration for the FPGA ClusterFinder Emulator component
209   // Some parameters can be later overwritten from the OCDB
210
211   fDoDeconvTime = 1;
212   fDoDeconvPad = 1;
213   fDoMC = 0;
214   fDoFlowControl = 0;
215   fDoSinglePadSuppression = 0;
216   fBypassMerger = 0;
217   fClusterLowerLimit = 10;
218   fSingleSeqLimit = 0;
219   fMergerDistance = 4;
220   fUseTimeBinWindow = 1;
221   fUseTimeFollow = 1;
222   fChargeFluctuation = 0;
223   fDebug = 0;
224   fBenchmark.Reset();
225   fBenchmark.SetTimer(0,"total");
226   fBenchmark.SetTimer(1,"reco");    
227 }
228
229 int AliHLTTPCHWCFEmulatorComponent::ReadConfigurationString(  const char* arguments )
230 {
231   // Set configuration parameters for the FPGA ClusterFinder Emulator component
232   // from the string
233
234   int iResult = 0;
235   if ( !arguments ) return iResult;
236   //cout<<"["<<arguments<<"]"<<endl;
237   TString allArgs = arguments;
238   TString argument;
239   int bMissingParam = 0;
240
241   TObjArray* pTokens = allArgs.Tokenize( " " );
242
243   int nArgs =  pTokens ? pTokens->GetEntries() : 0;
244
245   for ( int i = 0; i < nArgs; i++ ) {
246     argument = ( ( TObjString* )pTokens->At( i ) )->GetString();
247     if ( argument.IsNull() ) continue;
248
249     if ( argument.CompareTo( "-deconvolute-time" ) == 0 ) {
250       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
251       fDoDeconvTime  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
252       fUseTimeBinWindow = fDoDeconvTime;
253       HLTInfo( "Time deconvolution and using of TimeBin window are set to: %d", fDoDeconvTime );
254       continue;
255     }
256
257     if ( argument.CompareTo( "-deconvolute-pad" ) == 0 ) {
258       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
259       fDoDeconvPad  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
260       HLTInfo( "Pad deconvolution is set to: %d", fDoDeconvPad );
261       continue;
262     }
263
264     if ( argument.CompareTo( "-deconvolute" ) == 0 ) {
265       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
266       fDoDeconvTime  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
267       fUseTimeBinWindow = fDoDeconvTime;
268       fDoDeconvPad  = fDoDeconvTime;
269       HLTInfo( "Time and pad deconvolution and using of TimeBin window are set to: %d", fDoDeconvPad );
270       continue;
271     }
272  
273     if ( argument.CompareTo( "-do-mc" ) == 0 ) {
274       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
275       fDoMC  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
276       HLTInfo( "MC processing is set to: %d", fDoMC );
277       continue;
278     }
279
280     if ( argument.CompareTo( "-flow-control" ) == 0 ) {
281       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
282       fDoFlowControl  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
283       HLTInfo( "Flow control is set to: %d", fDoFlowControl );
284       continue;
285     }
286
287     if ( argument.CompareTo( "-single-pad-suppression" ) == 0 ) {
288       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
289       fDoSinglePadSuppression  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
290       HLTInfo( "Single pad suppression is set to: %d", fDoSinglePadSuppression );
291       continue;
292     }
293     
294     if ( argument.CompareTo( "-bypass-merger" ) == 0 ) {
295       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
296       fBypassMerger  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
297       HLTInfo( "Bypassing merger is set to: %d", fBypassMerger );
298       continue;
299     }
300
301     if ( argument.CompareTo( "-cluster-lower-limit" ) == 0 ) {
302       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
303       fClusterLowerLimit  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
304       HLTInfo( "Cluster lower limit is set to: %d", fClusterLowerLimit );
305       continue;
306     }
307
308     if ( argument.CompareTo( "-single-sequence-limit" ) == 0 ) {
309       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
310       fSingleSeqLimit  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
311       HLTInfo( "Single sequence limit is set to: %d", fSingleSeqLimit );
312       continue;
313     }
314     
315     if ( argument.CompareTo( "-merger-distance" ) == 0 ) {
316       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
317       fMergerDistance  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
318       HLTInfo( "Merger distance is set to: %d", fMergerDistance );
319       continue;
320     }
321  
322     if ( argument.CompareTo( "-use-timebin-window" ) == 0 ) {
323       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
324       fUseTimeBinWindow  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
325       fDoDeconvTime = fUseTimeBinWindow;
326       HLTInfo( "Using TimeBin window and Time deconvolution are set to: %d", fUseTimeBinWindow );
327       continue;
328     }
329    
330     if ( argument.CompareTo( "-charge-fluctuation" ) == 0 ) {
331       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
332       fChargeFluctuation  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
333       HLTInfo( "Charge fluctuation is set to: %d", fChargeFluctuation );
334       continue;
335     }
336     
337     if ( argument.CompareTo( "-use-time-follow" ) == 0 ) {
338       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
339       fUseTimeFollow  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
340       HLTInfo( "Use of time follow algorithm is set to: %d", fUseTimeFollow );
341       continue;
342     }
343    
344     if ( argument.CompareTo( "-debug-level" ) == 0 ) {
345       if ( ( bMissingParam = ( ++i >= pTokens->GetEntries() ) ) ) break;
346       fDebug  = ( ( TObjString* )pTokens->At( i ) )->GetString().Atoi();
347       HLTInfo( "Debug level is set to: %d", fDebug );
348       continue;
349     }
350     
351     HLTError( "Unknown option \"%s\"", argument.Data() );
352     iResult = -EINVAL;
353   }
354   delete pTokens;
355
356   if ( bMissingParam ) {
357     HLTError( "Specifier missed for parameter \"%s\"", argument.Data() );
358     iResult = -EINVAL;
359   }
360
361   return iResult;
362 }
363
364
365 int AliHLTTPCHWCFEmulatorComponent::ReadCDBEntry( const char* cdbEntry, const char* chainId )
366 {
367   // Read configuration from OCDB
368
369   const char* defaultNotify = "";
370
371   if ( !cdbEntry ){
372     cdbEntry = "HLT/ConfigTPC/TPCHWClusterFinder";
373     defaultNotify = " (default)";
374     chainId = 0;
375   }
376
377   HLTInfo( "configure from entry \"%s\"%s, chain id %s", cdbEntry, defaultNotify, ( chainId != NULL && chainId[0] != 0 ) ? chainId : "<none>" );
378   AliCDBEntry *pEntry = AliCDBManager::Instance()->Get( cdbEntry );//,GetRunNo());
379   
380   if ( !pEntry ) {
381     HLTError( "cannot fetch object \"%s\" from CDB", cdbEntry );
382     return -EINVAL;
383   }
384
385   TObjString* pString = dynamic_cast<TObjString*>( pEntry->GetObject() );
386
387   if ( !pString ) {
388     HLTError( "configuration object \"%s\" has wrong type, required TObjString", cdbEntry );
389     return -EINVAL;
390   }
391
392   HLTInfo( "received configuration object string: \"%s\"", pString->GetString().Data() );
393
394   return  ReadConfigurationString( pString->GetString().Data() );
395 }
396
397
398 int AliHLTTPCHWCFEmulatorComponent::Configure( const char* cdbEntry, const char* chainId, const char *commandLine )
399 {
400   // Configure the component
401   // There are few levels of configuration,
402   // parameters which are set on one step can be overwritten on the next step
403
404   //* read hard-coded values
405
406   SetDefaultConfiguration();
407
408   //* read the default CDB entry
409
410   int iResult1 = ReadCDBEntry( NULL, chainId );
411
412   //* read the actual CDB entry if required
413
414   int iResult2 = ( cdbEntry ) ? ReadCDBEntry( cdbEntry, chainId ) : 0;
415
416   //* read extra parameters from input (if they are)
417
418   int iResult3 = 0;
419
420   if ( commandLine && commandLine[0] != '\0' ) {
421     HLTInfo( "received configuration string from HLT framework: \"%s\"", commandLine );
422     iResult3 = ReadConfigurationString( commandLine );
423   }
424   
425   if( fDebug>1 ) fCFEmulator.SetDebugLevel( fDebug );
426   else fCFEmulator.SetDebugLevel(0);
427
428   return iResult1 ? iResult1 : ( iResult2 ? iResult2 : iResult3 );
429 }
430
431
432 int AliHLTTPCHWCFEmulatorComponent::DoDeinit()
433 {
434   // see header file for class documentation 
435   return 0;
436 }
437
438
439 int AliHLTTPCHWCFEmulatorComponent::DoEvent( const AliHLTComponentEventData& evtData, 
440                                                         const AliHLTComponentBlockData* blocks, 
441                                                         AliHLTComponentTriggerData& /*trigData*/, AliHLTUInt8_t* outputPtr, 
442                                                         AliHLTUInt32_t& size, 
443                                                         vector<AliHLTComponentBlockData>& outputBlocks )
444 {
445   // see header file for class documentation
446
447   int iResult=0;
448   AliHLTUInt32_t maxSize = size;
449   size = 0;
450   
451   if(!IsDataEvent()){
452     return 0;
453   }
454
455   fBenchmark.StartNewEvent();
456   fBenchmark.Start(0);
457
458   AliHLTUInt32_t configWord1=0, configWord2=0; 
459   AliHLTTPCHWCFEmulator::CreateConfiguration
460     ( fDoDeconvTime, fDoDeconvPad, fDoFlowControl, fDoSinglePadSuppression, fBypassMerger, fClusterLowerLimit, fSingleSeqLimit, fMergerDistance, fUseTimeBinWindow, fChargeFluctuation, fUseTimeFollow, configWord1, configWord2 );
461
462   for ( unsigned long ndx = 0; ndx < evtData.fBlockCnt; ndx++ )
463     {
464       const AliHLTComponentBlockData* iter = blocks+ndx;
465       
466       if (  iter->fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC) 
467             &&  iter->fDataType != AliHLTTPCDefinitions::fgkUnpackedRawDataType ) continue;
468
469       int slice = AliHLTTPCDefinitions::GetMinSliceNr( *iter );
470       int patch = AliHLTTPCDefinitions::GetMinPatchNr( *iter );
471  
472       fBenchmark.AddInput(iter->fSize);
473  
474       if (!iter->fPtr) continue;
475  
476       // create input block for the HW cluster finder
477
478       const AliHLTUInt32_t *rawEvent=0;
479       const AliHLTTPCClusterMCLabel *mcLabels = 0;
480       AliHLTUInt32_t rawEventSize32 = 0;
481       AliHLTUInt32_t nMCLabels = 0;
482
483       if( fCFSupport.CreateRawEvent( iter, rawEvent, rawEventSize32, mcLabels, nMCLabels )<0 ) continue; 
484       if( !fDoMC ){
485         mcLabels = 0;
486         nMCLabels = 0;
487       }
488
489       // book memory for the output
490       
491       AliHLTUInt32_t maxNClusters = rawEventSize32 + 1; // N 32-bit words in input
492       AliHLTUInt32_t clustersSize32 = maxNClusters*AliHLTTPCHWCFData::fgkAliHLTTPCHWClusterSize;
493       AliHLTUInt32_t nOutputMC = maxNClusters;
494
495       AliHLTUInt32_t headerSize = sizeof(AliRawDataHeader);
496       AliHLTUInt32_t outBlockSize=headerSize+clustersSize32*sizeof(AliHLTUInt32_t);
497       std::auto_ptr<AliHLTUInt8_t> outBlock(new AliHLTUInt8_t[outBlockSize]);
498       std::auto_ptr<AliHLTTPCClusterMCLabel> allocOutMC(new AliHLTTPCClusterMCLabel[nOutputMC+1]);
499       
500       if( !outBlock.get() || !allocOutMC.get() ){
501         iResult=-ENOMEM;
502         break;
503       }
504
505       memset(outBlock.get(), 0, outBlockSize*sizeof(AliHLTUInt8_t));
506       memset(allocOutMC.get(), 0, (nOutputMC+1)*sizeof(AliHLTTPCClusterMCLabel));
507       AliHLTTPCClusterMCData *outMC = reinterpret_cast<AliHLTTPCClusterMCData *>(allocOutMC.get());
508       
509       // fill CDH header here, since the HW clusterfinder does not receive it
510       
511       AliRawDataHeader *cdhHeader = reinterpret_cast<AliRawDataHeader*>(iter->fPtr);
512       AliRawDataHeader *outCDHHeader = reinterpret_cast<AliRawDataHeader*>(outBlock.get());      
513       *outCDHHeader = *cdhHeader;
514       outCDHHeader->fSize = 0xFFFFFFFF;
515
516       AliHLTUInt32_t *outClusters = reinterpret_cast<AliHLTUInt32_t*> (outBlock.get() + headerSize);
517      
518       fBenchmark.Start(1);
519       fCFEmulator.Init
520         ( fCFSupport.GetMapping(slice,patch), configWord1, configWord2 );
521       
522       int err = fCFEmulator.FindClusters( rawEvent, rawEventSize32, 
523                                           outClusters, clustersSize32, 
524                                           mcLabels, nMCLabels,
525                                           outMC );
526       fBenchmark.Stop(1);
527       if( err==-1 ){ HLTWarning("NULL input pointer (warning %d)",err);}
528       else if( err==-2 ){  HLTWarning("No space left in the output buffer (warning %d)",err); }
529       else if( err<0 ){ HLTWarning("HWCF emulator finished with error code %d",err); }
530       if( err<0 ){
531         continue;
532       }
533
534       if( fDebug ){
535         int elsize=AliHLTTPCHWCFData::fgkAliHLTTPCHWClusterSize;
536         printf("\nHWCF Emulator: output clusters for slice%d patch %d:\n",slice,patch);
537         for( AliHLTUInt32_t i=0; i<clustersSize32; i+=elsize ){
538           AliHLTUInt32_t *c = outClusters+i;
539           AliHLTUInt32_t flag = (c[0]>>30);       
540           if( flag == 0x3){ //beginning of a cluster
541             int padRow  = (c[0]>>24)&0x3f;
542             int q  = c[1];
543             double p   = *((AliHLTFloat32_t*)&c[2]);
544             double t  = *((AliHLTFloat32_t*)&c[3]);
545             AliHLTFloat32_t p2 = *((AliHLTFloat32_t*)&c[4]);
546             AliHLTFloat32_t t2 = *((AliHLTFloat32_t*)&c[5]);
547             printf("N: %3d    R: %3d    C: %4d    P:  %7.4f    T:  %8.4f    DP: %6.4f    DT: %6.4f\n", 
548                    i/elsize+1, padRow, q, p, t, sqrt(fabs(p2-p*p)), sqrt(fabs(t2-t*t)));
549
550             if( outMC && outMC->fCount>0 ){
551               printf("        MC: (%3d,%6.1f) (%3d,%6.1f) (%3d,%6.1f)\n",
552                      outMC->fLabels[i/elsize].fClusterID[0].fMCID,outMC->fLabels[i/elsize].fClusterID[0].fWeight,
553                      outMC->fLabels[i/elsize].fClusterID[1].fMCID,outMC->fLabels[i/elsize].fClusterID[1].fWeight,
554                      outMC->fLabels[i/elsize].fClusterID[2].fMCID,outMC->fLabels[i/elsize].fClusterID[2].fWeight
555                      );
556             }
557           }
558         }
559       }
560           
561
562       AliHLTUInt32_t outSize = headerSize + clustersSize32*sizeof(AliHLTUInt32_t);
563       
564       if( size + outSize <= maxSize ){
565         
566         memcpy( outputPtr, outBlock.get(), outSize );
567         
568         AliHLTComponentBlockData bd;
569         FillBlockData( bd );
570         bd.fOffset = size;
571         bd.fSize = outSize;
572         bd.fSpecification = iter->fSpecification;
573         bd.fDataType = AliHLTTPCDefinitions::fgkHWClustersDataType | kAliHLTDataOriginTPC;
574         outputBlocks.push_back( bd );
575         fBenchmark.AddOutput(bd.fSize);
576         size+= bd.fSize;
577         outputPtr+=bd.fSize;
578       } else {
579         HLTWarning( "Output buffer (%db) is too small, required %db", maxSize, size+outSize);
580         iResult=-ENOSPC;
581         break;
582       }
583
584       if( fDoMC && outMC && outMC->fCount>0 ){
585         int s = sizeof(AliHLTTPCClusterMCData) + outMC->fCount*sizeof(AliHLTTPCClusterMCLabel);
586         if( size + s <= maxSize ){
587           memcpy( outputPtr, outMC, s );                  
588           AliHLTComponentBlockData bdMCInfo;
589           FillBlockData( bdMCInfo );
590           bdMCInfo.fOffset = size;
591           bdMCInfo.fSize = s;
592           bdMCInfo.fSpecification = iter->fSpecification;
593           bdMCInfo.fDataType = AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo | kAliHLTDataOriginTPC;
594           outputBlocks.push_back( bdMCInfo );
595           fBenchmark.AddOutput(bdMCInfo.fSize);
596           size+=bdMCInfo.fSize;
597           outputPtr+=bdMCInfo.fSize; 
598         } else {        
599           HLTWarning( "Output buffer (%db) is too small, required %db", maxSize, size+s);
600           iResult=-ENOSPC;          
601           break;
602         }
603       }
604     }
605   
606   fBenchmark.Stop(0);  
607   HLTInfo(fBenchmark.GetStatistics());
608   fCFSupport.ReleaseEventMemory();
609   return iResult;
610 }
611
612