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