1 /**************************************************************************
2 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
19 // This class contains the implementation of the
20 // compression and decompression algorithms
21 // Compression is performed reading the Altro data block (called packet) backward.
22 // Similarly decompression is also done backward so that the final file is restored
23 // after the compression and decompression phase.
26 #include <TObjArray.h>
27 #include <Riostream.h>
30 #include "AliAltroBuffer.h"
31 #include "AliTPCHNode.h"
32 #include "AliTPCHTable.h"
33 #include "AliTPCCompression.h"
35 ClassImp(AliTPCCompression)
37 //////////////////////////////////////////////////////////////////////////////////////////////////
38 AliTPCCompression::AliTPCCompression(){
40 fDimBuffer=sizeof(UInt_t)*8;
41 fFreeBitsBuffer=fDimBuffer;
50 //////////////////////////////////////////////////////////////////////////////////////////////////
51 AliTPCCompression::AliTPCCompression(const AliTPCCompression &source)
54 this->fDimBuffer=source.fDimBuffer;
55 this->fFreeBitsBuffer=source.fFreeBitsBuffer;
56 this->fBitsToRead=source.fBitsToRead;
57 this->fPos=source.fPos;
58 this->fBuffer=source.fBuffer;
59 this->fVerbose=source.fVerbose;
60 this->fFillWords=source.fFillWords;
61 this->fPointBuffer=source.fPointBuffer;
64 //////////////////////////////////////////////////////////////////////////////////////////////////
65 AliTPCCompression& AliTPCCompression::operator=(const AliTPCCompression &source){
66 //Redefinition of the assignment operator
67 this->fDimBuffer=source.fDimBuffer;
68 this->fFreeBitsBuffer=source.fFreeBitsBuffer;
69 this->fBitsToRead=source.fBitsToRead;
70 this->fPos=source.fPos;
71 this->fBuffer=source.fBuffer;
72 this->fVerbose=source.fVerbose;
73 this->fFillWords=source.fFillWords;
74 this->fPointBuffer=source.fPointBuffer;
77 //////////////////////////////////////////////////////////////////////////////////////////////////
78 void AliTPCCompression::NextTable(Int_t Val,Int_t &NextTableType,Int_t &BunchLen,Int_t &Count)const{
79 //Depending on the data type (5 types of data) a specific table is called
82 0==> Bunch length value
88 switch (NextTableType){
95 if (BunchLen==1)NextTableType=2;
108 if (Count==(BunchLen-1)){
129 /////////////////////////////////////////////////////////////////////////////////////////////////////
131 Int_t AliTPCCompression::FillTables(const char* fSource,AliTPCHTable* table[],Int_t /*NumTables*/){
132 //This method is used to compute the frequencies of the symbols in the source file
133 AliAltroBuffer buff(fSource,0);
135 UInt_t countTrailer=0;
136 Int_t numWords,padNum,rowNum,secNum=0;
138 UInt_t stat[5]={0,0,0,0,0};
141 while(buff.ReadTrailerBackward(numWords,padNum,rowNum,secNum)){
143 endFill=buff.GetFillWordsNum();
148 fFillWords+=4-numWords%4;
149 for(Int_t j=0;j<(4-numWords%4);j++){
150 value=buff.GetNextBackWord();
157 for(Int_t i=0;i<345;i++)timePos[i]=0;
158 for(Int_t i=0;i<1024;i++)packet[i]=0;
160 Int_t nextTableType=0;
163 for(Int_t i=0;i<numWords;i++){
164 value=buff.GetNextBackWord();
166 if(nextTableType==1){
170 NextTable(value,nextTableType,bunchLen,count);
172 //computing the Time gap between two bunches
175 Int_t previousTime=packet[timePos[tp]];
176 for(Int_t i=tp-1;i>=0;i--){
177 Int_t timPos=timePos[i];
178 Int_t bunchLen=packet[timPos-1]-2;
180 packet[timPos]=packet[timPos]-previousTime-bunchLen;
186 for(Int_t i=0;i<numWords;i++){
188 table[nextTableType]->SetFrequency(value);
189 stat[nextTableType]++;
190 NextTable(value,nextTableType,bunchLen,count);
194 cout<<"Number of words: "<<countWords<<endl;
195 cout<<"Number of trailers: "<<countTrailer<<endl;
196 cout<<"Number of fill words "<<fFillWords+endFill<<endl;
197 cout<<"Total number of words: "<<countWords+countTrailer*4+fFillWords<<endl;
199 fStat.open("Statistics");
200 fStat<<"Number of words:..........................................."<<countWords<<endl;
201 fStat<<"Number of trailers (4 10 bits words in each one)..........."<<countTrailer<<endl;
202 fStat<<"Number of fill words:......................................"<<fFillWords+endFill<<endl;
203 fStat<<"Total number of words:....................................."<<countWords+countTrailer*4+fFillWords+endFill<<endl;
204 fStat<<"-----------------------------------------"<<endl;
205 fStat<<"Number of Bunches............."<<stat[0]<<endl;
206 fStat<<"Number of Time bin............"<<stat[1]<<endl;
207 fStat<<"Number of One Samples Bunch..."<<stat[2]<<endl;
208 fStat<<"Number of Central Samples....."<<stat[3]<<endl;
209 fStat<<"Number of Border Samples......"<<stat[4]<<endl;
210 fStat<<"-----------------------------------------"<<endl;
211 UInt_t fileDimension=(UInt_t)TMath::Ceil(double((countTrailer*4+countWords+fFillWords+endFill)*10/8));
212 fStat<<"Total file Size in bytes.."<<fileDimension<<endl;
213 Double_t percentage=TMath::Ceil((fFillWords+endFill)*125)/fileDimension;
214 fStat<<"Fill Words................"<<(UInt_t)TMath::Ceil((fFillWords+endFill)*10/8)<<" bytes "<<percentage<<"%"<<endl;
215 percentage=(Double_t)countTrailer*500/fileDimension;
216 fStat<<"Trailer..................."<<countTrailer*5<<" bytes "<<percentage<<"%"<<endl;
218 percentage=(Double_t)((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])) *125/fileDimension;
219 fStat<<"Data......................"<<(UInt_t)TMath::Ceil((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])*10/8)<<" bytes "<<percentage<<"%"<<endl;
221 percentage=(Double_t)(stat[0]*125)/fileDimension;
222 fStat<<"Bunch....................."<<(UInt_t)TMath::Ceil(stat[0]*10/8)<<" bytes "<<percentage<<"%"<<endl; //
223 percentage=(Double_t)(stat[1]*125)/fileDimension;
224 fStat<<"Time......................"<<(UInt_t)TMath::Ceil(stat[1]*10/8)<<" bytes "<<percentage<<"%"<<endl; //
227 percentage=(Double_t)((stat[2]+stat[3]+stat[4])) *125/fileDimension;
228 fStat<<"Amplitude values.........."<<(UInt_t)TMath::Ceil((stat[2]+stat[3]+stat[4])*10/8)<<" bytes "<<percentage<<"%"<<endl;
229 percentage=(Double_t)(stat[2]*125)/fileDimension;
230 fStat<<" One Samples..............."<<(UInt_t)TMath::Ceil(stat[2]*10/8)<<" bytes "<<percentage<<"%"<<endl; //
231 percentage=(Double_t)(stat[3]*125)/fileDimension;
232 fStat<<" Central Samples..........."<<(UInt_t)TMath::Ceil(stat[3]*10/8)<<" bytes "<<percentage<<"%"<<endl; //
233 percentage=(Double_t)(stat[4]*125)/fileDimension;
234 fStat<<" Border Samples............"<<(UInt_t)TMath::Ceil(stat[4]*10/8)<<" bytes "<<percentage<<"%"<<endl; //
238 ////////////////////////////////////////////////////////////////////////////////////////
239 Int_t AliTPCCompression::StoreTables(AliTPCHTable* table[],const Int_t NumTable){
240 //This method stores the tables in a sequence of binary file
243 for(Int_t k=0;k<NumTable;k++){
244 sprintf(filename,"Table%d.dat",k);
246 fTable.open(filename,ios::binary);
248 fTable.open(filename);
250 Int_t dim=table[k]->Size();
251 //Table dimension is written into a file
252 fTable.write((char*)(&dim),sizeof(Int_t));
253 //One table is written into a file
254 for(Int_t i=0;i<dim;i++){
255 UChar_t codeLen=table[k]->CodeLen()[i];
256 // UInt_t code=(UInt_t)table[k]->Code()[i];
257 Double_t code=table[k]->Code()[i];
258 fTable.write((char*)(&codeLen),sizeof(UChar_t));
259 //fTable.write((char*)(&code),sizeof(UInt_t));
260 fTable.write((char*)(&code),sizeof(Double_t));
266 ////////////////////////////////////////////////////////////////////////////////////////
267 Int_t AliTPCCompression::CreateTableFormula(Double_t beta,UInt_t M,Int_t dim,Int_t Type){
268 // Type = 0 for Bunch length
269 // Type = 1 for Time Gap
275 AliTPCHTable *table=new AliTPCHTable(dim);
278 Double_t freqArray[1024];
279 for(Int_t i=0;i<1024;i++){
282 alpha=M*0.000000602+0.0104;
284 cout<<"alpha "<<alpha<<endl;
285 for(Int_t x=0;x<dim;x++){
287 freqArray[x]=TMath::Power((x+1),-beta)*TMath::Exp(-alpha*(x+1));
289 freqArray[x]=TMath::Power((x+1),-beta);
291 if (freqArray[x]<min)min=freqArray[x];
294 cout<<"Minimun Value "<<min<<endl;
297 cout<<"a Value: "<<a<<endl;
298 for(Int_t x=0;x<dim;x++){
299 if (Type==0)//Bunch length
300 if (x>=3)//minimum bunch length
301 table->SetValFrequency(x,a*freqArray[x]*1000);
303 table->SetValFrequency(x,0);
305 table->SetValFrequency(x,a*freqArray[x]);
307 table->BuildHTable();
310 sprintf(filename,"Table%d.dat",Type);
312 fTable.open(filename,ios::binary);
314 fTable.open(filename);
316 Int_t dimTable=table->Size();
317 //Table dimension is written into a file
318 fTable.write((char*)(&dimTable),sizeof(Int_t));
319 //One table is written into a file
320 for(Int_t i=0;i<dimTable;i++){
321 UChar_t codeLen=table->CodeLen()[i];
322 Double_t code=table->Code()[i];
323 fTable.write((char*)(&codeLen),sizeof(UChar_t));
324 fTable.write((char*)(&code),sizeof(Double_t));
330 ////////////////////////////////////////////////////////////////////////////////////////
331 Int_t AliTPCCompression::CreateTables(const char* fSource,Int_t NumTables){
335 0==> Bunch length values
341 Int_t n=10;// 10 bits per symbol
342 AliTPCHTable ** table = new AliTPCHTable*[NumTables];
343 //The table is inizialized with the rigth number of rows
344 for(Int_t i=0;i<NumTables;i++){
345 table[i]=new AliTPCHTable((Int_t)(TMath::Power(2,n)));
346 table[i]->SetVerbose(fVerbose);
348 //The frequencies are calculated and the tables are filled
350 cout<<"Filling tables...\n";
351 //The get the frequencies
352 FillTables(fSource,table,NumTables);
354 //This part will be used in the table optimization phase
356 for(Int_t i=0;i<NumTables;i++){
357 table[i]->CompleteTable(i);
361 cout<<"Entropy of Bunch length table........."<<table[0]->GetEntropy()<<endl;
362 cout<<"Entropy of Time bin table............."<<table[1]->GetEntropy()<<endl;
363 cout<<"Entropy of one Sample bunch table....."<<table[2]->GetEntropy()<<endl;
364 cout<<"Entropy of Central Sample table......."<<table[3]->GetEntropy()<<endl;
365 cout<<"Entropy Border Samples table.........."<<table[4]->GetEntropy()<<endl;
367 fStat.open("Statistics",ios::app);
369 fStat<<"----------------- ENTROPY for castomized tables --------------------------"<<endl;
370 fStat<<"Entropy of Bunch length table......."<<table[0]->GetEntropy()<<endl;
371 fStat<<"Entropy of Time bin table..........."<<table[1]->GetEntropy()<<endl;
372 fStat<<"Entropy of one Sample bunch table..."<<table[2]->GetEntropy()<<endl;
373 fStat<<"Entropy of Central Sample table....."<<table[3]->GetEntropy()<<endl;
374 fStat<<"Entropy Border Samples table........"<<table[4]->GetEntropy()<<endl;
378 cout<<"Tables filled \n";
380 //Frequencies normalization
381 table[0]->NormalizeFrequencies();
382 table[1]->NormalizeFrequencies();
383 table[2]->NormalizeFrequencies();
384 table[3]->NormalizeFrequencies();
385 table[4]->NormalizeFrequencies();
387 //Tables are saved in a sequence of text file and using the macro Histo.C is it possible to get
388 //a series of histograms rappresenting the frequency distribution
389 table[0]->StoreFrequencies("BunchLenFreq.txt");
390 table[1]->StoreFrequencies("TimeFreq.txt");
391 table[2]->StoreFrequencies("Sample1Freq.txt");
392 table[3]->StoreFrequencies("SCentralFreq.txt");
393 table[4]->StoreFrequencies("SBorderFreq.txt");
395 cout<<"Creating Tables..\n";
396 //One Huffman tree is created for each table starting from the frequencies of the symbols
397 for(Int_t i=0;i<NumTables;i++){
398 table[i]->BuildHTable();
400 cout<<"Number of elements inside the table:"<<table[i]->GetWordsNumber();
403 cout<<" (Bunch Length)"<<endl;
407 cout<<" (Time Bin)"<<endl;
411 cout<<" (1 Samples Bunch)"<<endl;
415 cout<<" (Central Samples)"<<endl;
419 cout<<" (Border Samples)"<<endl;
423 table[i]->PrintTable();
426 //The tables are saved ad binary files
427 StoreTables(table,NumTables);
428 //The tables stored in memory are deleted;
429 for(Int_t i=0;i<NumTables;i++)delete table[i];
433 /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
434 Int_t AliTPCCompression::RetrieveTables(AliTPCHTable* table[],Int_t NumTable){
435 //This method retrieve the Huffman tables from a sequence of binary files
437 cout<<"Retrieving tables from files \n";
443 //The following for loop is used to generate the Huffman trees acording to the tables
444 for(Int_t k=0;k<NumTable;k++){
445 Int_t dim;//this variable contains the table dimension
446 sprintf(filename,"Table%d.dat",k);
448 fTable.open(filename,ios::binary);
450 fTable.open(filename);
452 if(!fTable && gSystem->Getenv("ALICE_ROOT")){
454 sprintf(filename,"%s/RAW/Table%d.dat",gSystem->Getenv("ALICE_ROOT"),k);
456 fTable.open(filename,ios::binary);
458 fTable.open(filename);
462 Error("RetrieveTables", "File doesn't exist: %s", filename);
465 fTable.read((char*)(&dim),sizeof(Int_t));
467 cout<<"Table dimension: "<<dim<<endl;
468 table[k]=new AliTPCHTable(dim);
469 for(Int_t i=0;i<dim;i++){
470 fTable.read((char*)(&codeLen),sizeof(UChar_t));
471 table[k]->SetCodeLen(codeLen,i);
472 // fTable.read((char*)(&code),sizeof(UInt_t));
473 fTable.read((char*)(&code),sizeof(Double_t));
474 table[k]->SetCode(Mirror((UInt_t)code,codeLen),i);
479 cout<<"Trees generated \n";
480 //At this point the trees are been built
484 Int_t AliTPCCompression::CreateTablesFromTxtFiles(Int_t NumTable){
485 //This method creates a set of binary tables, needed by the Huffman
486 //algorith, starting from a set of frequencies tables stored in form of
489 cout<<"Retrieving frequencies from txt files \n";
492 //Tables are read from the files (Each codeword has been "Mirrored")
493 AliTPCHTable **table = new AliTPCHTable*[NumTable];
494 for(Int_t k=0;k<NumTable;k++){
495 sprintf(filename,"Table%d.txt",k);
496 cout<<filename<<endl;
497 fTable.open(filename);
499 Error("CreateTablesFromTxtFiles", "File doesn't exist: %s", filename);
504 table[k]=new AliTPCHTable(1024);
505 while(!fTable.eof()){
509 cout<<"Frequency cannot be negative !!!\n";
512 table[k]->SetValFrequency(symbol,freq);
519 fStat.open("Statistics",ios::app);
521 fStat<<"----------------- ENTROPY for external txt tables --------------------------"<<endl;
522 fStat<<"Entropy of Bunch length table......."<<table[0]->GetEntropy()<<endl;
523 fStat<<"Entropy of Time bin table..........."<<table[1]->GetEntropy()<<endl;
524 fStat<<"Entropy of one Sample bunch table..."<<table[2]->GetEntropy()<<endl;
525 fStat<<"Entropy of Central Sample table....."<<table[3]->GetEntropy()<<endl;
526 fStat<<"Entropy Border Samples table........"<<table[4]->GetEntropy()<<endl;
528 for(Int_t k=0;k<NumTable;k++){
529 table[k]->BuildHTable();
531 //The tables are saved ad binary files
532 StoreTables(table,NumTable);
533 //The tables stored in memory are deleted;
534 for(Int_t i=0;i<NumTable;i++)delete table[i];
539 ////////////////////////////////////////////////////////////////////////////////////////
541 ////////////////////////////////////////////////////////////////////////////////////////
543 void AliTPCCompression::StoreValue(UInt_t val,UChar_t len){
544 //This method stores the value "val" of "len" bits into the internal buffer "fBuffer"
545 if (len<=fFreeBitsBuffer){ // val is not splitted in two buffer
546 fFreeBitsBuffer-=len;
547 fBuffer=fBuffer<<len;
549 if(!fFreeBitsBuffer){ // if the buffer is full it is written into a file
550 f.write((char*)(&fBuffer),sizeof(UInt_t));
551 fFreeBitsBuffer=fDimBuffer;
555 else{ //val has to be splitted in two buffers
556 fBuffer=fBuffer<<fFreeBitsBuffer;
559 temp=temp>>(len-fFreeBitsBuffer);
560 fBuffer=fBuffer|temp;
561 f.write((char*)(&fBuffer),sizeof(UInt_t));
562 fFreeBitsBuffer=fDimBuffer-(len-fFreeBitsBuffer);
563 val=val<<fFreeBitsBuffer;
564 val=val>>fFreeBitsBuffer;
569 //////////////////////////////////////////////////////////////////////////////////////////////////
570 void AliTPCCompression::Flush(){
571 //The last buffer cannot be completely full so to save it
572 //into the output file it is first necessary to fill it with an hexadecimal pattern
573 if(fFreeBitsBuffer<fDimBuffer){
574 fBuffer=fBuffer<<fFreeBitsBuffer;
575 f.write((char*)(&fBuffer),sizeof(UInt_t));
579 //////////////////////////////////////////////////////////////////////////////////////////////////
580 UInt_t AliTPCCompression::Mirror(UInt_t val,UChar_t len)const{
581 //This method inverts the digits of the number "val" and length "len"
582 //indicates the number of digits of the number considered in binary notation
586 for(Int_t i=0;i<len;i++){
589 specular=specular<<1;
590 specular=specular|bit;
596 /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
597 Int_t AliTPCCompression::CompressDataOptTables(Int_t NumTable,const char* fSource,const char* fDest){
598 //This method compress an Altro format file using a general set of tables stored as binary files to be provided
600 cout<<" BackWord COMPRESSION "<<endl;
601 cout<<"compression of the file "<<fSource<<" Output File: "<<fDest<<endl;
603 //Tables are read from the files (Each codeword has been "Mirrored")
604 AliTPCHTable **table = new AliTPCHTable*[NumTable];
605 for(Int_t i=0;i<NumTable;i++) table[i] = NULL;
606 if (RetrieveTables(table,NumTable) != 0) {
607 for(Int_t i=0;i<NumTable;i++) delete table[i];
611 //the output file is open
614 f.open(fDest,ios::binary|ios::out);
616 f.open(fDest,ios::out);
618 // Source file is open
619 AliAltroBuffer buff(fSource,0);
620 //coded words are written into a file
621 Int_t numWords,padNum,rowNum,secNum=0;
622 UInt_t storedWords=0;
625 Double_t stat[5]={0.,0.,0.,0.,0.};
626 UInt_t trailerNumbers=0;
627 Double_t numElem[5]={0,0,0,0,0};
628 Double_t fillWords=0.;
630 fStat.open("Statistics",ios::app);
632 fStat<<"-------------------COMPRESSION STATISTICS----------"<<endl;
634 while(buff.ReadTrailerBackward(numWords,padNum,rowNum,secNum)){
636 fillWords=buff.GetFillWordsNum();
642 fillWords+=4-numWords%4;
643 for(Int_t j=0;j<(4-numWords%4);j++){
644 value=buff.GetNextBackWord();
651 for(Int_t i=0;i<345;i++)timePos[i]=0;
652 for(Int_t i=0;i<1024;i++)packet[i]=0;
654 Int_t nextTableType=0;
657 for(Int_t i=0;i<numWords;i++){
658 value=buff.GetNextBackWord();
660 if(nextTableType==1){
664 NextTable(value,nextTableType,bunchLen,count);
666 //computing the Time gap between two bunches
669 Int_t previousTime=packet[timePos[tp]];
670 for(Int_t i=tp-1;i>=0;i--){
671 Int_t timPos=timePos[i];
672 Int_t bunchLen=packet[timPos-1]-2;
674 packet[timPos]=packet[timPos]-previousTime-bunchLen;
682 for(Int_t i=0;i<numWords;i++){
684 if(nextTableType==1)timeBin=value;
686 //UInt_t val=(UInt_t)table[nextTableType]->Code()[value]; // val is the code
687 Double_t val=table[nextTableType]->Code()[value]; // val is the code
688 UChar_t len=table[nextTableType]->CodeLen()[value]; // len is the length (number of bits)of val
689 stat[nextTableType]+=len;
690 numElem[nextTableType]++;
691 StoreValue((UInt_t)val,len);
694 NextTable(value,nextTableType,bunchLen,count);
695 if(nextTableType==0){
696 // UInt_t val=(UInt_t)table[1]->Code()[timeBin]; // val is the code
697 Double_t val=table[1]->Code()[timeBin]; // val is the code
698 UChar_t len=table[1]->CodeLen()[timeBin]; // len is the length (number of bits)of val
701 StoreValue((UInt_t)val,len);
702 // val=(UInt_t)table[nextTableType]->Code()[(bunchLen+2)]; // val is the code
703 val=table[nextTableType]->Code()[(bunchLen+2)]; // val is the code
704 len=table[nextTableType]->CodeLen()[(bunchLen+2)]; // len is the length (number of bits)of val
705 StoreValue((UInt_t)val,len);
706 stat[nextTableType]+=len;
707 numElem[nextTableType]++;
712 StoreValue(numWords,10);
713 StoreValue(padNum,10);
714 StoreValue(rowNum,10);
715 StoreValue(secNum,9);
720 StoreValue(numPackets,32);
722 cout<<"Number of strored packets: "<<numPackets<<endl;
724 //The last buffen cannot be completely full
727 cout<<"Number of stored words: "<<storedWords<<endl;
730 for(Int_t i=0;i<NumTable;i++){
734 Double_t dimension=(UInt_t)TMath::Ceil((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])/8)+trailerNumbers*5;
735 fStat<<"Trailer Dimension in bytes......"<<trailerNumbers*5<<endl;
736 fStat<<"Data Dimension in bytes........."<<(UInt_t)TMath::Ceil((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])/8)<<endl;
737 fStat<<"Compressed file dimension......."<<(UInt_t)dimension<<endl;
739 fStat<<(UInt_t)trailerNumbers<<endl;
740 fStat<<(UInt_t)fillWords<<endl;
741 fStat<<(UInt_t)numElem[0]<<endl;
742 fStat<<(UInt_t)numElem[1]<<endl;
743 fStat<<(UInt_t)numElem[2]<<endl;
744 fStat<<(UInt_t)numElem[3]<<endl;
745 fStat<<(UInt_t)numElem[4]<<endl;
747 fillWords=(fillWords+numElem[0]+numElem[1]+numElem[2]+numElem[3]+numElem[4]+trailerNumbers*4)*10/8;
748 fStat<<"Original file dimension........."<<(UInt_t)fillWords<<endl;
750 Double_t ratio=(dimension/fillWords)*100;
751 fStat<<"Compression ratio (Compressed/Uncompressed)..."<<ratio<<"%"<<endl;
754 fStat<<"Bunch length size in bytes......"<<(UInt_t)TMath::Ceil(stat[0]/8)<<" Comppression.."<<(stat[0]/numElem[0])*10<<"%"<<endl;
756 fStat<<"Time gap size in bytes.........."<<(UInt_t)TMath::Ceil(stat[1]/8)<<" Comppression.."<<(stat[1]/numElem[1])*10<<"%"<<endl;
757 if (numElem[2]+numElem[3]+numElem[4])
758 fStat<<"Amplitude values in bytes......."<<(UInt_t)TMath::Ceil((stat[2]+stat[3]+stat[4])/8)<<" Comppression.."<<
759 ((stat[2]+stat[3]+stat[4])/(numElem[2]+numElem[3]+numElem[4]))*10<<"%"<<endl;
761 fStat<<" One Samples in bytes............"<<(UInt_t)TMath::Ceil(stat[2]/8)<<" Comppression.."<<(stat[2]/numElem[2])*10<<"%"<<endl;
763 fStat<<" Central Samples size in bytes..."<<(UInt_t)TMath::Ceil(stat[3]/8)<<" Comppression.."<<(stat[3]/numElem[3])*10<<"%"<<endl;
765 fStat<<" Border Samples size in bytes...."<<(UInt_t)TMath::Ceil(stat[4]/8)<<" Comppression.."<<(stat[4]/numElem[4])*10<<"%"<<endl;
767 fStat<<"Average number of bits per word"<<endl;
769 fStat<<"Bunch length ......"<<stat[0]/numElem[0]<<endl;
771 fStat<<"Time gap .........."<<stat[1]/numElem[1]<<endl;
773 fStat<<"One Samples........"<<stat[2]/numElem[2]<<endl;
775 fStat<<"Central Samples ..."<<stat[3]/numElem[3]<<endl;
777 fStat<<"Border Samples....."<<stat[4]/numElem[4]<<endl;
782 ////////////////////////////////////////////////////////////////////////////////////////
784 ////////////////////////////////////////////////////////////////////////////////////////
786 ////////////////////////////////////////////////////////////////////////////////////////
788 Int_t AliTPCCompression::CreateTreesFromFile(AliTPCHNode *RootNode[],Int_t NumTables){
789 //For each table this method builds the associate Huffman tree starting from the codeword and
790 //the codelength of each symbol
792 cout<<"Creating the Huffman trees \n";
799 //The following for loop is used to generate the Huffman trees acording to the tables
800 //loop over the tables
801 for(Int_t k=0;k<NumTables;k++){
802 RootNode[k]=new AliTPCHNode(); //RootNode is the root of the tree
803 Int_t dim=0;//this variable contains the table dimension
804 sprintf(filename,"Table%d.dat",k);
806 fTable.open(filename,ios::binary);
808 fTable.open(filename);
810 if(!fTable && gSystem->Getenv("ALICE_ROOT")){
812 sprintf(filename,"%s/RAW/Table%d.dat",gSystem->Getenv("ALICE_ROOT"),k);
814 fTable.open(filename,ios::binary);
816 fTable.open(filename);
820 Error("CreateTreesFromFile", "File doesn't exist: %s", filename);
823 fTable.read((char*)(&dim),sizeof(Int_t));
825 cout<<"Table dimension: "<<dim<<endl;
826 //loop over the words of one table
827 for(Int_t i=0;i<dim;i++){
828 fTable.read((char*)(&codeLen),sizeof(UChar_t));
829 //fTable.read((char*)(&code),sizeof(UInt_t));
830 fTable.read((char*)(&code),sizeof(Double_t));
832 for(Int_t j=1;j<=codeLen;j++){
834 // val=(UInt_t)TMath::Power(2,codeLen-j);
835 val=(UInt_t)(1<<(codeLen-j));
836 bit=(UInt_t)code&val;
837 AliTPCHNode *temp=node;
839 node=node->GetRight();
841 node=new AliTPCHNode();
842 temp->SetRight(node);
846 node=node->GetLeft();
848 node=new AliTPCHNode();
855 node->SetFrequency(codeLen);
861 cout<<"Trees generated \n";
862 //At this point the trees are been built
866 Int_t AliTPCCompression::CreateLUTsFromTrees(AliTPCHNode *RootNode[],Int_t NumTables,UInt_t *LUTDimension[],AliTPCHNode **LUTNode[]){
867 //For each Huffman tree this method builds the associate look-up-table for fast
868 //decoding of compressed data.
869 //Should be called after CreateTreesFromFile.
871 cout<<"Creating the LUTs \n";
874 //loop over the tables
875 for(Int_t k=0;k<NumTables;k++){
876 UInt_t dim = 1<<fgkTableDimension;
877 LUTDimension[k] = new UInt_t[dim];
878 LUTNode[k] = new AliTPCHNode*[dim];
879 //loop over the words of one LUT
880 for(UInt_t i=0;i<dim;i++){
883 LUTDimension[k][i] = GetDecodedLUTBuffer(&node,&buffer);
884 LUTNode[k][i] = node;
885 // cout<<"LUT "<<k<<" "<<i<<" "<<LUTDimension[k][i]<<" "<<LUTNode[k][i]<<" "<<LUTNode[k][i]->GetSymbol()<<endl;
889 cout<<"LUTs generated \n";
890 //At this point the LUTs are built
893 //////////////////////////////////////////////////////////////////////////////////////////////////
894 void AliTPCCompression::DeleteHuffmanTree(AliTPCHNode* node){
895 //This function deletes all the nodes of an Huffman tree
896 //In an Huffman tree any internal node has always two children
898 DeleteHuffmanTree(node->GetLeft());
899 DeleteHuffmanTree(node->GetRight());
900 // cout<<node->GetSymbol()<<" "<<(Int_t)node->GetFrequency()<<endl;
904 //////////////////////////////////////////////////////////////////////////////////////////////////
905 void AliTPCCompression::VisitHuffmanTree(AliTPCHNode* node){
906 //This function realizes an in order visit of a binary tree
908 cout<<node->GetSymbol()<<" "<<node->GetFrequency()<<endl;
909 VisitHuffmanTree(node->GetLeft());
910 VisitHuffmanTree(node->GetRight());
913 //////////////////////////////////////////////////////////////////////////////////////////////////
914 UInt_t AliTPCCompression::ReadWord(UInt_t NumberOfBit){
915 //This method retrieves a word of a specific number of bits from the file through the internal buffer
918 for (UInt_t i=0;i<NumberOfBit;i++){
920 fPos-=sizeof(UInt_t);
922 f.read((char*)(&fBuffer),sizeof(UInt_t));
925 UInt_t mask=(UInt_t)(1<<(32-fBitsToRead));
927 bit=bit>>(32-fBitsToRead);
934 //////////////////////////////////////////////////////////////////////////////////////////////////
935 UInt_t AliTPCCompression::ReadWordBuffer(UInt_t NumberOfBit){
936 //This method retrieves a word of a specific number of bits from the file through the buffer
939 if(NumberOfBit<=fBitsToRead) {
940 result = fBuffer&((1<<NumberOfBit)-1);
941 fBuffer=fBuffer>>NumberOfBit;
942 fBitsToRead-=NumberOfBit;
946 UInt_t tempbuffer = *fPointBuffer;
947 if((NumberOfBit-fBitsToRead) != 32)
948 tempbuffer=tempbuffer&((1<<(NumberOfBit-fBitsToRead))-1);
949 tempbuffer=tempbuffer<<fBitsToRead;
950 result = fBuffer|tempbuffer;
951 fBuffer=*fPointBuffer;
952 fBitsToRead=(32+fBitsToRead)-NumberOfBit;
953 fBuffer=fBuffer>>(32-fBitsToRead);
959 //////////////////////////////////////////////////////////////////////////////////////////////////
960 inline void AliTPCCompression::AdjustWordBuffer(UInt_t NumberOfBit){
961 //This method retrieves a word of a specific number of bits from the file through the buffer
962 //The method is used together with LUTs for fast decoding
963 if(NumberOfBit<=fBitsToRead) {
964 fBuffer=fBuffer>>NumberOfBit;
965 fBitsToRead-=NumberOfBit;
969 fBuffer=*fPointBuffer;
970 fBitsToRead=(32+fBitsToRead)-NumberOfBit;
971 fBuffer=fBuffer>>(32-fBitsToRead);
975 //////////////////////////////////////////////////////////////////////////////////////////////////
976 inline UInt_t AliTPCCompression::ReadWordBufferWithLUTs(){
977 //This method retrieves a word of a specific number of bits from the file through the buffer
978 //The method is used together with LUTs for fast decoding
979 if(fgkTableDimension<=fBitsToRead)
980 return fBuffer&((1<<fgkTableDimension)-1);
982 UInt_t tempbuffer = *(fPointBuffer-1);
983 tempbuffer=tempbuffer&((1<<(fgkTableDimension-fBitsToRead))-1);
984 tempbuffer=tempbuffer<<fBitsToRead;
985 return fBuffer|tempbuffer;
989 //////////////////////////////////////////////////////////////////////////////////////////////////
990 inline UInt_t AliTPCCompression::ReadBitFromWordBuffer(){
991 //This method retrieves a bit from the file through the buffer
996 fBuffer=*fPointBuffer;
1005 //////////////////////////////////////////////////////////////////////////////////////////////////
1006 inline UInt_t AliTPCCompression::ReadBitFromLUTBuffer(UInt_t *buffer){
1007 //This method retrieves a word of a bit out of the LUT buffer
1015 //////////////////////////////////////////////////////////////////////////////////////////////////
1016 void AliTPCCompression::ReadTrailer(Int_t &WordsNumber,Int_t &PadNumber,Int_t &RowNumber,Int_t &SecNumber,Bool_t Memory){
1017 //It retrieves a trailer
1019 ReadBitFromWordBuffer();
1020 SecNumber=ReadWordBuffer(9);
1021 RowNumber=ReadWordBuffer(10);
1022 PadNumber=ReadWordBuffer(10);
1023 WordsNumber=ReadWordBuffer(10);
1027 SecNumber=ReadWord(9);
1028 RowNumber=ReadWord(10);
1029 PadNumber=ReadWord(10);
1030 WordsNumber=ReadWord(10);
1034 //////////////////////////////////////////////////////////////////////////////////////////////////
1035 inline UInt_t AliTPCCompression::GetDecodedWordBuffer(AliTPCHNode* root){
1036 //This method retrieves a decoded word.
1037 AliTPCHNode *node=root;
1042 UInt_t bit=ReadBitFromWordBuffer();
1044 node=node->GetRight();
1046 node=node->GetLeft();
1047 if (!(node->GetLeft())){
1048 symbol=node->GetSymbol();
1055 //////////////////////////////////////////////////////////////////////////////////////////////////
1056 inline UInt_t AliTPCCompression::GetDecodedWordBufferWithLUTs(UInt_t* LUTDimension,AliTPCHNode** LUTNode){
1057 //This method retrieves a decoded word.
1062 buffer = ReadWordBufferWithLUTs();
1063 AliTPCHNode *node=LUTNode[buffer];
1064 if (!(node->GetLeft())){
1065 symbol=node->GetSymbol();
1068 AdjustWordBuffer(LUTDimension[buffer]);
1070 UInt_t bit=ReadBitFromWordBuffer();
1072 node=node->GetRight();
1074 node=node->GetLeft();
1075 if (!(node->GetLeft())){
1076 symbol=node->GetSymbol();
1083 //////////////////////////////////////////////////////////////////////////////////////////////////
1084 inline UInt_t AliTPCCompression::GetDecodedLUTBuffer(AliTPCHNode** node,UInt_t *buffer){
1085 //This method retrieves a decoded word for the LUTs.
1089 while(!decoded && counter<fgkTableDimension){
1090 UInt_t bit=ReadBitFromLUTBuffer(buffer);
1093 *node=(*node)->GetRight();
1095 *node=(*node)->GetLeft();
1096 if (!((*node)->GetLeft())){
1103 inline UInt_t AliTPCCompression::GetDecodedWord(AliTPCHNode* root){
1104 //This method retrieves a decoded word.
1105 AliTPCHNode *node=root;
1110 UInt_t bit=ReadWord(1);
1112 node=node->GetRight();
1114 node=node->GetLeft();
1115 if (!(node->GetLeft())){
1116 symbol=node->GetSymbol();
1123 //////////////////////////////////////////////////////////////////////////////////////////////////
1125 Int_t AliTPCCompression::DecompressDataOptTables(Int_t NumTables,const char* fname, const char* fDest){
1126 //This method decompress a file using separate Huffman tables
1128 cout<<" DECOMPRESSION:"<<endl;
1129 cout<<"Source File "<<fname<<" Destination File "<<fDest<<endl;
1131 AliTPCHNode ** rootNode = new AliTPCHNode*[NumTables];
1132 for(Int_t i=0;i<NumTables;i++) rootNode[i] = NULL;
1133 //Creation of the Huffman trees
1134 if (CreateTreesFromFile(rootNode,NumTables) != 0) {
1135 for(Int_t i=0;i<NumTables;i++) {
1136 if (rootNode[i]) DeleteHuffmanTree(rootNode[i]);
1143 f.open(fname,ios::binary|ios::in);
1145 f.open(fname,ios::in);
1148 Error("DecompressDataOptTables", "File doesn't exist:",fname);
1151 //to go to the end of the file
1152 f.seekg(0,ios::end);
1153 //to get the file dimension in byte
1155 fPos-=sizeof(UInt_t);
1159 f.read((char*)(&fBuffer),sizeof(UInt_t));
1167 UInt_t packetNumber=ReadWord(sizeof(UInt_t)*8);
1169 cout<<"Number of Packect: "<<packetNumber<<endl;
1171 AliAltroBuffer bufferFile(fDest,1);
1173 UInt_t wordsRead=0; //number of read coded words
1174 while(k<packetNumber){
1175 Int_t numWords,padNumber,rowNumber,secNumber=0;
1176 ReadTrailer(numWords,padNumber,rowNumber,secNumber,kFALSE);
1179 Int_t previousTime=-1;
1181 Int_t nextTableType=0;
1184 for(Int_t i=0;i<numWords;i++){
1185 UInt_t symbol=GetDecodedWord(rootNode[nextTableType]);
1187 //Time reconstruction
1188 if (nextTableType==1){
1189 if (previousTime!=-1){
1190 previousTime=symbol+previousTime+bunchLen;
1192 else previousTime=symbol;
1196 bufferFile.FillBuffer(symbol);
1197 NextTable(symbol,nextTableType,bunchLen,count);
1198 if(nextTableType==0){
1199 bufferFile.FillBuffer(time);
1200 bufferFile.FillBuffer(bunchLen+2);
1204 bufferFile.WriteTrailer(numWords,padNumber,rowNumber,secNumber);
1207 cout<<"Number of decoded words:"<<wordsRead<<endl;
1210 //The trees are deleted
1211 for(Int_t j=0;j<NumTables;j++){
1212 DeleteHuffmanTree(rootNode[j]);
1218 //////////////////////////////////////////////////////////////////////////////////////////////////
1219 Int_t AliTPCCompression::Decompress(AliTPCHNode *RootNode[],Int_t /*NumTables*/,char* PointBuffer,UInt_t BufferSize,UShort_t out[],UInt_t &dim){
1220 //This method decompress a file using separate Huffman tables
1222 // fPointBuffer=((UInt_t *)PointBuffer)+(UInt_t)(BufferSize/4)-1;
1223 fPointBuffer=(UInt_t *)(PointBuffer+BufferSize-4);
1227 fBuffer=*fPointBuffer;
1234 UInt_t packetNumber=ReadWordBuffer(sizeof(UInt_t)*8); //32 bits
1236 cout<<"First one has been found "<<endl;
1237 cout<<"Number of packets:"<<packetNumber<<endl;
1240 UInt_t wordsRead=0; //number of read coded words
1241 while(k<packetNumber){
1242 Int_t numWords,padNumber,rowNumber,secNumber=0;
1243 ReadTrailer(numWords,padNumber,rowNumber,secNumber,kTRUE);
1252 //ftxt<<"S:"<<secNumber<<" R:"<<rowNumber<<" P:"<<padNumber<<" W:"<<numWords<<endl;
1253 // padDigits->SetPadID(padNumber,rowNumber,secNumber,DDL);
1256 Int_t previousTime=-1;
1258 Int_t nextTableType=0;
1262 for(Int_t i=0;i<numWords;i++){
1263 UInt_t symbol=GetDecodedWordBuffer(RootNode[nextTableType]);
1265 //Time reconstruction
1266 if (nextTableType==1){
1267 if (previousTime!=-1){
1268 previousTime=symbol+previousTime+bunchLen;
1270 else previousTime=symbol;
1272 out[dim]=bunchLen+2;
1276 timeDigit=time-bunchLen;
1278 if(nextTableType>1){
1280 //ftxt<<symbol<<endl;
1284 //padDigits->SetDigits(symbol,timeDigit);
1286 NextTable(symbol,nextTableType,bunchLen,count);
1287 if(nextTableType==0){
1290 // ftxt<<(bunchLen+2)<<endl;
1298 //////////////////////////////////////////////////////////////////////////////////////////////////
1299 Int_t AliTPCCompression::DecompressWithLUTs(AliTPCHNode *RootNode[],UInt_t *LUTDimension[],AliTPCHNode **LUTNode[],Int_t /*NumTables*/,char* PointBuffer,UInt_t BufferSize,UShort_t out[],UInt_t &dim){
1300 //This method decompress a file using separate Huffman tables and already prepared LUTs for fast decoding
1302 // fPointBuffer=((UInt_t *)PointBuffer)+(UInt_t)(BufferSize/4)-1;
1303 fPointBuffer=(UInt_t *)(PointBuffer+BufferSize-4);
1307 fBuffer=*fPointBuffer;
1314 UInt_t packetNumber=ReadWordBuffer(sizeof(UInt_t)*8); //32 bits
1316 cout<<"First one has been found "<<endl;
1317 cout<<"Number of packets:"<<packetNumber<<endl;
1320 UInt_t wordsRead=0; //number of read coded words
1321 while(k<packetNumber){
1322 Int_t numWords,padNumber,rowNumber,secNumber=0;
1323 ReadTrailer(numWords,padNumber,rowNumber,secNumber,kTRUE);
1332 //ftxt<<"S:"<<secNumber<<" R:"<<rowNumber<<" P:"<<padNumber<<" W:"<<numWords<<endl;
1333 // padDigits->SetPadID(padNumber,rowNumber,secNumber,DDL);
1336 Int_t previousTime=-1;
1338 Int_t nextTableType=0;
1342 for(Int_t i=0;i<numWords;i++){
1344 if(i == (numWords-1))
1345 symbol = GetDecodedWordBuffer(RootNode[nextTableType]);
1347 symbol=GetDecodedWordBufferWithLUTs(LUTDimension[nextTableType],LUTNode[nextTableType]);
1349 //Time reconstruction
1350 if (nextTableType==1){
1351 if (previousTime!=-1){
1352 previousTime=symbol+previousTime+bunchLen;
1354 else previousTime=symbol;
1356 out[dim]=bunchLen+2;
1360 timeDigit=time-bunchLen;
1362 if(nextTableType>1){
1364 //ftxt<<symbol<<endl;
1368 //padDigits->SetDigits(symbol,timeDigit);
1370 NextTable(symbol,nextTableType,bunchLen,count);
1371 if(nextTableType==0){
1374 // ftxt<<(bunchLen+2)<<endl;
1382 //////////////////////////////////////////////////////////////////////////////////////////////////
1383 Int_t AliTPCCompression::DestroyTables(AliTPCHNode *RootNode[],Int_t NumTables){
1384 //The trees are deleted
1385 for(Int_t j=0;j<NumTables;j++){
1386 DeleteHuffmanTree(RootNode[j]);
1389 cout<<"Huffman trees destroyed"<<endl;
1392 //////////////////////////////////////////////////////////////////////////////////////////////////
1394 void AliTPCCompression::ReadAltroFormat(char* fileOut,char* fileIn)const{
1395 //This method creates a text file containing the same information stored in
1396 //an Altro file. The information in the text file is organized pad by pad and
1397 //and for each pad it consists in a sequence of bunches (Bunch length +2,
1398 //Time bin of the last amplitude sample in the bunch, amplitude values)
1399 //It is used mainly for debugging
1400 ofstream ftxt(fileOut);
1401 AliAltroBuffer buff(fileIn,0);
1402 Int_t numWords,padNum,rowNum,secNum=0;
1404 if (fVerbose) cout<<"Creating a txt file from an Altro Format file"<<endl;
1405 while(buff.ReadTrailerBackward(numWords,padNum,rowNum,secNum)){
1406 ftxt<<"S:"<<secNum<<" R:"<<rowNum<<" P:"<<padNum<<" W:"<<numWords<<endl;
1408 for(Int_t j=0;j<(4-numWords%4);j++){
1409 value=buff.GetNextBackWord();
1412 for(Int_t i=0;i<numWords;i++){
1413 value=buff.GetNextBackWord();
1421 //////////////////////////////////////////////////////////////////////////////////////////