]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTDataDeflater.cxx
bugfix: correct range of DDL for specified detector
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTDataDeflater.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /// @file   AliHLTDataDeflater.cxx
20 /// @author Matthias Richter, Timm Steinbeck
21 /// @date   2011-08-10
22 /// @brief  Data deflater class storing only necessary bits
23 /// @note   Code original from AliHLTTPCCompModelDeflater
24
25 #include "AliHLTDataDeflater.h"
26 #include "AliHLTErrorGuard.h"
27 #include "TObjArray.h"
28 #include "TH1I.h"
29 #include "TH2D.h"
30 #include "TFile.h"
31 #include "TMath.h"
32 #include <memory>
33 #include <algorithm>
34 #include <iostream>
35
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTDataDeflater)
38
39 AliHLTDataDeflater::AliHLTDataDeflater()
40   : AliHLTLogging()
41   , fBitDataCurrentWord(0)
42   , fBitDataCurrentPosInWord(0)
43   , fBitDataCurrentOutput(NULL)
44   , fBitDataCurrentOutputStart(NULL)
45   , fBitDataCurrentOutputEnd(NULL)
46   , fHistograms(NULL)
47   , fParameterCompression(NULL)
48   , fParameterSize(NULL)
49 {
50   // see header file for class documentation
51   // or
52   // refer to README to build package
53   // or
54   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
55   if (fHistograms) fHistograms->SetOwner(kTRUE);
56 }
57
58 AliHLTDataDeflater::~AliHLTDataDeflater()
59 {
60   // destructor
61   Clear();
62
63   if (fHistograms)
64     delete fHistograms;
65   fHistograms=NULL;
66   if (fParameterCompression)
67     delete fParameterCompression;
68   fParameterCompression=NULL;
69   if (fParameterSize)
70     delete fParameterSize;
71   fParameterSize=NULL;
72 }
73
74 int AliHLTDataDeflater::InitBitDataOutput( AliHLTUInt8_t* output, UInt_t outputSize)
75 {
76   // init the target buffer
77   fBitDataCurrentWord = 0;
78   fBitDataCurrentPosInWord = 7;
79   fBitDataCurrentOutput = fBitDataCurrentOutputStart = output;
80   fBitDataCurrentOutputEnd = output+outputSize;
81
82   return 0;
83 }
84
85 void AliHLTDataDeflater::CloseBitDataOutput()
86 {
87   // pad to full byte and clear internal pointer references
88   Pad8Bits();
89   fBitDataCurrentWord=0;
90   fBitDataCurrentPosInWord=0;
91   fBitDataCurrentOutput=NULL;
92   fBitDataCurrentOutputStart=NULL;
93   fBitDataCurrentOutputEnd=NULL;
94 }
95
96 AliHLTUInt8_t AliHLTDataDeflater::GetCurrentOutputByte( Int_t offset ) const
97 {
98   // get the current byte
99   if ( !offset )
100     return fBitDataCurrentWord;
101   else
102     return *(fBitDataCurrentOutput+offset);
103 }
104
105 bool AliHLTDataDeflater::OutputBit( AliHLTUInt32_t const & value )
106 {
107   // write one bit to the current byte and position
108   if ( fBitDataCurrentOutput>=fBitDataCurrentOutputEnd )
109     return false;
110   fBitDataCurrentWord |= (value & 1) << fBitDataCurrentPosInWord;
111   if ( fBitDataCurrentPosInWord )
112     fBitDataCurrentPosInWord--;
113   else {
114     *fBitDataCurrentOutput = fBitDataCurrentWord;
115     fBitDataCurrentPosInWord = 7;
116     fBitDataCurrentOutput++;
117     fBitDataCurrentWord = 0;
118   }
119   return true;
120 }
121
122 bool AliHLTDataDeflater::OutputBits( AliHLTUInt64_t const & value, UInt_t const & bitCount )
123 {
124   // write bit pattern to the current byte and position
125   if ( bitCount>64 ) {
126     HLTFatal( "Internal error: Attempt to write more than 64 bits (%u)", (unsigned)bitCount );
127     return false;
128   }
129   UInt_t bitsToWrite=bitCount;
130   UInt_t curBitCount;
131   while ( bitsToWrite>0 ) {
132     if ( fBitDataCurrentOutput>=fBitDataCurrentOutputEnd )
133       return false;
134 #if 1
135     if ( bitsToWrite >= fBitDataCurrentPosInWord+1 )
136       curBitCount = fBitDataCurrentPosInWord+1;
137     else
138       curBitCount = bitsToWrite;
139     fBitDataCurrentWord |= ( (value >> (bitsToWrite-curBitCount)) & ((1<<curBitCount)-1) ) << (fBitDataCurrentPosInWord+1-curBitCount);
140     if ( fBitDataCurrentPosInWord < curBitCount )
141       {
142         *fBitDataCurrentOutput = fBitDataCurrentWord;
143         fBitDataCurrentPosInWord = 7;
144         fBitDataCurrentOutput++;
145         fBitDataCurrentWord = 0;
146       }
147     else
148       fBitDataCurrentPosInWord -= curBitCount;
149     bitsToWrite -= curBitCount;
150
151 #else
152     AliHLTUInt8_t curValue;
153     if ( bitsToWrite>=8 )
154       {
155         curBitCount=8;
156         curValue = (value >> bitsToWrite-8) & 0xFF;
157         bitsToWrite -= 8;
158       }
159     else
160       {
161         curBitCount=bitsToWrite;
162         curValue = value & ( (1<<bitsToWrite)-1 );
163         bitsToWrite = 0;
164       }
165     if ( fBitDataCurrentPosInWord+1>curBitCount )
166       {
167         fBitDataCurrentWord |= curValue << (fBitDataCurrentPosInWord-curBitCount+1);
168         fBitDataCurrentPosInWord -= curBitCount;
169       }
170     else if ( fBitDataCurrentPosInWord+1==curBitCount )
171       {
172         fBitDataCurrentWord |= curValue;
173         *fBitDataCurrentOutput = fBitDataCurrentWord;
174         fBitDataCurrentPosInWord = 7;
175         fBitDataCurrentOutput++;
176         fBitDataCurrentWord = 0;
177       }
178     else
179       {
180         const UInt_t first = fBitDataCurrentPosInWord+1; // Number of bits for first block
181         const UInt_t second = curBitCount-first; // Number of bits for second block
182         fBitDataCurrentWord |= ( curValue >> second ) & ((1<<first)-1);
183         *fBitDataCurrentOutput = fBitDataCurrentWord;
184         fBitDataCurrentOutput++;
185         if ( fBitDataCurrentOutput>=fBitDataCurrentOutputEnd )
186           return false;
187         fBitDataCurrentWord = curValue & ((1<<second)-1) << (8-second);
188         fBitDataCurrentPosInWord = 7-second;
189       }
190 #endif
191   }
192   return true;
193 }
194
195 bool AliHLTDataDeflater::OutputBits( std::bitset<64> const & value, UInt_t const & bitCount )
196 {
197   // write bit pattern to the current byte and position
198   if ( bitCount>64 ) {
199     HLTFatal( "Internal error: Attempt to write more than 64 bits (%u)", (unsigned)bitCount );
200     return false;
201   }
202   static const std::bitset<64> mask8bit(255ul);
203   UInt_t bitsToWrite=bitCount;
204   UInt_t curBitCount;
205   while ( bitsToWrite>0 ) {
206     if ( fBitDataCurrentOutput>=fBitDataCurrentOutputEnd )
207       return false;
208     if ( bitsToWrite >= fBitDataCurrentPosInWord+1 )
209       curBitCount = fBitDataCurrentPosInWord+1;
210     else
211       curBitCount = bitsToWrite;
212     std::bitset<64> valwrite=(value >> (bitsToWrite-curBitCount)) & mask8bit;
213     fBitDataCurrentWord |= ( valwrite.to_ulong() & ((1<<curBitCount)-1) ) << (fBitDataCurrentPosInWord+1-curBitCount);
214     if ( fBitDataCurrentPosInWord < curBitCount )
215       {
216         *fBitDataCurrentOutput = fBitDataCurrentWord;
217         fBitDataCurrentPosInWord = 7;
218         fBitDataCurrentOutput++;
219         fBitDataCurrentWord = 0;
220       }
221     else
222       fBitDataCurrentPosInWord -= curBitCount;
223     bitsToWrite -= curBitCount;
224   }
225   return true;
226 }
227
228 void AliHLTDataDeflater::Pad8Bits()
229 {
230   // finish the current word
231   if ( fBitDataCurrentPosInWord==7 )
232     return;
233   *fBitDataCurrentOutput = fBitDataCurrentWord;
234   fBitDataCurrentPosInWord = 7;
235   fBitDataCurrentOutput++;
236   fBitDataCurrentWord = 0;
237 }
238
239 bool AliHLTDataDeflater::OutputBytes( AliHLTUInt8_t const * data, UInt_t const & byteCount )
240 {
241   // write sequence of bytes
242   Pad8Bits();
243   if ( fBitDataCurrentOutput+byteCount>fBitDataCurrentOutputEnd )
244     return false;
245   memcpy( fBitDataCurrentOutput, data, byteCount );
246   fBitDataCurrentOutput += byteCount;
247   return true;
248 }
249
250 bool AliHLTDataDeflater::OutputParameterBits( int parameterId, AliHLTUInt64_t const & value )
251 {
252   // write bit pattern of a member to the current byte and position
253   return OutputParameterBits(parameterId, value, 0);
254 }
255
256 bool AliHLTDataDeflater::OutputParameterBits( int /*(parameterId*/, AliHLTUInt64_t const & /*value*/ , int /*lengthOffset*/)
257 {
258   // write bit pattern of a member to the current byte and position
259   ALIHLTERRORGUARD(1,"method needs to be implemented in child class");
260   return false;
261 }
262
263 void AliHLTDataDeflater::Clear(Option_t * /*option*/)
264 {
265   // internal cleanup
266 }
267
268 void AliHLTDataDeflater::Print(Option_t *option) const
269 {
270   // print info
271   Print(cout, option);
272 }
273
274 void AliHLTDataDeflater::Print(ostream& out, Option_t */*option*/) const
275 {
276   // print to stream
277   out << "AliHLTDataDeflater: " << endl;
278 }
279
280 ostream& operator<<(ostream &out, const AliHLTDataDeflater& me)
281 {
282   me.Print(out);
283   return out;
284 }
285
286 int AliHLTDataDeflater::EnableStatistics()
287 {
288   /// enable statistics accounting
289   if (!fHistograms) {
290     fHistograms=new TObjArray;
291     if (!fHistograms) return -ENOMEM;
292     fHistograms->SetOwner(kTRUE);
293   }
294   return 0;
295 }
296
297 int AliHLTDataDeflater::AddHistogram(int id, const char* name, int bitWidth, TH1* h)
298 {
299   /// add a histogram for deflater statistic of the corresponding parameter
300   /// a histogram is created with default settings if h is not specified; if
301   /// provided, the ownership goes over to the base class
302   if (!fHistograms) {
303     fHistograms=new TObjArray;
304     if (!fHistograms) return -ENOMEM;
305     fHistograms->SetOwner(kTRUE);
306   }
307   if (id>=0 && fHistograms->GetEntriesFast()>id && fHistograms->At(id)!=NULL) {
308     HLTWarning("parameter with id %d has existing object (%s), skipping histogram %s",
309                id, fHistograms->At(id)->GetName(), h?h->GetName():name);
310     return -EEXIST;
311   }
312   if (id<0 && h!=NULL && fHistograms->FindObject(h->GetName())) {
313     HLTWarning("parameter with name %s already existing, skipping histogram", h->GetName());
314     return -EEXIST;
315   }
316   if (!h)
317     h=new TH1I(name, name, 100, 0, 1<<bitWidth);
318   if (!h) return -ENOMEM;
319   fHistograms->Add(h);
320   return 0;
321 }
322
323 int AliHLTDataDeflater::FillStatistics(int id, AliHLTUInt64_t value, unsigned length, float codingWeight)
324 {
325   /// fill statistics for a parameter
326   if (!fHistograms ||
327       fHistograms->GetEntriesFast()<=id ||
328       id<0) return 0;
329
330   TObject* o=fHistograms->At(id);
331   if (o) {
332     TH1* h=dynamic_cast<TH1*>(o);
333     if (h) {
334       h->Fill(value);
335     }
336   }
337
338   if (!fParameterCompression) {
339     int bins=fHistograms->GetEntriesFast();
340     fParameterCompression=new TH2D("ParameterCompression", "ParameterCompression", bins, 0, bins, 1000, 0., 5.0);
341   }
342   if (fParameterCompression) {
343     fParameterCompression->Fill(id, codingWeight);
344   }
345   if (!fParameterSize) {
346     int bins=fHistograms->GetEntriesFast();
347     fParameterSize=new TH2D("ParameterSize", "ParameterSize", bins, 0, bins, 1000, 0., 64.0);
348   }
349   if (fParameterSize) {
350     fParameterSize->Fill(id, length);
351   }
352
353   return 0;
354 }
355
356 void AliHLTDataDeflater::SaveAs(const char *filename,Option_t */*option*/) const
357 {
358   // safe histograms to file
359   std::auto_ptr<TFile> file(TFile::Open(filename, "RECREATE"));
360   if (!file.get() || file->IsZombie()) {
361     HLTError("can not open file %s", filename);;
362     return;
363   }
364   file->cd();
365   if (fHistograms) {
366     for (int i=0; i<fHistograms->GetEntries(); i++) {
367       if (fHistograms->At(i)==NULL || 
368           !fHistograms->At(i)->InheritsFrom("TH1") ||
369           fHistograms->At(i)->InheritsFrom("TH2") ||
370           fHistograms->At(i)->InheritsFrom("TH3")
371           ) continue; // only TH1 objects in the list
372       TH1* h=reinterpret_cast<TH1*>(fHistograms->At(i));
373       if (!h) continue;
374       float entropy=CalcEntropy(h);
375       if (entropy<0) continue;
376       TString title=h->GetTitle();
377       title+=Form(" entropy %.2f", entropy);
378       h->SetTitle(title);
379     }
380     fHistograms->Write();
381     if (fParameterCompression)
382       fParameterCompression->Write();
383     if (fParameterSize)
384       fParameterSize->Write();
385   }
386
387   file->Close();
388 }
389
390 float AliHLTDataDeflater::CalcEntropy(TH1* histo, const char* /*option*/, int mode)
391 {
392
393   if (!histo) return -1000.;
394
395   float l2=TMath::Log(2.0);
396   float integral=histo->Integral(0,histo->GetNbinsX());
397   int centerbin=mode*histo->GetNbinsX()/2;
398   int nofBins=histo->GetNbinsX()-centerbin;
399   float entropy=0.0;
400   for (int offset=0; offset<nofBins; offset++) {
401     float abundance=histo->GetBinContent(offset);
402     if (abundance<1.0) continue;
403     entropy += (- (Double_t) abundance / (Double_t) integral ) * log( ( (Double_t) abundance / (Double_t) integral )) / (l2);
404   }
405
406   return entropy;
407 }