]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTDataDeflaterSimple.cxx
Update master to aliroot
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTDataDeflaterSimple.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   AliHLTDataDeflaterSimple.cxx
20 /// @author Matthias Richter
21 /// @date   2011-08-10
22 /// @brief  Simple deflater implementation storing frequent values below a
23 ///         maximum value with a reduced bit number and others with the full
24 ///         number of bits.
25
26 #include "AliHLTDataDeflaterSimple.h"
27 #include <memory>
28 #include <algorithm>
29 #include <iostream>
30
31 /** ROOT macro for the implementation of ROOT specific class methods */
32 ClassImp(AliHLTDataDeflaterSimple)
33
34 AliHLTDataDeflaterSimple::AliHLTDataDeflaterSimple()
35   : AliHLTDataDeflater()
36   , fParameterDefinitions()
37 {
38   // see header file for class documentation
39   // or
40   // refer to README to build package
41   // or
42   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
43 }
44
45 AliHLTDataDeflaterSimple::~AliHLTDataDeflaterSimple()
46 {
47   // destructor
48   Clear();
49
50 }
51
52 int AliHLTDataDeflaterSimple::AddParameterDefinition(const char* name, int bitLength, int reducedBitLength)
53 {
54   /// add a parameter definition to the configuration, return reference id
55   fParameterDefinitions.push_back(AliHLTDataDeflaterParameter(name, bitLength, reducedBitLength));
56   int memberId=fParameterDefinitions.size()-1;
57   if (DoStatistics()) {
58     AddHistogram(memberId, name, bitLength);
59   }
60   return memberId;
61 }
62
63 bool AliHLTDataDeflaterSimple::OutputParameterBits( int memberId, AliHLTUInt64_t const & value )
64 {
65   // write bit pattern of a member to the current byte and position
66   if (memberId>=(int)fParameterDefinitions.size()) return false;
67
68   AliHLTUInt32_t switchBit=fParameterDefinitions[memberId].SwitchBit(value); // 0 -> reduced, 1 -> full
69   AliHLTUInt64_t v=fParameterDefinitions[memberId].Value(value);
70   AliHLTUInt32_t length=fParameterDefinitions[memberId].ValueLength(value);
71   fParameterDefinitions[memberId].IncrementBitCount(value);
72
73   if (DoStatistics())
74     FillStatistics(memberId, length, value);
75
76   if (!OutputBit(switchBit)) return false;
77   return OutputBits(v, length);
78 }
79
80 void AliHLTDataDeflaterSimple::Clear(Option_t * option)
81 {
82   // internal cleanup
83   unsigned i=0;
84   for (vector<AliHLTDataDeflaterParameter>::iterator m=fParameterDefinitions.begin();
85        m!=fParameterDefinitions.end(); m++, i++) {
86     m->ResetBitCount();
87   }
88   AliHLTDataDeflater::Clear(option);
89 }
90
91 void AliHLTDataDeflaterSimple::Print(Option_t *option) const
92 {
93   // print info
94   Print(cout, option);
95 }
96
97 void AliHLTDataDeflaterSimple::Print(ostream& out, Option_t *option) const
98 {
99   // print to stream
100   out << "AliHLTDataDeflaterSimple:" << endl;
101   AliHLTUInt64_t bitCount=0;
102   AliHLTUInt64_t fullSize=0;
103   for (vector<AliHLTDataDeflaterParameter>::const_iterator m=fParameterDefinitions.begin();
104        m!=fParameterDefinitions.end(); m++) {
105     cout << "   "; m->Print(option);
106     bitCount+=m->GetBitCount();
107     fullSize+=m->GetValueCount()*m->GetBitLength();
108   }
109   out << " total: " << bitCount << "/" << fullSize << " " << (fullSize>0?float(bitCount)/fullSize:0.0) << endl;
110 }
111
112 void AliHLTDataDeflaterSimple::AliHLTDataDeflaterParameter::Print(const char* /*option*/) const
113 {
114   // print info
115   cout << fName << " (" << fFullBitLength << "," << fReducedBitLength << "): "
116        << fValueCount << " entries  "
117        << fBitCount << "/" << fFullBitLength*fValueCount;
118   if (fFullBitLength && fValueCount) {
119     cout << " " << float(fBitCount)/(fValueCount*fFullBitLength);
120   }
121   cout << endl;
122 }
123
124 ostream& operator<<(ostream &out, const AliHLTDataDeflaterSimple& me)
125 {
126   me.Print(out);
127   return out;
128 }