]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/comp/AliHLTTPCDataCompressionComponent.cxx
moving decoding of compressed TPC data to generic class
[u/mrichter/AliRoot.git] / HLT / TPCLib / comp / AliHLTTPCDataCompressionComponent.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: Matthias Richter <Matthias.Richter@ift.uib.no>        *
7 //*                  for The ALICE HLT Project.                            *
8 //*                                                                        *
9 //* Permission to use, copy, modify and distribute this software and its   *
10 //* documentation strictly for non-commercial purposes is hereby granted   *
11 //* without fee, provided that the above copyright notice appears in all   *
12 //* copies and that both the copyright notice and this permission notice   *
13 //* appear in the supporting documentation. The authors make no claims     *
14 //* about the suitability of this software for any purpose. It is          *
15 //* provided "as is" without express or implied warranty.                  *
16 //**************************************************************************
17
18 /// @file   AliHLTTPCDataCompressionComponent.cxx
19 /// @author Matthias Richter
20 /// @date   2011-08-08
21 /// @brief  TPC component for data compression
22 ///
23
24 #include "AliHLTTPCDataCompressionComponent.h"
25 #include "AliHLTTPCDefinitions.h"
26 #include "AliHLTTPCTrackGeometry.h"
27 #include "AliHLTTPCSpacePointContainer.h"
28 #include "AliHLTTPCHWCFSpacePointContainer.h"
29 #include "AliHLTGlobalBarrelTrack.h"
30 #include "AliHLTComponentBenchmark.h"
31 #include "AliHLTDataDeflaterSimple.h"
32 #include "AliHLTDataDeflaterHuffman.h"
33 #include "AliHLTTPCTransform.h"
34 #include "AliHLTTPCClusterMCData.h"
35 #include "AliHLTTPCClusterTransformation.h"
36 #include "AliRawDataHeader.h"
37 #include "AliCDBManager.h"
38 #include "AliCDBPath.h"
39 #include "AliCDBId.h"
40 #include "AliCDBMetaData.h"
41 #include "AliCDBEntry.h"
42 #include "TH1F.h"
43 #include "TFile.h"
44 #include <memory>
45
46 ClassImp(AliHLTTPCDataCompressionComponent)
47
48 AliHLTTPCDataCompressionComponent::AliHLTTPCDataCompressionComponent()
49   : AliHLTProcessor()
50   , fMode(0)
51   , fDeflaterMode(0)
52   , fVerificationMode(0)
53   , fMaxDeltaPad(AliHLTTPCDefinitions::GetMaxClusterDeltaPad())
54   , fMaxDeltaTime(AliHLTTPCDefinitions::GetMaxClusterDeltaTime())
55   , fRawInputClusters(NULL)
56   , fInputClusters(NULL)
57   , fTrackGrid(NULL)
58   , fSpacePointGrid(NULL)
59   , fpDataDeflater(NULL)
60   , fHistoCompFactor(NULL)
61   , fHistoResidualPad(NULL)
62   , fHistoResidualTime(NULL)
63   , fHistoClustersOnTracks(NULL)
64   , fHistoClusterRatio(NULL)
65   , fHistoTrackClusterRatio(NULL)
66   , fHistogramFile()
67   , fTrainingTableOutput()
68   , fpBenchmark(NULL)
69   , fpWrittenAssociatedClusterIds(NULL)
70   , fDriftTimeFactorA(1.)
71   , fDriftTimeOffsetA(0.)
72   , fDriftTimeFactorC(1.)
73   , fDriftTimeOffsetC(0.)
74   , fVerbosity(0)
75 {
76 }
77
78 AliHLTTPCDataCompressionComponent::~AliHLTTPCDataCompressionComponent()
79 {
80   /// destructor
81   if (fpWrittenAssociatedClusterIds) delete fpWrittenAssociatedClusterIds;
82 }
83
84
85 const char* AliHLTTPCDataCompressionComponent::GetComponentID()
86 {
87   /// inherited from AliHLTComponent: id of the component
88   return "TPCDataCompressor";
89 }
90
91
92 void AliHLTTPCDataCompressionComponent::GetInputDataTypes( AliHLTComponentDataTypeList& tgtList)
93 {
94   /// inherited from AliHLTComponent: list of data types in the vector reference
95   tgtList.clear();
96   tgtList.push_back(AliHLTTPCDefinitions::fgkHWClustersDataType);
97   tgtList.push_back(AliHLTTPCDefinitions::fgkClustersDataType);
98   tgtList.push_back(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
99 }
100
101 AliHLTComponentDataType AliHLTTPCDataCompressionComponent::GetOutputDataType()
102 {
103   /// inherited from AliHLTComponent: output data type of the component.
104   return kAliHLTMultipleDataType;
105 }
106
107 int AliHLTTPCDataCompressionComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
108 {
109   /// inherited from AliHLTComponent: multiple output data types of the component.
110   tgtList.clear();
111   tgtList.push_back(AliHLTTPCDefinitions::RawClustersDataType());
112   tgtList.push_back(AliHLTTPCDefinitions::RemainingClustersCompressedDataType());
113   tgtList.push_back(AliHLTTPCDefinitions::RemainingClusterIdsDataType());
114   tgtList.push_back(AliHLTTPCDefinitions::ClusterTracksCompressedDataType());
115   tgtList.push_back(AliHLTTPCDefinitions::ClusterIdTracksDataType());
116   return tgtList.size();
117 }
118
119 void AliHLTTPCDataCompressionComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
120 {
121   /// inherited from AliHLTComponent: output data size estimator
122   constBase=0;
123   inputMultiplier=1.;  // there should not be more data than input
124   inputMultiplier+=.3; // slightly more data when using the old HWCF data with 20 Byte and raw clusters 22 Byte
125   if (fpWrittenAssociatedClusterIds) inputMultiplier+=.3; // space for optional cluster id array
126 }
127
128 AliHLTComponent* AliHLTTPCDataCompressionComponent::Spawn()
129 {
130   /// inherited from AliHLTComponent: spawn function.
131   return new AliHLTTPCDataCompressionComponent;
132 }
133
134 int AliHLTTPCDataCompressionComponent::DoEvent( const AliHLTComponentEventData& evtData, 
135                                                 const AliHLTComponentBlockData* inputBlocks, 
136                                                 AliHLTComponentTriggerData& /*trigData*/,
137                                                 AliHLTUInt8_t* outputPtr,
138                                                 AliHLTUInt32_t& size,
139                                                 AliHLTComponentBlockDataList& outputBlocks )
140 {
141   /// inherited from AliHLTProcessor: data processing
142   int iResult=0;
143   AliHLTUInt32_t capacity=size;
144   size=0;
145
146   if (!IsDataEvent()) return 0;
147
148   if (!fRawInputClusters) {
149     return -ENODEV;
150   }
151
152   if (GetBenchmarkInstance()) {
153     GetBenchmarkInstance()->StartNewEvent();
154     GetBenchmarkInstance()->Start(0);
155   }
156
157   // Process an event
158   // Loop over all input blocks in the event
159   bool bHaveMC=(GetFirstInputBlock(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo | kAliHLTDataOriginTPC))!=NULL;
160   if ((bHaveMC || fVerificationMode>0) && fpWrittenAssociatedClusterIds==NULL) {
161     fpWrittenAssociatedClusterIds=new vector<AliHLTUInt32_t>;
162   }
163
164   const AliHLTComponentBlockData* pDesc=NULL;
165
166   AliHLTUInt8_t minSlice=0xFF, maxSlice=0xFF, minPatch=0xFF, maxPatch=0xFF;
167   AliHLTUInt32_t inputRawClusterSize=0;
168   AliHLTUInt32_t outputDataSize=0;
169   int allClusters=0;
170   int associatedClusters=0;
171
172   /// input track array
173   vector<AliHLTGlobalBarrelTrack> inputTrackArray;
174
175   if (GetBenchmarkInstance()) {
176     GetBenchmarkInstance()->Start(2);
177   }
178
179   // transformed clusters
180   if (fMode==10) { // FIXME: condition to be adjusted
181     for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType);
182          pDesc!=NULL; pDesc=GetNextInputBlock()) {
183       if (GetBenchmarkInstance()) {
184         GetBenchmarkInstance()->AddInput(pDesc->fSize);
185       }
186       AliHLTUInt8_t slice = 0;
187       AliHLTUInt8_t patch = 0;
188       slice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification );
189       patch = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification );
190       if ( minSlice==0xFF || slice<minSlice )   minSlice = slice;
191       if ( maxSlice==0xFF || slice>maxSlice )   maxSlice = slice;
192       if ( minPatch==0xFF || patch<minPatch )   minPatch = patch;
193       if ( maxPatch==0xFF || patch>maxPatch )   maxPatch = patch;
194       if (fInputClusters) {
195         fInputClusters->AddInputBlock(pDesc);
196       }
197     }
198     if (GetBenchmarkInstance()) {
199       GetBenchmarkInstance()->Stop(2);
200       GetBenchmarkInstance()->Start(3);
201     }
202   }
203
204   // track data input
205   if (fMode==2) { // FIXME: condition to be adjusted
206     for (pDesc=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
207          pDesc!=NULL; pDesc=GetNextInputBlock()) {
208       if (GetBenchmarkInstance()) {
209         GetBenchmarkInstance()->AddInput(pDesc->fSize);
210       }
211       AliHLTUInt8_t slice = 0;
212       AliHLTUInt8_t patch = 0;
213       slice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification );
214       patch = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification );
215       if ( minSlice==0xFF || slice<minSlice )   minSlice = slice;
216       if ( maxSlice==0xFF || slice>maxSlice )   maxSlice = slice;
217       if ( minPatch==0xFF || patch<minPatch )   minPatch = patch;
218       if ( maxPatch==0xFF || patch>maxPatch )   maxPatch = patch;
219       const AliHLTTracksData* pTracks=reinterpret_cast<const AliHLTTracksData*>(pDesc->fPtr);
220       if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(pTracks, pDesc->fSize, inputTrackArray))<0) {
221         return iResult;
222       }
223     }
224   }
225
226   if (GetBenchmarkInstance()) {
227     GetBenchmarkInstance()->Stop(3);
228     GetBenchmarkInstance()->Start(4);
229   }
230
231   // processing
232   for (vector<AliHLTGlobalBarrelTrack>::iterator track=inputTrackArray.begin();
233        track!=inputTrackArray.end();
234        track++) {
235     int trackID=track->GetID();
236     if (trackID<0) {
237       // FIXME: error guard
238       HLTError("invalid track ID");
239       continue;
240     }
241
242     if (fVerbosity>0) {
243       UInt_t nofPoints=track->GetNumberOfPoints();
244       const UInt_t* points=track->GetPoints();
245       for (unsigned i=0; i<nofPoints; i++) {
246         int slice=AliHLTTPCSpacePointData::GetSlice(points[i]);
247         int partition=AliHLTTPCSpacePointData::GetPatch(points[i]);
248         int number=AliHLTTPCSpacePointData::GetNumber(points[i]);
249         HLTInfo("track %d point %d id 0x%08x slice %d partition %d number %d", track->GetID(), i, points[i], slice, partition, number);
250       }
251     }
252
253     AliHLTTPCTrackGeometry* trackpoints=new AliHLTTPCTrackGeometry;
254     if (!trackpoints) continue;
255     trackpoints->InitDriftTimeTransformation(fDriftTimeFactorA, fDriftTimeOffsetA, fDriftTimeFactorC, fDriftTimeOffsetC);
256     trackpoints->SetTrackId(trackID);
257     trackpoints->CalculateTrackPoints(*track);
258     trackpoints->RegisterTrackPoints(fTrackGrid);
259     track->SetTrackGeometry(trackpoints);
260   }
261
262   for (vector<AliHLTGlobalBarrelTrack>::const_iterator track=inputTrackArray.begin();
263        track!=inputTrackArray.end();
264        track++) {
265     AliHLTTrackGeometry* trackpoints=track->GetTrackGeometry();
266     if (!trackpoints) continue;
267     trackpoints->FillTrackPoints(fTrackGrid);
268   }
269   if (fVerbosity>0) {
270     fTrackGrid->Print();
271   }
272
273   if (GetBenchmarkInstance()) {
274     GetBenchmarkInstance()->Stop(4);
275     GetBenchmarkInstance()->Start(5);
276   }
277
278   // loop over raw cluster blocks, assign to tracks and write
279   // unassigned clusters
280   for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkHWClustersDataType);
281        pDesc!=NULL; pDesc=GetNextInputBlock()) {
282     if (pDesc->fSize<=sizeof(AliRawDataHeader)) continue;
283     if (GetBenchmarkInstance()) {
284       GetBenchmarkInstance()->Start(1);
285       GetBenchmarkInstance()->AddInput(pDesc->fSize);
286     }
287     AliHLTUInt8_t slice = 0;
288     AliHLTUInt8_t patch = 0;
289     slice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification );
290     patch = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification );
291     if ( minSlice==0xFF || slice<minSlice )     minSlice = slice;
292     if ( maxSlice==0xFF || slice>maxSlice )     maxSlice = slice;
293     if ( minPatch==0xFF || patch<minPatch )     minPatch = patch;
294     if ( maxPatch==0xFF || patch>maxPatch )     maxPatch = patch;
295     inputRawClusterSize+=pDesc->fSize;
296
297     // add the data and populate the index grid
298     fRawInputClusters->AddInputBlock(pDesc);
299     fRawInputClusters->PopulateAccessGrid(fSpacePointGrid, pDesc->fSpecification);
300     if (fVerbosity>0 && fSpacePointGrid->GetNumberOfSpacePoints()>0) {
301       HLTInfo("index grid slice %d partition %d", slice, patch);
302       fSpacePointGrid->Print();
303       for (AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid::iterator& cl=fSpacePointGrid->begin();
304            cl!=fSpacePointGrid->end(); cl++) {
305         AliHLTUInt32_t id=cl.Data().fId;
306         float row=fRawInputClusters->GetX(id);
307         float pad=fRawInputClusters->GetY(id);
308         float time=fRawInputClusters->GetZ(id);
309         HLTInfo("    cluster id 0x%08x: row %f  pad %f  time %f", id, row, pad, time);
310       }
311     }
312     if (GetBenchmarkInstance()) {
313       GetBenchmarkInstance()->Stop(1);
314       GetBenchmarkInstance()->Start(4);
315     }
316
317     // process the clusters per padrow and check the track grid
318     // for tracks crossing that particular padrow
319     if (GetBenchmarkInstance()) {
320       GetBenchmarkInstance()->Stop(4);
321       GetBenchmarkInstance()->Start(5);
322     }
323     allClusters+=fSpacePointGrid->GetNumberOfSpacePoints();
324     iResult=ProcessTrackClusters(&inputTrackArray[0], inputTrackArray.size(), fTrackGrid, fSpacePointGrid, fRawInputClusters, slice, patch);
325     int assignedInThisPartition=0;
326     if (iResult>=0) {
327       assignedInThisPartition=iResult;
328       associatedClusters+=iResult;
329     }
330     iResult=ProcessRemainingClusters(&inputTrackArray[0], inputTrackArray.size(), fTrackGrid, fSpacePointGrid, fRawInputClusters, slice, patch);
331     if (iResult>=0) {
332       if (fSpacePointGrid->GetNumberOfSpacePoints()>0) {
333         if (fVerbosity>0) HLTInfo("associated %d (%d) of %d clusters in slice %d partition %d", iResult+assignedInThisPartition, assignedInThisPartition, fSpacePointGrid->GetNumberOfSpacePoints(), slice, patch);
334       }
335       associatedClusters+=iResult;
336     }
337
338     // write all remaining clusters not yet assigned to tracks
339     // the index grid is used to write sorted in padrow
340     // FIXME: decoder index instead of data specification to be used
341     // use an external access grid to reduce allocated memory
342     // set to NULL after writing the clusters
343     const char* writeoptions="";
344     if (fpWrittenAssociatedClusterIds) {
345       writeoptions="write-cluster-ids";
346     }
347     fRawInputClusters->SetSpacePointPropertyGrid(pDesc->fSpecification, fSpacePointGrid);
348     iResult=fRawInputClusters->Write(outputPtr+size, capacity-size, outputBlocks, fpDataDeflater, writeoptions);
349     fRawInputClusters->SetSpacePointPropertyGrid(pDesc->fSpecification, NULL);
350     if (iResult>=0) {
351       size+=iResult;
352       outputDataSize+=iResult;
353       if (GetBenchmarkInstance()) GetBenchmarkInstance()->AddOutput(iResult);
354     }
355     if (GetBenchmarkInstance()) {
356       GetBenchmarkInstance()->Stop(5);
357     }
358
359     // forward MC labels
360     if (bHaveMC) {
361       // loop over input blocks and find MC block of current specification
362       unsigned mcBlock=0;
363       for (; mcBlock<evtData.fBlockCnt; mcBlock++) {
364         if (inputBlocks[mcBlock].fDataType!=(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo | kAliHLTDataOriginTPC) ||
365             inputBlocks[mcBlock].fSpecification!=pDesc->fSpecification) {
366           continue;
367         }
368         if (size+inputBlocks[mcBlock].fSize>capacity) {
369           iResult=-ENOSPC;
370           break;
371         }
372         iResult=ForwardMCLabels(inputBlocks[mcBlock], fSpacePointGrid, outputPtr+size, capacity-size, size, outputBlocks);
373         if (iResult>0) {
374           size+=iResult;
375         }
376
377         HLTDebug("forwarded MC data block of slice %d partition %d", slice, patch);
378         break;
379       }
380       if (mcBlock==evtData.fBlockCnt) {
381         HLTWarning("no mc data found for slice %d partition %d", slice, patch);
382       }
383     }
384
385     fSpacePointGrid->Clear();
386   }
387   if (fHistoClusterRatio && allClusters>0) {
388     if (fVerbosity>0) HLTInfo("associated %d of %d clusters to tracks", associatedClusters, allClusters);
389     float ratio=associatedClusters; ratio/=allClusters;
390     fHistoClusterRatio->Fill(ratio);
391   }
392
393   // output of track model clusters
394   if (iResult>=0) do {
395     AliHLTUInt32_t tracksBufferOffset=sizeof(AliHLTTPCTrackModelBlock);
396     if (capacity-size<tracksBufferOffset) {
397       iResult=-ENOSPC;
398       break;
399     }
400     if (fpWrittenAssociatedClusterIds) fpWrittenAssociatedClusterIds->clear();
401     AliHLTTPCTrackModelBlock* trackModelBlock=reinterpret_cast<AliHLTTPCTrackModelBlock*>(outputPtr+size);
402     trackModelBlock->fVersion=1;
403     trackModelBlock->fDeflaterMode=fpDataDeflater?fpDataDeflater->GetDeflaterVersion():0;
404     trackModelBlock->fTrackCount=inputTrackArray.size();
405     trackModelBlock->fClusterCount=0;
406     trackModelBlock->fGlobalParameterCnt=5;
407     tracksBufferOffset+=trackModelBlock->fGlobalParameterCnt*sizeof(trackModelBlock->fGlobalParameters);
408     if (capacity-size<tracksBufferOffset) {
409       iResult=-ENOSPC;
410       break;
411     }
412
413     AliHLTUInt32_t parameterIndex=0;
414     trackModelBlock->fGlobalParameters[parameterIndex++]=GetBz();
415     trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeFactorA;
416     trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeOffsetA;
417     trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeFactorC;
418     trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeOffsetC;
419     if (parameterIndex!=trackModelBlock->fGlobalParameterCnt) {
420       HLTError("internal error, size of parameter array has changed without providing all values");
421       iResult=-EFAULT;
422       break;
423     }
424
425     iResult=WriteTrackClusters(inputTrackArray, fRawInputClusters, fpDataDeflater, outputPtr+size+tracksBufferOffset, capacity-size-tracksBufferOffset);
426     if (iResult>=0) {
427       AliHLTComponent_BlockData bd;
428       FillBlockData(bd);
429       bd.fOffset        = size;
430       bd.fSize          = tracksBufferOffset+iResult;
431       bd.fDataType      = AliHLTTPCDefinitions::ClusterTracksCompressedDataType();
432       bd.fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(minSlice, maxSlice, minPatch, maxPatch);
433       outputBlocks.push_back(bd);
434       size += bd.fSize;
435       outputDataSize+=bd.fSize;
436       HLTBenchmark("track data block of %d tracks: size %d", inputTrackArray.size(), bd.fSize);
437
438       if (fpWrittenAssociatedClusterIds && fpWrittenAssociatedClusterIds->size()>0) {
439         AliHLTComponent::FillBlockData(bd);
440         bd.fOffset        = size;
441         bd.fSize        = fpWrittenAssociatedClusterIds->size()*sizeof(vector<AliHLTUInt32_t>::value_type);
442         if (capacity-size>bd.fSize) {
443           memcpy(outputPtr+bd.fOffset, &(*fpWrittenAssociatedClusterIds)[0], bd.fSize);
444           bd.fDataType    = AliHLTTPCDefinitions::ClusterIdTracksDataType();
445           bd.fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(minSlice, maxSlice, minPatch, maxPatch);
446           outputBlocks.push_back(bd);    
447           size += bd.fSize;
448         } else {
449           iResult=-ENOSPC;
450         }
451         
452         fpWrittenAssociatedClusterIds->clear();
453       }
454     }
455   } while (0);
456
457   fRawInputClusters->Clear();
458
459   float compressionFactor=(float)inputRawClusterSize;
460   if ((outputDataSize)>0) compressionFactor/=outputDataSize;
461   else compressionFactor=0.;
462   if (fHistoCompFactor) fHistoCompFactor->Fill(compressionFactor);
463
464   if (GetBenchmarkInstance()) {
465     GetBenchmarkInstance()->Stop(0);
466     if (fDeflaterMode!=3) {
467       HLTBenchmark("%s - compression factor %.2f", GetBenchmarkInstance()->GetStatistics(), compressionFactor);
468     } else {
469       HLTBenchmark("%s", GetBenchmarkInstance()->GetStatistics());
470     }
471   }
472
473   if (fInputClusters) {
474     fInputClusters->Clear();
475   }
476   if (fRawInputClusters) {
477     fRawInputClusters->Clear();
478   }
479   if (fTrackGrid) {
480     fTrackGrid->Clear();
481   }
482
483   return iResult;
484 }
485
486 int AliHLTTPCDataCompressionComponent::ProcessTrackClusters(AliHLTGlobalBarrelTrack* pTracks, unsigned nofTracks,
487                                                             AliHLTTrackGeometry::AliHLTTrackGrid* pTrackIndex,
488                                                             AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid* pClusterIndex,
489                                                             AliHLTSpacePointContainer* pClusters,
490                                                             int slice, int partition) const
491 {
492   // process to assigned track clusters
493   int assignedClusters=0;
494   if (!pTracks || nofTracks==0) return 0;
495
496   vector<AliHLTUInt32_t> processedTracks;
497   for (AliHLTTrackGeometry::AliHLTTrackGrid::iterator& trackId=pTrackIndex->begin(slice, partition, -1);
498        trackId!=pTrackIndex->end(); trackId++) {
499     if (find(processedTracks.begin(), processedTracks.end(), trackId.Data())!=processedTracks.end()) {
500       continue;
501     }
502     unsigned trackindex=0;
503     for (; trackindex<nofTracks; trackindex++) {
504       if ((unsigned)pTracks[trackindex].GetID()==trackId.Data()) break;
505     }
506     if (trackindex>=nofTracks) {
507       HLTError("can not find track of id %d", trackId.Data());
508       continue;
509     }
510     processedTracks.push_back(trackId.Data());
511     AliHLTGlobalBarrelTrack& track=pTracks[trackindex];
512     if (!track.GetTrackGeometry()) {
513       HLTError("can not find track geometry for track %d", trackId.Data());
514       continue;
515     }
516     AliHLTTPCTrackGeometry* pTrackPoints=dynamic_cast<AliHLTTPCTrackGeometry*>(track.GetTrackGeometry());
517     if (!pTrackPoints) {
518       HLTError("invalid track geometry type for track %d, expecting AliHLTTPCTrackGeometry", trackId.Data());
519       continue; 
520     }
521
522     UInt_t nofTrackPoints=track.GetNumberOfPoints();
523     const UInt_t* trackPoints=track.GetPoints();
524     for (unsigned i=0; i<nofTrackPoints; i++) {
525       const AliHLTUInt32_t& clusterId=trackPoints[i];
526       if (AliHLTTPCSpacePointData::GetSlice(clusterId)!=(unsigned)slice ||
527           AliHLTTPCSpacePointData::GetPatch(clusterId)!=(unsigned)partition) {
528         // not in the current partition;
529         continue;
530       }
531           
532       AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid::iterator& cl=pClusterIndex->find(AliHLTSpacePointContainer::AliHLTSpacePointProperties(clusterId));
533       if (cl==pClusterIndex->end()) {
534         HLTError("can not find cluster no 0x%08x of track %d in index grid", clusterId, track.GetID());
535         continue;
536       }
537
538       int clusterrow=(int)pClusters->GetX(clusterId);
539       float clusterpad=pClusters->GetY(clusterId);
540       float clustertime=pClusters->GetZ(clusterId);
541
542       AliHLTUInt32_t pointId=AliHLTTPCSpacePointData::GetID(slice, partition, clusterrow);
543       AliHLTTrackGeometry::AliHLTTrackPoint* point=pTrackPoints->GetRawTrackPoint(pointId);
544       if (!point) {
545         //HLTError("can not find track point slice %d partition %d padrow %d (0x%08x) of track %d", slice, partition, clusterrow, pointId, trackId.Data());
546         continue;
547       }
548       float pad=point->GetU();
549       float time=point->GetV();
550       if (TMath::Abs(clusterpad-pad)<fMaxDeltaPad &&
551           TMath::Abs(clustertime-time)<fMaxDeltaTime) {
552         // add this cluster to the track point and mark in the index grid
553         cl.Data().fTrackId=track.GetID();
554         point->AddAssociatedSpacePoint(clusterId, clusterpad-pad, clustertime-time);
555         assignedClusters++;
556       }
557     }
558   }
559   return assignedClusters;
560 }
561
562 int AliHLTTPCDataCompressionComponent::ProcessRemainingClusters(AliHLTGlobalBarrelTrack* pTracks, unsigned nofTracks,
563                                                                 AliHLTTrackGeometry::AliHLTTrackGrid* pTrackIndex,
564                                                                 AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid* pClusterIndex,
565                                                                 AliHLTSpacePointContainer* pClusters,
566                                                                 int slice, int partition) const
567 {
568   // assign remaining clusters to tracks
569   int iResult=0;
570   int associatedClusters=0;
571   if (!pTracks || nofTracks==0) return 0;
572
573   for (int padrow=0; padrow<AliHLTTPCTransform::GetNRows(partition); padrow++) {
574     for (AliHLTTrackGeometry::AliHLTTrackGrid::iterator& trackId=pTrackIndex->begin(slice, partition, padrow);
575          trackId!=pTrackIndex->end(); trackId++) {
576       unsigned i=0;
577       for (; i<nofTracks; i++) {
578         if ((unsigned)pTracks[i].GetID()==trackId.Data()) break;
579       }
580       if (i>=nofTracks) {
581         HLTError("can not find track of id %d", trackId.Data());
582         continue;
583       }
584       AliHLTGlobalBarrelTrack& track=pTracks[i];
585       if (!track.GetTrackGeometry()) {
586         HLTError("can not find track geometry for track %d", trackId.Data());
587         continue;
588       }
589       AliHLTTPCTrackGeometry* pTrackPoints=dynamic_cast<AliHLTTPCTrackGeometry*>(track.GetTrackGeometry());
590       if (!pTrackPoints) {
591         HLTError("invalid track geometry type for track %d, expecting AliHLTTPCTrackGeometry", trackId.Data());
592         continue;       
593       }
594       AliHLTUInt32_t pointId=AliHLTTPCSpacePointData::GetID(slice, partition, padrow);
595       AliHLTTrackGeometry::AliHLTTrackPoint* point=pTrackPoints->GetRawTrackPoint(pointId);
596       if (!point) {
597         //HLTError("can not find track point slice %d partition %d padrow %d (0x%08x) of track %d", slice, partition, padrow, pointId, trackId.Data());
598         continue;
599       }
600       float pad=point->GetU();
601       float time=point->GetV();
602
603       iResult=FindCellClusters(trackId.Data(), padrow, pad, time, pClusterIndex, pClusters, point);
604       if (iResult>0) associatedClusters+=iResult;
605       if (fVerbosity>0) {
606         HLTInfo("trackpoint track %d slice %d partition %d padrow %d: %.3f \t%.3f - associated %d", track.GetID(), slice, partition, padrow, pad, time, iResult);
607       }
608     }
609   }
610   if (iResult<0) return iResult;
611   return associatedClusters;
612 }
613
614 int AliHLTTPCDataCompressionComponent::FindCellClusters(int trackId, int padrow, float pad, float time,
615                                                         AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid* pClusterIndex,
616                                                         AliHLTSpacePointContainer* pClusters,
617                                                         AliHLTTrackGeometry::AliHLTTrackPoint* pTrackPoint) const
618 {
619   // check index cell for entries and assign to track
620   int count=0;
621   // search a 4x4 matrix out of the 9x9 matrix around the cell addressed by
622   // pad and time
623   int rowindex=pClusterIndex->GetYIndex((float)padrow);
624   int padindex=pClusterIndex->GetYIndex(pad);
625   int timeindex=pClusterIndex->GetZIndex(time);
626   if (pClusterIndex->GetCenterY(padindex)>pad) padindex--;
627   if (pClusterIndex->GetCenterZ(timeindex)>pad) timeindex--;
628   for (int padcount=0; padcount<2; padcount++, padindex++) {
629     if (padindex<0) continue;
630     if (padindex>=pClusterIndex->GetDimensionY()) break;
631     for (int timecount=0; timecount<2; timecount++, timeindex++) {
632       if (timeindex<0) continue;
633       if (timeindex>=pClusterIndex->GetDimensionZ()) break;
634       int cellindex=pClusterIndex->Index(rowindex, padindex, timeindex);
635       pad=pClusterIndex->GetCenterY(cellindex);
636       time=pClusterIndex->GetCenterZ(cellindex);
637       for (AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid::iterator& cl=pClusterIndex->begin((float)padrow, pad, time);
638            cl!=pClusterIndex->end(); cl++) {
639         if (cl.Data().fTrackId>=0) continue;
640         float clusterpad=pClusters->GetY(cl.Data().fId);
641         float clustertime=pClusters->GetZ(cl.Data().fId);
642         if (TMath::Abs(clusterpad-pad)<fMaxDeltaPad &&
643             TMath::Abs(clustertime-time)<fMaxDeltaTime) {
644           // add this cluster to the track point and mark in the index grid
645           cl.Data().fTrackId=trackId;
646           pTrackPoint->AddAssociatedSpacePoint(cl.Data().fId, clusterpad-pad, clustertime-time);
647           count++;
648         }
649       }
650     }
651   }
652   return count;
653 }
654
655 int AliHLTTPCDataCompressionComponent::WriteTrackClusters(const vector<AliHLTGlobalBarrelTrack>& tracks,
656                                                           AliHLTSpacePointContainer* pSpacePoints,
657                                                           AliHLTDataDeflater* pDeflater,
658                                                           AliHLTUInt8_t* outputPtr,
659                                                           AliHLTUInt32_t capacity) const
660 {
661   // write the track data block including all associated clusters
662   AliHLTUInt32_t size=0;
663   for (vector<AliHLTGlobalBarrelTrack>::const_iterator track=tracks.begin();
664        track!=tracks.end();
665        track++) {
666     if (!track->GetTrackGeometry()) {
667       HLTError("can not find track geometry for track %d", track->GetID());
668       return -EBADF;
669     }
670     AliHLTTPCTrackGeometry* pTrackPoints=dynamic_cast<AliHLTTPCTrackGeometry*>(track->GetTrackGeometry());
671     if (!pTrackPoints) {
672       HLTError("invalid track geometry type for track %d, expecting AliHLTTPCTrackGeometry", track->GetID());
673       return -EBADF;
674     }
675
676     int result=pTrackPoints->Write(*track, pSpacePoints, pDeflater, outputPtr+size, capacity-size, fpWrittenAssociatedClusterIds);
677     if (result<0) return result;
678     size+=result;
679
680     UInt_t nofTrackPoints=track->GetNumberOfPoints();
681     const UInt_t* trackPoints=track->GetPoints();
682
683     int assignedPoints=0;
684     int assignedTrackPoints=0;
685     const vector<AliHLTTrackGeometry::AliHLTTrackPoint>& rawPoints=pTrackPoints->GetRawPoints();
686     for (vector<AliHLTTrackGeometry::AliHLTTrackPoint>::const_iterator point=rawPoints.begin();
687          point!=rawPoints.end(); point++) {
688       const vector<AliHLTTrackGeometry::AliHLTTrackSpacepoint>& spacePoints=point->GetSpacepoints();
689       for (vector<AliHLTTrackGeometry::AliHLTTrackSpacepoint>::const_iterator spacePoint=spacePoints.begin();
690            spacePoint!=spacePoints.end(); spacePoint++) {
691         float dpad=spacePoint->GetResidual(0);
692         float dtime=spacePoint->GetResidual(1);
693         if (dpad>-1000 && dtime>-1000 && fHistoResidualPad && fHistoResidualTime) {
694           fHistoResidualPad->Fill(dpad);
695           fHistoResidualTime->Fill(dtime);
696         }
697         assignedPoints++;
698         for (unsigned i=0; i<nofTrackPoints; i++) {
699           if (trackPoints[i]==spacePoint->fId) {
700             assignedTrackPoints++;
701             break;
702           }
703         }
704       }
705     }
706     if (fHistoClustersOnTracks) {
707       fHistoClustersOnTracks->Fill(assignedPoints);
708     }
709     if (fHistoTrackClusterRatio && nofTrackPoints>0) {
710       float ratio=assignedTrackPoints; ratio/=nofTrackPoints;
711       fHistoTrackClusterRatio->Fill(ratio);
712     }
713   }
714   return size;
715 }
716
717 int AliHLTTPCDataCompressionComponent::DoInit( int argc, const char** argv )
718 {
719   /// inherited from AliHLTComponent: component initialisation and argument scan.
720   int iResult=0;
721
722   // component configuration
723   //Stage 1: default initialization.
724   //Default values.
725
726   //Stage 2: OCDB.
727   TString cdbPath("HLT/ConfigTPC/");
728   cdbPath += GetComponentID();
729   //
730   iResult = ConfigureFromCDBTObjString(cdbPath);
731   if (iResult < 0) 
732     return iResult;
733
734   //Stage 3: command line arguments.
735   if (argc && (iResult = ConfigureFromArgumentString(argc, argv)) < 0)
736     return iResult;
737
738   std::auto_ptr<AliHLTComponentBenchmark> benchmark(new AliHLTComponentBenchmark);
739   if (benchmark.get()) {
740     benchmark->SetTimer(0,"total");
741     benchmark->SetTimer(1,"rawclusterinput");
742     benchmark->SetTimer(2,"clusterinput");
743     benchmark->SetTimer(3,"trackinput");
744     benchmark->SetTimer(4,"processing");
745     benchmark->SetTimer(5,"output");
746   } else {
747     return -ENOMEM;
748   }
749
750   unsigned spacePointContainerMode=(fMode==2)?AliHLTTPCHWCFSpacePointContainer::kModeCreateMap:0;
751   std::auto_ptr<AliHLTTPCHWCFSpacePointContainer> rawInputClusters(new AliHLTTPCHWCFSpacePointContainer(spacePointContainerMode));
752   std::auto_ptr<AliHLTTPCSpacePointContainer> inputClusters(new AliHLTTPCSpacePointContainer);
753
754   std::auto_ptr<TH1F> histoCompFactor(new TH1F("CompressionFactor",
755                                                "HLT TPC data compression factor",
756                                                100, 0., 10.));
757   std::auto_ptr<TH1F> histoResidualPad(new TH1F("PadResidual",
758                                                 "HLT TPC pad residual",
759                                                 100, -fMaxDeltaPad, fMaxDeltaPad));
760   std::auto_ptr<TH1F> histoResidualTime(new TH1F("TimeResidual",
761                                                  "HLT TPC time residual",
762                                                  100, -fMaxDeltaTime, fMaxDeltaTime));
763   std::auto_ptr<TH1F> histoClustersOnTracks(new TH1F("ClustersOnTracks",
764                                                  "Clusters in track model compression",
765                                                  200, 0., 600));
766   std::auto_ptr<TH1F> histoClusterRatio(new TH1F("ClusterRatio",
767                                                  "Fraction of clusters in track model compression",
768                                                  100, 0., 1.));
769   std::auto_ptr<TH1F> histoTrackClusterRatio(new TH1F("UsedTrackClusters",
770                                                  "Fraction of track clusters in track model compression",
771                                                  100, 0., 1.));
772
773   // track grid: 36 slices, each 6 partitions with max 33 rows
774   fTrackGrid=new AliHLTTrackGeometry::AliHLTTrackGrid(36, 1, 6, 1, 33, 1, 20000);
775   fSpacePointGrid=AliHLTTPCHWCFSpacePointContainer::AllocateIndexGrid();
776
777   if (!rawInputClusters.get() ||
778       !inputClusters.get() ||
779       !fTrackGrid ||
780       !fSpacePointGrid) {
781     if (fTrackGrid) delete fTrackGrid; fTrackGrid=NULL;
782     if (fSpacePointGrid) delete fSpacePointGrid; fSpacePointGrid=NULL;
783     return -ENOMEM;
784   }
785
786   if (fDeflaterMode>0 && (iResult=InitDeflater(fDeflaterMode))<0)
787     return iResult;
788
789   fpBenchmark=benchmark.release();
790   fRawInputClusters=rawInputClusters.release();
791   fInputClusters=inputClusters.release();
792
793   // initialize the histograms if stored at the end
794   // condition might be extended
795   if (!fHistogramFile.IsNull()) {
796     fHistoCompFactor=histoCompFactor.release();
797     fHistoResidualPad=histoResidualPad.release();
798     fHistoResidualTime=histoResidualTime.release();
799     fHistoClustersOnTracks=histoClustersOnTracks.release();
800     fHistoClusterRatio=histoClusterRatio.release();
801     fHistoTrackClusterRatio=histoTrackClusterRatio.release();
802   }
803
804   if (iResult>=0 && (iResult=InitDriftTimeTransformation())<0) return iResult;
805
806   return iResult;
807 }
808
809 int AliHLTTPCDataCompressionComponent::InitDeflater(int mode)
810 {
811   /// init the data deflater
812   int iResult=0;
813   if (mode==2 || mode==3) {
814     // huffman deflater
815     std::auto_ptr<AliHLTDataDeflaterHuffman> deflater(new AliHLTDataDeflaterHuffman(mode==3));
816     if (!deflater.get()) return -ENOMEM;
817
818     if (!deflater->IsTrainingMode()) {
819       TString cdbPath("HLT/ConfigTPC/");
820       cdbPath += GetComponentID();
821       cdbPath += "HuffmanTables";
822       TObject* pConf=LoadAndExtractOCDBObject(cdbPath);
823       if (!pConf) return -ENOENT;
824       if (dynamic_cast<TList*>(pConf)==NULL) {
825         HLTError("huffman table configuration object of inconsistent type");
826         return -EINVAL;
827       }
828       iResult=deflater->InitDecoders(dynamic_cast<TList*>(pConf));
829       if (iResult<0) return iResult;
830     }
831     
832     unsigned nofParameters=AliHLTTPCDefinitions::GetNumberOfClusterParameterDefinitions();
833     unsigned p=0;
834     for (; p<nofParameters; p++) {
835       const AliHLTTPCDefinitions::AliClusterParameter& parameter=AliHLTTPCDefinitions::fgkClusterParameterDefinitions[p];
836       if (deflater->AddParameterDefinition(parameter.fName,
837                                            parameter.fBitLength)!=(int)parameter.fId) {
838         // for performance reason the parameter id is simply used as index in the array of
839         // definitions, the position must match the id
840         HLTFatal("mismatch between parameter id and position in array for parameter %s, rearrange definitions!", parameter.fName);
841         return -EFAULT;
842       }
843     }
844     fpDataDeflater=deflater.release();
845     return 0;
846   }
847   if (mode==1) {
848     std::auto_ptr<AliHLTDataDeflaterSimple> deflater(new AliHLTDataDeflaterSimple);
849     if (!deflater.get()) return -ENOMEM;
850
851     unsigned nofParameters=AliHLTTPCDefinitions::GetNumberOfClusterParameterDefinitions();
852     unsigned p=0;
853     for (; p<nofParameters; p++) {
854       const AliHLTTPCDefinitions::AliClusterParameter& parameter=AliHLTTPCDefinitions::fgkClusterParameterDefinitions[p];
855       if (deflater->AddParameterDefinition(parameter.fName,
856                                            parameter.fBitLength,
857                                            parameter.fOptional)!=(int)parameter.fId) {
858         // for performance reason the parameter id is simply used as index in the array of
859         // definitions, the position must match the id
860         HLTFatal("mismatch between parameter id and position in array for parameter %s, rearrange definitions!", parameter.fName);
861         return -EFAULT;
862       }
863     }
864     fpDataDeflater=deflater.release();
865     return 0;
866   }
867   HLTError("invalid deflater mode %d, allowed 1=simple 2=huffman", mode);
868   return -EINVAL;
869 }
870
871 int AliHLTTPCDataCompressionComponent::DoDeinit()
872 {
873   /// inherited from AliHLTComponent: component cleanup
874   int iResult=0;
875   if (fpBenchmark) delete fpBenchmark; fpBenchmark=NULL;
876   if (fRawInputClusters) delete fRawInputClusters; fRawInputClusters=NULL;
877   if (fInputClusters) delete fInputClusters; fInputClusters=NULL;
878   if (!fHistogramFile.IsNull()) {
879     TFile out(fHistogramFile, "RECREATE");
880     if (!out.IsZombie()) {
881       out.cd();
882       if (fHistoCompFactor) fHistoCompFactor->Write();
883       if (fHistoResidualPad) fHistoResidualPad->Write();
884       if (fHistoResidualTime) fHistoResidualTime->Write();
885       if (fHistoClusterRatio) fHistoClusterRatio->Write();
886       if (fHistoClustersOnTracks) fHistoClustersOnTracks->Write();
887       if (fHistoTrackClusterRatio) fHistoTrackClusterRatio->Write();
888       out.Close();
889     }
890   }
891   if (fHistoCompFactor) delete fHistoCompFactor;
892   fHistoCompFactor=NULL;
893   if (fHistoResidualPad) delete fHistoResidualPad;
894   fHistoResidualPad=NULL;
895   if (fHistoResidualTime) delete fHistoResidualTime;
896   fHistoResidualTime=NULL;
897   if (fHistoClustersOnTracks) delete fHistoClustersOnTracks;
898   fHistoClustersOnTracks=NULL;
899   if (fHistoClusterRatio) delete fHistoClusterRatio;
900   fHistoClusterRatio=NULL;
901   if (fHistoTrackClusterRatio) delete fHistoTrackClusterRatio;
902   fHistoTrackClusterRatio=NULL;
903
904   if (fpDataDeflater) {
905     if (!fHistogramFile.IsNull()) {
906       TString filename=fHistogramFile;
907       filename.ReplaceAll(".root", "-deflater.root");
908       fpDataDeflater->SaveAs(filename);
909     }
910     if (fDeflaterMode==3) {
911       if (fTrainingTableOutput.IsNull()) {
912         fTrainingTableOutput=GetComponentID();
913         fTrainingTableOutput+="-huffman.root";
914       }
915       // TODO: currently, the code tables are also calculated in FindObject
916       // check if a different function is more appropriate
917       TObject* pConf=fpDataDeflater->FindObject("DeflaterConfiguration");
918       if (pConf) {
919         TString cdbEntryPath("HLT/ConfigTPC/");
920         cdbEntryPath += GetComponentID();
921         cdbEntryPath += "HuffmanTables";
922         AliCDBPath cdbPath(cdbEntryPath);
923         AliCDBId cdbId(cdbPath, AliCDBManager::Instance()->GetRun(), AliCDBRunRange::Infinity(), 0, 0);
924         AliCDBMetaData* cdbMetaData=new AliCDBMetaData;
925         cdbMetaData->SetResponsible("ALICE HLT Matthias.Richter@cern.ch");
926         cdbMetaData->SetComment("Huffman encoder configuration");
927         AliCDBEntry* entry=new AliCDBEntry(pConf, cdbId, cdbMetaData, kTRUE);
928
929         entry->SaveAs(fTrainingTableOutput);
930         // this is a small memory leak
931         // seg fault in ROOT object handling if the two objects are deleted
932         // investigate later
933         //delete entry;
934         //delete cdbMetaData;
935       }
936     }
937     delete fpDataDeflater;
938   }
939   fpDataDeflater=NULL;
940
941
942   if (fTrackGrid) delete fTrackGrid; fTrackGrid=NULL;
943   if (fSpacePointGrid) delete fSpacePointGrid; fSpacePointGrid=NULL;
944
945   return iResult;
946 }
947
948 int AliHLTTPCDataCompressionComponent::ScanConfigurationArgument(int argc, const char** argv)
949 {
950   /// inherited from AliHLTComponent: argument scan
951   int iResult=0;
952   if (argc<1) return 0;
953   int bMissingParam=0;
954   int i=0;
955   TString argument=argv[i];
956
957   do {
958     // -mode
959     if (argument.CompareTo("-mode")==0) {
960       if ((bMissingParam=(++i>=argc))) break;
961       TString parameter=argv[i];
962       if (parameter.IsDigit()) {
963         fMode=parameter.Atoi();
964         return 2;
965       } else {
966         HLTError("invalid parameter for argument %s, expecting number instead of %s", argument.Data(), parameter.Data());
967         return -EPROTO;
968       }
969     }
970
971     // -deflater-mode
972     if (argument.CompareTo("-deflater-mode")==0) {
973       if ((bMissingParam=(++i>=argc))) break;
974       TString parameter=argv[i];
975       if (parameter.IsDigit()) {
976         fDeflaterMode=parameter.Atoi();
977         return 2;
978       } else {
979         HLTError("invalid parameter for argument %s, expecting number instead of %s", argument.Data(), parameter.Data());
980         return -EPROTO;
981       }
982     }
983
984     // -histogram-file
985     if (argument.CompareTo("-histogram-file")==0) {
986       if ((bMissingParam=(++i>=argc))) break;
987       fHistogramFile=argv[i++];
988       return 2;
989     }
990     // -save-histogram-table
991     if (argument.CompareTo("-save-huffman-table")==0) {
992       if ((bMissingParam=(++i>=argc))) break;
993       fTrainingTableOutput=argv[i++];
994       return 2;
995     }
996     // -cluster-verification
997     if (argument.CompareTo("-cluster-verification")==0) {
998       if ((bMissingParam=(++i>=argc))) break;
999       TString parameter=argv[i];
1000       if (parameter.IsDigit()) {
1001         fVerificationMode=parameter.Atoi();
1002         return 2;
1003       } else {
1004         HLTError("invalid parameter for argument %s, expecting number instead of %s", argument.Data(), parameter.Data());
1005         return -EPROTO;
1006       }
1007     }
1008   } while (0); // using do-while only to have break available
1009
1010   if (bMissingParam) {
1011     HLTError("missing parameter for argument %s", argument.Data());
1012     iResult=-EPROTO;
1013   }
1014
1015   return iResult;
1016 }
1017
1018 int AliHLTTPCDataCompressionComponent::ForwardMCLabels(const AliHLTComponentBlockData& desc,
1019                                                        AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid* pIndex,
1020                                                        AliHLTUInt8_t* outputPtr,
1021                                                        AliHLTUInt32_t size, AliHLTUInt32_t offset,
1022                                                        vector<AliHLTComponentBlockData>& outputBlocks) const
1023 {
1024   // forward the mc labels in the same order as the clusters
1025   // sorted according to the index grid
1026   if (!outputPtr || !pIndex) return -EINVAL;
1027   if (!desc.fPtr) return -ENODATA;
1028   if (size<desc.fSize) return -ENOSPC;
1029
1030   int slice=AliHLTTPCDefinitions::GetMinSliceNr(desc.fSpecification);
1031   int part=AliHLTTPCDefinitions::GetMinPatchNr(desc.fSpecification);
1032
1033   const AliHLTTPCClusterMCData* pInput = reinterpret_cast<const AliHLTTPCClusterMCData*>(desc.fPtr);
1034   Int_t nLabels = (Int_t) pInput->fCount;
1035   if (nLabels*sizeof(AliHLTTPCClusterMCLabel) + sizeof(AliHLTTPCClusterMCData) != desc.fSize) {
1036     HLTError("inconsistent cluster mc data block size, skipping block");
1037     return -EBADF;
1038   }
1039   const AliHLTTPCClusterMCLabel *pInputLabels = pInput->fLabels;
1040
1041   AliHLTTPCClusterMCData* pOutput = reinterpret_cast<AliHLTTPCClusterMCData*>(outputPtr);
1042   AliHLTTPCClusterMCLabel *pOutputLabels = pOutput->fLabels;
1043
1044   unsigned outIndex=0;
1045   for (AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid::iterator clusterID=pIndex->begin();
1046        clusterID!=pIndex->end(); clusterID++, outIndex++) {
1047       if (clusterID.Data().fTrackId>-1) {
1048         // this is an assigned cluster, skip
1049         // TODO: introduce selectors into AliHLTIndexGrid::begin to loop
1050         // consistently over entries, e.g. this check has to be done also
1051         // in the forwarding of MC labels in
1052         // AliHLTTPCHWCFSpacePointContainer::WriteSorted
1053         continue;
1054       }
1055       if ((unsigned)slice!=AliHLTTPCSpacePointData::GetSlice(clusterID.Data().fId) ||
1056           (unsigned)part!=AliHLTTPCSpacePointData::GetPatch(clusterID.Data().fId)) {
1057         HLTError("spacepoint index 0x%08x out of slice %d partition %d", clusterID.Data().fId, slice, part);
1058       }
1059       int index=AliHLTTPCSpacePointData::GetNumber(clusterID.Data().fId);
1060       pOutputLabels[outIndex]=pInputLabels[index];
1061   }
1062   if (outIndex==pInput->fCount) {
1063     pOutput->fCount=outIndex;
1064   } else {
1065     HLTError("failed to copy MC label data block 0x%08x: expecting %d, copied %d entries", desc.fSpecification, pInput->fCount, outIndex);
1066     return -EBADMSG;
1067   }
1068
1069   AliHLTComponent_BlockData bd;
1070   AliHLTComponent::FillBlockData(bd);
1071   bd.fSize          = desc.fSize;
1072   bd.fOffset        = offset;
1073   bd.fSpecification = desc.fSpecification;
1074   bd.fDataType      = desc.fDataType;
1075   outputBlocks.push_back(bd);
1076
1077   return bd.fSize;
1078 }
1079
1080 int AliHLTTPCDataCompressionComponent::InitDriftTimeTransformation()
1081 {
1082   /// calculate correction factor and offset for a linear approximation of the
1083   /// drift time transformation, separately for A and C side
1084   int iResult=0;
1085   AliHLTTPCClusterTransformation transform;
1086   if ((iResult=transform.Init( GetBz(), GetTimeStamp()))<0) {
1087     HLTError("failed to init AliHLTTPCClusterTransformation: %d", iResult);
1088     return iResult;
1089   }
1090
1091   if ((iResult=CalculateDriftTimeTransformation(transform, 0, 0, fDriftTimeFactorA, fDriftTimeOffsetA))<0) return iResult;
1092   if (fVerbosity>0) HLTInfo("drift time transformation A side: m=%f n=%f", fDriftTimeFactorA, fDriftTimeOffsetA);
1093   if ((iResult=CalculateDriftTimeTransformation(transform, 18, 0, fDriftTimeFactorC, fDriftTimeOffsetC))<0) return iResult;
1094   if (fVerbosity>0) HLTInfo("drift time transformation C side: m=%f n=%f", fDriftTimeFactorC, fDriftTimeOffsetC);
1095
1096   return 0;
1097 }
1098
1099 int AliHLTTPCDataCompressionComponent::CalculateDriftTimeTransformation(AliHLTTPCClusterTransformation& transform,
1100                                                                         int slice, int padrow,
1101                                                                         float& m, float& n) const
1102 {
1103   /// calculate correction factor and offset for a linear approximation of the
1104   /// drift time transformation by just probing the range of timebins with
1105   /// AliHLTTPCClusterTransformation
1106   const int nofSteps=100;
1107   vector<float> zvalues;
1108
1109   int nofTimebins=AliHLTTPCTransform::GetNTimeBins();
1110   int stepWidth=nofTimebins/nofSteps;
1111   int time=0;
1112   int count=0;
1113   float meanT=0.;
1114   float meanZ=0.;
1115   for (time=0; time<nofTimebins; time+=stepWidth, count++) {
1116     Float_t xyz[3];
1117     transform.Transform(slice, padrow, 0, time, xyz);
1118     zvalues.push_back(xyz[2]);
1119     meanT+=time;
1120     meanZ+=xyz[2];
1121   }
1122   meanT/=count;
1123   meanZ/=count;
1124   float sumTZ=.0;
1125   float sumT2=.0;
1126   time=0;
1127   for (vector<float>::const_iterator z=zvalues.begin();
1128        z!=zvalues.end(); z++, time+=stepWidth) {
1129     sumTZ+=(time-meanT)*((*z)-meanZ);
1130     sumT2+=(time-meanT)*(time-meanT);
1131   }
1132   m=sumTZ/sumT2;
1133   n=meanZ-m*meanT;
1134
1135   return 0;
1136 }