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