]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/comp/AliHLTTPCDataCompressionMonitorComponent.cxx
adding publishing of histograms, by default all histograms are sent as separate objec...
[u/mrichter/AliRoot.git] / HLT / TPCLib / comp / AliHLTTPCDataCompressionMonitorComponent.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   AliHLTTPCDataCompressionMonitorComponent.cxx
19 /// @author Matthias Richter
20 /// @date   2011-09-12
21 /// @brief  TPC component for monitoring of data compression
22 ///
23
24 #include "AliHLTTPCDataCompressionMonitorComponent.h"
25 #include "AliHLTTPCDataCompressionComponent.h"
26 #include "AliHLTTPCDataCompressionDecoder.h"
27 #include "AliHLTTPCDefinitions.h"
28 #include "AliHLTTPCHWCFData.h"
29 #include "AliHLTTPCDefinitions.h"
30 #include "AliHLTTPCClusterDataFormat.h"
31 #include "AliHLTTPCRawCluster.h"
32 #include "AliHLTTPCTransform.h"
33 #include "AliHLTTPCTrackGeometry.h"
34 #include "AliHLTTPCHWCFSpacePointContainer.h"
35 #include "AliHLTErrorGuard.h"
36 #include "AliRawDataHeader.h"
37 #include "AliTPCclusterMI.h"
38 #include "TH1I.h"
39 #include "TH1F.h"
40 #include "TH2I.h"
41 #include "TH2F.h"
42 #include "TH3I.h"
43 #include "TH3F.h"
44 #include "TFile.h"
45 #include "TObjArray.h"
46 #include "TList.h"
47 #include <memory>
48
49 ClassImp(AliHLTTPCDataCompressionMonitorComponent)
50
51 AliHLTTPCDataCompressionMonitorComponent::AliHLTTPCDataCompressionMonitorComponent()
52   : AliHLTProcessor()
53   , fpHWClusterDecoder(NULL)
54   , fHistoHWCFDataSize(NULL)
55   , fHistoHWCFReductionFactor(NULL)
56   , fHistoNofClusters(NULL)
57   , fHistogramFile()
58   , fMonitoringContainer(NULL)
59   , fVerbosity(0)
60   , fFlags(0)
61   , fPublishingMode(kPublishSeparate)
62 {
63 }
64
65 AliHLTTPCDataCompressionMonitorComponent::~AliHLTTPCDataCompressionMonitorComponent()
66 {
67   /// destructor
68 }
69
70
71 const char* AliHLTTPCDataCompressionMonitorComponent::GetComponentID()
72 {
73   /// inherited from AliHLTComponent: id of the component
74   return "TPCDataCompressorMonitor";
75 }
76
77
78 void AliHLTTPCDataCompressionMonitorComponent::GetInputDataTypes( AliHLTComponentDataTypeList& tgtList)
79 {
80   /// inherited from AliHLTComponent: list of data types in the vector reference
81   tgtList.clear();
82   tgtList.push_back(AliHLTTPCDefinitions::fgkHWClustersDataType);
83   tgtList.push_back(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
84   tgtList.push_back(AliHLTTPCDefinitions::fgkRawClustersDataType);
85   tgtList.push_back(AliHLTTPCDefinitions::RemainingClustersCompressedDataType());
86   tgtList.push_back(AliHLTTPCDefinitions::ClusterTracksCompressedDataType());  
87 }
88
89 AliHLTComponentDataType AliHLTTPCDataCompressionMonitorComponent::GetOutputDataType()
90 {
91   /// inherited from AliHLTComponent: output data type of the component.
92   return kAliHLTMultipleDataType;
93 }
94
95 int AliHLTTPCDataCompressionMonitorComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
96 {
97   /// inherited from AliHLTComponent: multiple output data types of the component.
98   tgtList.clear();
99   tgtList.push_back(kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
100   return tgtList.size();
101 }
102
103 void AliHLTTPCDataCompressionMonitorComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
104 {
105   /// inherited from AliHLTComponent: output data size estimator
106   constBase=100000;
107   inputMultiplier=1.0;
108 }
109
110 AliHLTComponent* AliHLTTPCDataCompressionMonitorComponent::Spawn()
111 {
112   /// inherited from AliHLTComponent: spawn function.
113   return new AliHLTTPCDataCompressionMonitorComponent;
114 }
115
116 int AliHLTTPCDataCompressionMonitorComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/, 
117                                                        const AliHLTComponentBlockData* /*inputBlocks*/, 
118                                                        AliHLTComponentTriggerData& /*trigData*/,
119                                                        AliHLTUInt8_t* /*outputPtr*/,
120                                                        AliHLTUInt32_t& /*size*/,
121                                                        AliHLTComponentBlockDataList& /*outputBlocks*/ )
122 {
123   /// inherited from AliHLTProcessor: data processing
124   int iResult=0;
125
126   if (!IsDataEvent()) return 0;
127
128   const AliHLTComponentBlockData* pDesc=NULL;
129   unsigned rawDataSize=0;
130   unsigned hwclustersDataSize=0;
131   unsigned nofClusters=0;
132   unsigned compDataSize=0; 
133   
134   // check size of TPC raw data
135   for (pDesc=GetFirstInputBlock(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC);
136        pDesc!=NULL; pDesc=GetNextInputBlock()) {
137     fFlags|=kHaveRawData;
138     rawDataSize+=pDesc->fSize;
139   }
140
141   // check size of HWCF data and add to the MonitoringContainer
142   for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkHWClustersDataType);
143        pDesc!=NULL; pDesc=GetNextInputBlock()) {
144     fFlags|=kHaveHWClusters;
145     // FIXME: the decoding can now be handled via the data container
146     if (pDesc->fSize<=sizeof(AliRawDataHeader)) continue;
147     if (fpHWClusterDecoder) {
148       hwclustersDataSize+=pDesc->fSize;
149       AliHLTUInt8_t* pData=reinterpret_cast<AliHLTUInt8_t*>(pDesc->fPtr);
150       pData+=sizeof(AliRawDataHeader);
151       if (fpHWClusterDecoder->Init(pData, pDesc->fSize-sizeof(AliRawDataHeader))<0 ||
152           (fpHWClusterDecoder->CheckVersion()<0 && (int)pDesc->fSize>fpHWClusterDecoder->GetRCUTrailerSize())) {
153         HLTError("data block of type %s corrupted: can not decode format",
154                  AliHLTComponent::DataType2Text(pDesc->fDataType).c_str());
155       } else {
156         nofClusters+=fpHWClusterDecoder->GetNumberOfClusters();
157       }
158     }
159     if (fMonitoringContainer) {
160       fMonitoringContainer->AddRawData(pDesc);
161     }
162   }
163
164   if (fMonitoringContainer) {
165     for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::RemainingClusterIdsDataType());
166          pDesc!=NULL; pDesc=GetNextInputBlock()) {
167       iResult=fMonitoringContainer->AddClusterIds(pDesc);
168     }
169
170     for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::ClusterIdTracksDataType());
171          pDesc!=NULL; pDesc=GetNextInputBlock()) {
172       iResult=fMonitoringContainer->AddClusterIds(pDesc);
173     }
174
175     // read data
176     AliHLTTPCDataCompressionDecoder decoder;
177     bool bHaveRawClusters=false;
178     for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::RawClustersDataType());
179          pDesc!=NULL; pDesc=GetNextInputBlock()) {
180       // Note: until r51411 and v5-01-Rev-03 the compressed cluster format was sent with data
181       // type {CLUSTRAW,TPC }, the version member indicated the actual type of data
182       // transparently handled in the decoder
183       bHaveRawClusters=true;
184       iResult=decoder.ReadClustersPartition(fMonitoringContainer->BeginRemainingClusterBlock(0, pDesc->fSpecification),
185                                             reinterpret_cast<AliHLTUInt8_t*>(pDesc->fPtr),
186                                             pDesc->fSize,
187                                             pDesc->fSpecification);
188       if (iResult<0) {
189         HLTError("reading of partition clusters failed with error %d", iResult);
190       }
191     }
192
193     if (!bHaveRawClusters) {
194     for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::RemainingClustersCompressedDataType());
195          pDesc!=NULL; pDesc=GetNextInputBlock()) {
196       iResult=decoder.ReadClustersPartition(fMonitoringContainer->BeginRemainingClusterBlock(0, pDesc->fSpecification),
197                                             reinterpret_cast<AliHLTUInt8_t*>(pDesc->fPtr),
198                                             pDesc->fSize,
199                                             pDesc->fSpecification);
200       compDataSize+=pDesc->fSize;
201     }
202
203     for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::ClusterTracksCompressedDataType());
204          pDesc!=NULL; pDesc=GetNextInputBlock()) {
205       iResult=decoder.ReadTrackModelClustersCompressed(fMonitoringContainer->BeginTrackModelClusterBlock(0),
206                                                reinterpret_cast<AliHLTUInt8_t*>(pDesc->fPtr),
207                                                pDesc->fSize,
208                                                pDesc->fSpecification);
209       compDataSize+=pDesc->fSize;
210     }
211     } else {
212       if (GetFirstInputBlock(AliHLTTPCDefinitions::RemainingClustersCompressedDataType()) ||
213           GetFirstInputBlock(AliHLTTPCDefinitions::ClusterTracksCompressedDataType())) {
214         ALIHLTERRORGUARD(5, "conflicting data blocks, monitoring histograms already filled from raw cluster data, ignoring blocks of compressed partition and track clusters");
215       }              
216     }
217
218     fMonitoringContainer->Clear();
219   }
220
221   float ratio=0;
222   if (compDataSize) {ratio=(float)hwclustersDataSize; ratio/=compDataSize;}
223   if (fHistoHWCFDataSize)        fHistoHWCFDataSize       ->Fill(hwclustersDataSize/1024, compDataSize/1024);
224   if (fHistoHWCFReductionFactor) fHistoHWCFReductionFactor->Fill(hwclustersDataSize/1024, ratio);
225   if (fHistoNofClusters)         fHistoNofClusters        ->Fill(hwclustersDataSize/1024, nofClusters);
226   HLTInfo("raw data %d, hwcf data %d, comp data %d, ratio %f, %d clusters\n", rawDataSize, hwclustersDataSize, compDataSize, ratio, nofClusters);
227
228   if (iResult>=0 && fPublishingMode!=kPublishOff) {
229     iResult=Publish(fPublishingMode);
230   }
231
232   return iResult;
233 }
234
235 int AliHLTTPCDataCompressionMonitorComponent::Publish(int mode)
236 {
237   /// publish to output
238
239   // FIXME: code needs to be optimized, maybe a bit to much new and delete for the
240   // moment, the data type might need adjustment
241   int iResult=0;
242   TObjArray* pArray=mode==kPublishArray?(new TObjArray):NULL;
243   TList* pList=mode==kPublishList?(new TList):NULL;
244   if (mode==kPublishSeparate) {
245     if (fHistoHWCFDataSize)        PushBack(fHistoHWCFDataSize       , kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
246     if (fHistoHWCFReductionFactor) PushBack(fHistoHWCFReductionFactor, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
247     if (fHistoNofClusters)         PushBack(fHistoNofClusters        , kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
248   } else if (pList) {
249     if (fHistoHWCFDataSize)        pList->Add(fHistoHWCFDataSize->Clone());
250     if (fHistoHWCFReductionFactor) pList->Add(fHistoHWCFReductionFactor->Clone());
251     if (fHistoNofClusters)         pList->Add(fHistoNofClusters->Clone());
252   } else if (pArray) {
253     if (fHistoHWCFDataSize)        pArray->Add(fHistoHWCFDataSize->Clone());
254     if (fHistoHWCFReductionFactor) pArray->Add(fHistoHWCFReductionFactor->Clone());
255     if (fHistoNofClusters)         pArray->Add(fHistoNofClusters->Clone());
256   }
257
258   if (fMonitoringContainer) {
259     static const char* searchIds[] = {"fHistograms", "fHistograms2D", "fHistograms3D", NULL};
260     const char** searchId=searchIds;
261     while (*searchId && iResult>=0) {
262       const TObject* o=fMonitoringContainer->FindObject(*searchId);
263       if (o) {
264         const TObjArray* histograms=dynamic_cast<const TObjArray*>(o);
265         if (histograms) {
266           for (int i=0; i<histograms->GetEntriesFast() && iResult>=0; i++) {
267             if (!histograms->At(i)) continue;
268             if (mode==kPublishSeparate) {
269               iResult=PushBack(histograms->At(i), kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
270             } else if (pList) {
271               pList->Add(histograms->At(i)->Clone());
272             } else if (pArray) {
273               pArray->Add(histograms->At(i)->Clone());
274             }
275           }
276         }
277       } else {
278         HLTError("failed to find object \"%s\"", *searchId);
279       }
280       searchId++;
281     }
282   }
283
284   if (pArray) {
285     iResult=PushBack(pArray, kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC);
286     pArray->SetOwner(kTRUE);
287     delete pArray;
288     pArray=NULL;
289   }
290   if (pList) {
291     iResult=PushBack(pList, kAliHLTDataTypeTObject|kAliHLTDataOriginTPC);
292     pList->SetOwner(kTRUE);
293     delete pList;
294     pList=NULL;
295   }
296   return iResult;
297 }
298
299 int AliHLTTPCDataCompressionMonitorComponent::DoInit( int argc, const char** argv )
300 {
301   /// inherited from AliHLTComponent: component initialisation and argument scan.
302   int iResult=0;
303
304   // component configuration
305   //Stage 1: default initialization.
306   //Default values.
307   fFlags=0;
308
309   //Stage 2: OCDB.
310   TString cdbPath("HLT/ConfigTPC/");
311   cdbPath += GetComponentID();
312   //
313   // iResult = ConfigureFromCDBTObjString(cdbPath);
314   // if (iResult < 0) 
315   //   return iResult;
316
317   //Stage 3: command line arguments.
318   if (argc && (iResult = ConfigureFromArgumentString(argc, argv)) < 0)
319     return iResult;
320
321   std::auto_ptr<AliHLTTPCHWCFData> hwClusterDecoder(new AliHLTTPCHWCFData);
322   std::auto_ptr<AliDataContainer> dataContainer(new AliDataContainer);
323
324   std::auto_ptr<TH2I> histoHWCFDataSize(new TH2I("HWCFDataSize",
325                                                  "HW ClusterFinder Size",
326                                                  100, 0., 80000., 100, 0., 80000.));
327   if (histoHWCFDataSize.get()) {
328     TAxis* xaxis=histoHWCFDataSize->GetXaxis();
329     if (xaxis) xaxis->SetTitle("hwcf size [kB]");
330     TAxis* yaxis=histoHWCFDataSize->GetYaxis();
331     if (yaxis) yaxis->SetTitle("compressed data size [kb]");
332   }
333
334   std::auto_ptr<TH2I> histoHWCFReductionFactor(new TH2I("HWCFReductionFactor",
335                                                         "Data reduction HW ClusterFinder",
336                                                         100, 0., 80000., 100, 0., 10.));
337   if (histoHWCFReductionFactor.get()) {
338     TAxis* xaxis=histoHWCFReductionFactor->GetXaxis();
339     if (xaxis) xaxis->SetTitle("hwcf size [kB]");
340     TAxis* yaxis=histoHWCFReductionFactor->GetYaxis();
341     if (yaxis) yaxis->SetTitle("reduction factor");
342   }
343
344   std::auto_ptr<TH2I> histoNofClusters(new TH2I("NofClusters",
345                                                "Number of HLT TPC clusters",
346                                                100, 0., 80000., 100, 0., 3000000.));
347   if (histoNofClusters.get()) {
348     TAxis* xaxis=histoNofClusters->GetXaxis();
349     if (xaxis) xaxis->SetTitle("hwcf size [kB]");
350     TAxis* yaxis=histoNofClusters->GetYaxis();
351     if (yaxis) yaxis->SetTitle("N. of clusters");
352   }
353
354   fHistoHWCFDataSize=histoHWCFDataSize.release();
355   fHistoHWCFReductionFactor=histoHWCFReductionFactor.release();
356   fHistoNofClusters=histoNofClusters.release();
357
358   fpHWClusterDecoder=hwClusterDecoder.release();
359   fMonitoringContainer=dataContainer.release();
360
361   return iResult;
362 }
363
364 int AliHLTTPCDataCompressionMonitorComponent::DoDeinit()
365 {
366   /// inherited from AliHLTComponent: component cleanup
367   int iResult=0;
368
369   if (fpHWClusterDecoder) delete fpHWClusterDecoder;
370   fpHWClusterDecoder=NULL;
371
372   if (!fHistogramFile.IsNull()) {
373     TFile out(fHistogramFile, "RECREATE");
374     if (!out.IsZombie()) {
375       out.cd();
376       if (fHistoHWCFDataSize) fHistoHWCFDataSize->Write();
377       if (fHistoHWCFReductionFactor) fHistoHWCFReductionFactor->Write();
378       if (fHistoNofClusters) fHistoNofClusters->Write();
379       if (fMonitoringContainer) {
380         const TObject* o1=fMonitoringContainer->FindObject("fHistograms");
381         const TObject* o2=fMonitoringContainer->FindObject("fHistograms2D");
382         const TObject* o3=fMonitoringContainer->FindObject("fHistograms3D");
383         if (o1) o1->Write();
384         if (o2) o2->Write();
385         if (o3) o3->Write();
386       }
387       out.Close();
388     }
389   }
390   if (fHistoHWCFDataSize) delete fHistoHWCFDataSize;
391   fHistoHWCFDataSize=NULL;
392   if (fHistoHWCFReductionFactor) delete fHistoHWCFReductionFactor;
393   fHistoHWCFReductionFactor=NULL;
394   if (fHistoNofClusters) delete fHistoNofClusters;
395   fHistoNofClusters=NULL;
396   if (fMonitoringContainer) {
397     fMonitoringContainer->Clear();
398     delete fMonitoringContainer;
399   }
400   fMonitoringContainer=NULL;
401
402
403   return iResult;
404 }
405
406 int AliHLTTPCDataCompressionMonitorComponent::ScanConfigurationArgument(int argc, const char** argv)
407 {
408   /// inherited from AliHLTComponent: argument scan
409   int iResult=0;
410   if (argc<1) return 0;
411   int bMissingParam=0;
412   int i=0;
413   TString argument=argv[i];
414
415   do {
416     // -histogram-file
417     if (argument.CompareTo("-histogram-file")==0) {
418       if ((bMissingParam=(++i>=argc))) break;
419       fHistogramFile=argv[i++];
420       return i;
421     }
422     // -publishing-mode
423     if (argument.CompareTo("-publishing-mode")==0) {
424       if ((bMissingParam=(++i>=argc))) break;
425       TString option=argv[i++];
426       if (option.CompareTo("off")==0)           fPublishingMode=kPublishOff     ;
427       else if (option.CompareTo("separate")==0) fPublishingMode=kPublishSeparate;
428       else if (option.CompareTo("list")==0)     fPublishingMode=kPublishList    ;
429       else if (option.CompareTo("array")==0)    fPublishingMode=kPublishArray   ;
430       else {
431         HLTError("invalid option \"%s\" for argument \"%s\", expecting 'off', 'separate', 'list', or 'array'", option.Data(), argument.Data());
432         return -EPROTO;
433       }
434       return i;
435     }
436   } while (0); // using do-while only to have break available
437
438   if (bMissingParam) {
439     HLTError("missing parameter for argument %s", argument.Data());
440     iResult=-EPROTO;
441   }
442
443   return iResult;
444 }
445
446 AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::AliDataContainer()
447   : fHistograms(new TObjArray)  
448   , fHistograms2D(new TObjArray)  
449   , fHistograms3D(new TObjArray)    
450   , fHistogramPointers()
451   , fHistogram2DPointers()
452   , fHistogram3DPointers()
453   , fRemainingClusterIds()
454   , fTrackModelClusterIds()
455   , fCurrentClusterIds(NULL)
456   , fRawData(NULL)
457   , fLastPadRow(-1)
458   , fSector(-1)
459   , fBegin()
460 {
461   /// constructor
462   if (fHistograms) {
463     fHistograms->SetOwner(kTRUE);
464     fHistogramPointers.resize(kNumberOfHistograms, NULL);
465     for (const AliHistogramDefinition* definition=fgkHistogramDefinitions;
466          definition->fName!=NULL; definition++) {
467       fHistogramPointers[definition->fId]=new TH1F(definition->fName,
468                                                   definition->fTitle,
469                                                   definition->fBins,
470                                                   definition->fLowerBound,
471                                                   definition->fUpperBound
472                                                   );
473       fHistograms->AddAt(fHistogramPointers[definition->fId], definition->fId);
474     }
475   }
476   ///
477   if (fHistograms2D) {
478     fHistograms2D->SetOwner(kTRUE);
479     fHistogram2DPointers.resize(kNumberOfHistograms2D, NULL);
480     for (const AliHistogramDefinition2D* definition=fgkHistogramDefinitions2D;
481          definition->fName!=NULL; definition++) {
482       fHistogram2DPointers[definition->fId]=new TH2F(definition->fName,
483                                                      definition->fTitle,
484                                                      definition->fBinsX,
485                                                      definition->fLowerBoundX,
486                                                      definition->fUpperBoundX,
487                                                      definition->fBinsY,
488                                                      definition->fLowerBoundY,
489                                                      definition->fUpperBoundY
490                                                      );
491       fHistograms2D->AddAt(fHistogram2DPointers[definition->fId], definition->fId);
492     }
493   }
494   ///
495   if (fHistograms3D) {
496     fHistograms3D->SetOwner(kTRUE);
497     fHistogram3DPointers.resize(kNumberOfHistograms3D, NULL);
498     for (const AliHistogramDefinition3D* definition=fgkHistogramDefinitions3D;
499          definition->fName!=NULL; definition++) {
500       fHistogram3DPointers[definition->fId]=new TH3F(definition->fName,
501                                                      definition->fTitle,
502                                                      definition->fBinsX,
503                                                      definition->fLowerBoundX,
504                                                      definition->fUpperBoundX,
505                                                      definition->fBinsY,
506                                                      definition->fLowerBoundY,
507                                                      definition->fUpperBoundY,
508                                                      definition->fBinsZ,
509                                                      definition->fLowerBoundZ,
510                                                      definition->fUpperBoundZ
511                                                      );
512       fHistograms3D->AddAt(fHistogram3DPointers[definition->fId], definition->fId);
513     }
514   }
515   
516 }
517
518 const AliHLTTPCDataCompressionMonitorComponent::AliHistogramDefinition AliHLTTPCDataCompressionMonitorComponent::fgkHistogramDefinitions[] = {
519   {kHistogramPadrow,        "padrow"   , "padrow; padrow; counts"                  ,  159,   0.,   158.},
520   {kHistogramPad,           "pad"      , "pad; pad; counts"                        ,  140,   0.,   139.},
521   {kHistogramTime,          "timebin"  , "timebin; time; counts"                   , 1024,   0.,  1023.},
522   {kHistogramSigmaY2,       "sigmaY2"  , "sigmaY2; #sigma_{Y}^{2}; counts"         ,  100,   0.,     1.},
523   {kHistogramSigmaZ2,       "sigmaZ2"  , "sigmaZ2; #sigma_{Z}^{2}; counts"         ,  100,   0.,     1.},
524   {kHistogramCharge,        "charge"   , "charge; charge; counts"                  , 1024,   0., 65535.},
525   {kHistogramQMax,          "qmax"     , "qmax; Q_{max}; counts"                   ,  128,   0.,  1023.},
526   {kHistogramDeltaPadrow,   "d_padrow" , "d_padrow; #Delta padrow; counts"         , 1000,  -1.,     1.},
527   {kHistogramDeltaPad,      "d_pad"    , "d_pad; #Delta pad; counts"               , 1000,  -.1,     .1},
528   {kHistogramDeltaTime,     "d_time"   , "d_time; #Delta time; counts"             , 1000,  -.1,     .1},
529   {kHistogramDeltaSigmaY2,  "d_sigmaY2", "d_sigmaY2; #Delta #sigma_{Y}^{2}; counts", 1000,  -1.,     1.},
530   {kHistogramDeltaSigmaZ2,  "d_sigmaZ2", "d_sigmaZ2; #Delta #sigma_{Z}^{2}; counts", 1000,  -1.,     1.},
531   {kHistogramDeltaCharge,   "d_charge" , "d_charge; #Delta charge"                 , 1000,  -1.,     1.},
532   {kHistogramDeltaQMax,     "d_qmax"   , "d_qmax; #Delta Q_{max}"                  , 1000,  -1.,     1.},
533   {kHistogramOutOfRange,    "OutOfR"   , "OutOfR"                                  ,  159,   0.,   158.},
534   {kNumberOfHistograms, NULL, NULL, 0,0.,0.}
535 };
536
537  const AliHLTTPCDataCompressionMonitorComponent::AliHistogramDefinition2D AliHLTTPCDataCompressionMonitorComponent::fgkHistogramDefinitions2D[] = {
538    {kHistogramQMaxSector,    "qmaxsector"   , "qmaxsector; sector; Q_{max}"           , 72,0.,71., 128,0.,1023.},
539    {kHistogramSigmaY2Sector, "sigmaY2sector", "sigmaY2sector; sector; #sigma_{Y}^{2}" , 72,0.,71., 100,0.,1.},
540    {kHistogramSigmaZ2Sector, "sigmaZ2sector", "sigmaZ2sector; sector; #sigma_{Z}^{2}" , 72,0.,71., 100,0.,1.},
541    {kNumberOfHistograms2D, NULL, NULL, 0,0.,0., 0,0.,0.}
542  };
543
544  const AliHLTTPCDataCompressionMonitorComponent::AliHistogramDefinition3D AliHLTTPCDataCompressionMonitorComponent::fgkHistogramDefinitions3D[] = {
545    {kHistogramPadrowPadSector,"padrowpadsector","padrowpadsector; sector; pad;padrow", 72,0.,71., 140,0.,139., 159,0.,158.},
546    {kNumberOfHistograms3D, NULL, NULL, 0,0.,0., 0,0.,0., 0,0.,0.}
547  };
548
549 AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::~AliDataContainer()
550 {
551   /// dectructor
552   if (fRawData) delete fRawData;
553   if (fHistograms) delete fHistograms;
554   if (fHistograms2D) delete fHistograms2D;
555   if (fHistograms3D) delete fHistograms3D;
556 }
557
558 AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::iterator& AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::BeginRemainingClusterBlock(int /*count*/, AliHLTUInt32_t specification)
559 {
560   /// iterator of remaining clusters block of specification
561   AliHLTUInt8_t slice=AliHLTTPCDefinitions::GetMinSliceNr(specification);
562   AliHLTUInt8_t partition=AliHLTTPCDefinitions::GetMinPatchNr(specification);
563   unsigned index=slice*AliHLTTPCTransform::GetNumberOfPatches()+partition;
564   if (index<fRemainingClusterIds.size())
565     fCurrentClusterIds=&fRemainingClusterIds[index];
566   else
567     fCurrentClusterIds=NULL;
568   fBegin=iterator(this);
569   return fBegin;
570 }
571
572 AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::iterator& AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::BeginTrackModelClusterBlock(int /*count*/)
573 {
574   /// iterator of track model clusters
575   if (fTrackModelClusterIds.fIds && fTrackModelClusterIds.fSize>0)
576     fCurrentClusterIds=&fTrackModelClusterIds;
577   else
578     fCurrentClusterIds=NULL;
579   fBegin=iterator(this);
580   return fBegin;
581 }
582
583 int AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::AddRawData(const AliHLTComponentBlockData* pDesc)
584 {
585     /// add raw data bloack
586   if (pDesc->fDataType==AliHLTTPCDefinitions::HWClustersDataType()) {
587     if (!fRawData) fRawData=new AliHLTTPCHWCFSpacePointContainer(AliHLTTPCHWCFSpacePointContainer::kModeCreateMap);
588     if (!fRawData) return -ENOMEM;
589     return fRawData->AddInputBlock(pDesc);
590   }
591   return -ENODATA;  
592 }
593
594 int AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::AddClusterIds(const AliHLTComponentBlockData* pDesc)
595 {
596   /// add cluster id block for remaining or track model clusters
597   if (!pDesc) return -EINVAL;
598   if (pDesc->fDataType==AliHLTTPCDefinitions::ClusterIdTracksDataType()) {
599     fTrackModelClusterIds.fIds=reinterpret_cast<AliHLTUInt32_t*>(pDesc->fPtr);
600     fTrackModelClusterIds.fSize=pDesc->fSize/sizeof(AliHLTUInt32_t);
601     return 0;
602   }
603   if (pDesc->fDataType==AliHLTTPCDefinitions::RemainingClusterIdsDataType()) {
604     AliHLTUInt8_t slice=AliHLTTPCDefinitions::GetMinSliceNr(pDesc->fSpecification);
605     AliHLTUInt8_t partition=AliHLTTPCDefinitions::GetMinPatchNr(pDesc->fSpecification);
606     unsigned index=slice*AliHLTTPCTransform::GetNumberOfPatches()+partition;
607     if (fRemainingClusterIds.size()<=index) {
608       if ((int)fRemainingClusterIds.size()<AliHLTTPCTransform::GetNSlice()*AliHLTTPCTransform::GetNumberOfPatches()) {
609         fRemainingClusterIds.resize(AliHLTTPCTransform::GetNSlice()*AliHLTTPCTransform::GetNumberOfPatches());
610       } else {
611         fRemainingClusterIds.resize(index+1);
612       }
613     }
614     fRemainingClusterIds[index].fIds=reinterpret_cast<AliHLTUInt32_t*>(pDesc->fPtr);
615     fRemainingClusterIds[index].fSize=pDesc->fSize/sizeof(AliHLTUInt32_t);
616     return 0;
617   }
618   return -ENODATA;
619 }
620
621 AliHLTUInt32_t AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::GetClusterId(int clusterNo) const
622 {
623   /// get the cluster id from the current cluster id block (optional)
624   if (!fCurrentClusterIds ||
625       clusterNo<0 ||
626       (int)fCurrentClusterIds->fSize<=clusterNo)
627     return kAliHLTVoidDataSpec;
628   return fCurrentClusterIds->fIds[clusterNo];
629 }
630
631 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FillPadRow(int row, int slice, AliHLTUInt32_t clusterId)
632 {
633   /// fill padrow histogram
634   unsigned index=kHistogramPadrow;
635   fLastPadRow=row;
636   // the inner sectors consist of readout partitions 0 and 1, if the row
637   // is smaller than first row of readout partition 2, its an inner sector
638   if (row<AliHLTTPCTransform::GetFirstRow(2)) {
639     fSector = slice;
640   } else {
641     fSector = slice+AliHLTTPCTransform::GetNSlice();
642   }
643   if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL)
644     fHistogramPointers[index]->Fill(row);
645   if (clusterId!=kAliHLTVoidDataSpec) {
646     index=kHistogramDeltaPadrow;
647     if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
648       if (fRawData->Check(clusterId)) {
649         fHistogramPointers[index]->Fill(row-fRawData->GetX(clusterId));
650       }
651     }
652   }
653 }
654
655 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FillPad(float pad, AliHLTUInt32_t clusterId)
656 {
657   /// fill pad histogram
658   unsigned index=kHistogramPad;
659   if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL)
660     fHistogramPointers[index]->Fill(pad);
661
662   index=kHistogramPadrowPadSector;
663   if (index<fHistogram3DPointers.size() && fHistogram3DPointers[index]!=NULL)
664     fHistogram3DPointers[index]->Fill(fSector,pad,fLastPadRow);
665
666   if (clusterId!=kAliHLTVoidDataSpec) {
667     index=kHistogramDeltaPad;
668     if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
669       if (fRawData->Check(clusterId)) {
670         float dPad=pad-fRawData->GetY(clusterId);
671         fHistogramPointers[index]->Fill(dPad);
672         static const float maxdPad=0.015; // better 100um for 4 and 6mm pad width
673         if (TMath::Abs(dPad)>maxdPad) {
674           AliHLTUInt8_t slice = AliHLTTPCSpacePointData::GetSlice(clusterId);
675           AliHLTUInt8_t partition = AliHLTTPCSpacePointData::GetPatch(clusterId);
676           HLTError("cluster 0x%08x slice %d partition %d: pad difference %f - max %f", clusterId, slice, partition, dPad, maxdPad);
677           index=kHistogramOutOfRange;
678           if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
679             fHistogramPointers[index]->Fill(fLastPadRow>=0?fLastPadRow:0);
680           }
681         }
682       }
683     }
684   }
685 }
686
687 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FillTime(float time, AliHLTUInt32_t clusterId)
688 {
689   /// fill pad histogram
690   unsigned index=kHistogramTime;
691   if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL)
692     fHistogramPointers[index]->Fill(time);
693   if (clusterId!=kAliHLTVoidDataSpec) {
694     index=kHistogramDeltaTime;
695     if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
696       if (fRawData->Check(clusterId)) {
697         float dTime=time-fRawData->GetZ(clusterId);
698         fHistogramPointers[index]->Fill(dTime);
699         static const float maxdTime=0.04; // corresponds to 100um
700         if (TMath::Abs(dTime)>maxdTime) {
701           AliHLTUInt8_t slice = AliHLTTPCSpacePointData::GetSlice(clusterId);
702           AliHLTUInt8_t partition = AliHLTTPCSpacePointData::GetPatch(clusterId);
703           HLTError("cluster 0x%08x slice %d partition %d: time difference %f - max %f", clusterId, slice, partition, dTime, maxdTime);
704           index=kHistogramOutOfRange;
705           if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
706             fHistogramPointers[index]->Fill(fLastPadRow>=0?fLastPadRow:0);
707           }
708         }
709       }
710     }
711   }
712 }
713
714 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FillSigmaY2(float sigmaY2, AliHLTUInt32_t clusterId, int partition)
715 {
716   /// fill sigmaY2 histogram
717   unsigned index=kHistogramSigmaY2;
718   /// take account for different pad widths
719   float weight=AliHLTTPCTransform::GetPadPitchWidth(partition);
720   if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL)
721     fHistogramPointers[index]->Fill(sigmaY2*weight*weight);
722
723   index=kHistogramSigmaY2Sector;
724   if (index<fHistogram2DPointers.size() && fHistogram2DPointers[index]!=NULL)
725     fHistogram2DPointers[index]->Fill(fSector,sigmaY2*weight*weight);
726
727   if (clusterId!=kAliHLTVoidDataSpec) {
728     index=kHistogramDeltaSigmaY2;
729     if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
730       if (fRawData->Check(clusterId)) {
731         fHistogramPointers[index]->Fill(sigmaY2-fRawData->GetYWidth(clusterId));
732       }
733     }
734   }
735 }
736
737 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FillSigmaZ2(float sigmaZ2, AliHLTUInt32_t clusterId)
738 {
739   /// fill sigmaZ2 histogram
740   unsigned index=kHistogramSigmaZ2;
741   // FIXME: this is just a fixed value, to be correct the values from the global
742   // parameter block has to be used
743   float weight=AliHLTTPCTransform::GetZWidth();
744   if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL)
745     fHistogramPointers[index]->Fill(sigmaZ2*weight*weight);
746
747   index=kHistogramSigmaZ2Sector;
748   if (index<fHistogram2DPointers.size() && fHistogram2DPointers[index]!=NULL)
749     fHistogram2DPointers[index]->Fill(fSector,sigmaZ2*weight*weight);
750
751   if (clusterId!=kAliHLTVoidDataSpec) {
752     index=kHistogramDeltaSigmaZ2;
753     if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
754       if (fRawData->Check(clusterId)) {
755         fHistogramPointers[index]->Fill(sigmaZ2-fRawData->GetZWidth(clusterId));
756       }
757     }
758   }
759 }
760
761 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FillCharge(unsigned charge, AliHLTUInt32_t clusterId)
762 {
763   /// fill charge histogram
764   unsigned index=kHistogramCharge;
765   if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL)
766     fHistogramPointers[index]->Fill(charge);
767   if (clusterId!=kAliHLTVoidDataSpec) {
768     index=kHistogramDeltaCharge;
769     if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
770       if (fRawData->Check(clusterId)) {
771         fHistogramPointers[index]->Fill(charge-fRawData->GetCharge(clusterId));
772       }
773     }
774   }
775 }
776
777 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FillQMax(unsigned qmax, AliHLTUInt32_t clusterId)
778 {
779   /// fill qmax histogram
780   unsigned index=kHistogramQMax;
781   if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL)
782     fHistogramPointers[index]->Fill(qmax);
783
784   index=kHistogramQMaxSector;
785   if (index<fHistogram2DPointers.size() && fHistogram2DPointers[index]!=NULL)
786     fHistogram2DPointers[index]->Fill(fSector,qmax);
787
788   if (clusterId!=kAliHLTVoidDataSpec) {
789     index=kHistogramDeltaQMax;
790     if (index<fHistogramPointers.size() && fHistogramPointers[index]!=NULL && fRawData) {
791       if (fRawData->Check(clusterId)) {
792         fHistogramPointers[index]->Fill(qmax-fRawData->GetQMax(clusterId));
793       }
794     }
795   }
796 }
797
798 void AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::Clear(Option_t * option)
799 {
800   /// internal cleanup
801   if (fRawData) fRawData->Clear(option);
802 }
803
804 TObject* AliHLTTPCDataCompressionMonitorComponent::AliDataContainer::FindObject(const char *name) const
805 {
806   /// get histogram object  
807   if (!name) return NULL;
808   if ( strcmp(name,"fHistograms")   == 0 )
809     return fHistograms;
810   if ( strcmp(name,"fHistograms2D") == 0 )
811     return fHistograms2D;
812   if ( strcmp(name,"fHistograms3D") == 0 )
813     return fHistograms3D;
814
815   return NULL;
816 }