]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCompression.cxx
Removing unused include
[u/mrichter/AliRoot.git] / TPC / AliTPCCompression.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18
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.
24
25 #include <stdlib.h>
26 #include <TObjArray.h>
27 #include <Riostream.h>
28 #include <TMath.h>
29 #include "AliTPCBuffer160.h"
30 #include "AliTPCHuffman.h"
31 #include "AliTPCCompression.h"
32
33 ClassImp(AliTPCCompression)
34 //////////////////////////////////////////////////////////////////////////////////////////////////
35 AliTPCCompression::AliTPCCompression(){
36   //Defaul constructor
37   fDimBuffer=sizeof(ULong_t)*8;
38   fFreeBitsBuffer=fDimBuffer;
39   fReadBits=0;
40   fPos=0;
41   fBuffer=0;
42   fVerbose=0;
43   fFillWords=0;
44   fPointBuffer=0;
45   return;
46 }
47 //////////////////////////////////////////////////////////////////////////////////////////////////
48 AliTPCCompression::AliTPCCompression(const AliTPCCompression &source){
49   //Constructor
50   this->fDimBuffer=source.fDimBuffer;
51   this->fFreeBitsBuffer=source.fFreeBitsBuffer;
52   this->fReadBits=source.fReadBits;
53   this->fPos=source.fPos;
54   this->fBuffer=source.fBuffer;
55   this->fVerbose=source.fVerbose;
56   this->fFillWords=source.fFillWords;
57   this->fPointBuffer=source.fPointBuffer;
58   return;
59 }
60 //////////////////////////////////////////////////////////////////////////////////////////////////
61 AliTPCCompression&  AliTPCCompression::operator=(const AliTPCCompression &source){
62   //Redefinition of the assignment operator
63   this->fDimBuffer=source.fDimBuffer;
64   this->fFreeBitsBuffer=source.fFreeBitsBuffer;
65   this->fReadBits=source.fReadBits;
66   this->fPos=source.fPos;
67   this->fBuffer=source.fBuffer;
68   this->fVerbose=source.fVerbose;
69   this->fFillWords=source.fFillWords;
70   this->fPointBuffer=source.fPointBuffer;
71   return *this;
72
73 //////////////////////////////////////////////////////////////////////////////////////////////////
74 void AliTPCCompression::NextTable(Int_t Val,Int_t &NextTableType,Int_t &BunchLen,Int_t &Count)const{
75   //Depending on the data type (5 types of data) a specific table is called
76   /*
77     Table index:
78     0==> Bunch length value     
79     1==> Time Bin value 
80     2==> 1-samples bunch
81     3==> Central samples
82     4==> Border samples
83   */  
84   switch (NextTableType){
85   case 0:{
86     BunchLen=Val-2;
87     NextTableType=1;
88     break;
89   }//end case 0
90   case 1:{
91     if (BunchLen==1)NextTableType=2;
92     else{
93       NextTableType=4;
94       Count=1;
95     }
96     break;
97   }//end case 1
98   case 2:{
99     NextTableType=0;
100     break;
101   }//end case 2
102   case 3:{
103     Count++;
104     if (Count==(BunchLen-1)){
105       NextTableType=4;
106     }
107     break;
108   }//end case 3
109   case 4:{
110     if (Count==1){
111       if (BunchLen>2)
112         NextTableType=3;
113       else
114         Count++;
115     }
116     else
117       NextTableType=0;
118     break;
119   }//end case 4
120   }//end switch
121   return;
122 }
123
124
125 /////////////////////////////////////////////////////////////////////////////////////////////////////
126
127 Int_t AliTPCCompression::FillTables(const char* fSource,AliTPCHTable* table[],const Int_t NumTables){
128   //This method is used to compute the frequencies of the symbols in the source file
129   AliTPCBuffer160 buff(fSource,0);
130   ULong_t countWords=0;
131   ULong_t countTrailer=0;
132   Int_t numWords,padNum,rowNum,secNum=0;
133   Int_t value=0;
134   ULong_t stat[5]={0,0,0,0,0};
135   Int_t endFill=0;
136   Int_t end=1;
137   while(buff.ReadTrailerBackward(numWords,padNum,rowNum,secNum) !=-1 ){
138     if(end){
139       endFill=buff.GetFillWordsNum();
140       end=0;
141     }//endif
142     countTrailer++;
143     if (numWords%4){
144       fFillWords+=4-numWords%4;
145       for(Int_t j=0;j<(4-numWords%4);j++){
146         value=buff.GetNextBackWord();
147       }//end for
148     }//end if
149     
150     Int_t packet[1024];
151     Int_t timePos[345];
152     Int_t tp=0;
153     for(Int_t i=0;i<345;i++)timePos[i]=0;
154     for(Int_t i=0;i<1024;i++)packet[i]=0;
155     
156     Int_t nextTableType=0;
157     Int_t bunchLen=0;
158     Int_t count=0;
159     for(Int_t i=0;i<numWords;i++){
160       value=buff.GetNextBackWord();
161       packet[i]=value;
162       if(nextTableType==1){
163         timePos[tp]=i;
164         tp++;
165       }
166       NextTable(value,nextTableType,bunchLen,count);
167     }//end for
168     //computing the Time gap between two bunches
169     Int_t temp=0;
170     tp--;
171     Int_t previousTime=packet[timePos[tp]];
172     for(Int_t i=tp-1;i>=0;i--){
173       Int_t timPos=timePos[i];
174       Int_t bunchLen=packet[timPos-1]-2;
175       temp=packet[timPos];
176       packet[timPos]=packet[timPos]-previousTime-bunchLen;
177       previousTime=temp;
178     }
179     nextTableType=0;
180     count=0;
181     bunchLen=0;
182     for(Int_t i=0;i<numWords;i++){
183       value=packet[i];
184       table[nextTableType]->SetFrequency(value);
185       stat[nextTableType]++;
186       NextTable(value,nextTableType,bunchLen,count);
187       countWords++;
188     }//end for
189   }//end while
190   cout<<"Number of words:       "<<countWords<<endl;
191   cout<<"Number of trailers:    "<<countTrailer<<endl;
192   cout<<"Number of fill words   "<<fFillWords+endFill<<endl;
193   cout<<"Total number of words: "<<countWords+countTrailer*4+fFillWords<<endl;
194   //STATISTICS  
195   fStat.open("Statistics");
196   fStat<<"Number of words:..........................................."<<countWords<<endl;
197   fStat<<"Number of trailers (4 10 bits words in each one)..........."<<countTrailer<<endl;
198   fStat<<"Number of fill words:......................................"<<fFillWords+endFill<<endl;
199   fStat<<"Total number of words:....................................."<<countWords+countTrailer*4+fFillWords+endFill<<endl;
200   fStat<<"-----------------------------------------"<<endl;
201   fStat<<"Number of Bunches............."<<stat[0]<<endl;
202   fStat<<"Number of Time bin............"<<stat[1]<<endl;
203   fStat<<"Number of One Samples Bunch..."<<stat[2]<<endl;
204   fStat<<"Number of Central Samples....."<<stat[3]<<endl;
205   fStat<<"Number of Border Samples......"<<stat[4]<<endl;
206   fStat<<"-----------------------------------------"<<endl;
207   ULong_t fileDimension=(ULong_t)TMath::Ceil(double((countTrailer*4+countWords+fFillWords+endFill)*10/8));
208   fStat<<"Total file Size in bytes.."<<fileDimension<<endl;
209   Double_t percentage=TMath::Ceil((fFillWords+endFill)*125)/fileDimension;
210   fStat<<"Fill Words................"<<(ULong_t)TMath::Ceil((fFillWords+endFill)*10/8)<<" bytes   "<<percentage<<"%"<<endl;  
211   percentage=(Double_t)countTrailer*500/fileDimension;
212   fStat<<"Trailer..................."<<countTrailer*5<<" bytes   "<<percentage<<"%"<<endl;
213
214   percentage=(Double_t)((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])) *125/fileDimension;
215   fStat<<"Data......................"<<(ULong_t)TMath::Ceil((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])*10/8)<<" bytes   "<<percentage<<"%"<<endl;
216
217   percentage=(Double_t)(stat[0]*125)/fileDimension;
218   fStat<<"Bunch....................."<<(ULong_t)TMath::Ceil(stat[0]*10/8)<<" bytes  "<<percentage<<"%"<<endl;  //  
219   percentage=(Double_t)(stat[1]*125)/fileDimension;
220   fStat<<"Time......................"<<(ULong_t)TMath::Ceil(stat[1]*10/8)<<" bytes  "<<percentage<<"%"<<endl;  //  
221
222
223   percentage=(Double_t)((stat[2]+stat[3]+stat[4])) *125/fileDimension;
224   fStat<<"Amplitude values.........."<<(ULong_t)TMath::Ceil((stat[2]+stat[3]+stat[4])*10/8)<<" bytes  "<<percentage<<"%"<<endl;
225   percentage=(Double_t)(stat[2]*125)/fileDimension;
226   fStat<<"     One Samples..............."<<(ULong_t)TMath::Ceil(stat[2]*10/8)<<" bytes  "<<percentage<<"%"<<endl;  //  
227   percentage=(Double_t)(stat[3]*125)/fileDimension;
228   fStat<<"     Central Samples..........."<<(ULong_t)TMath::Ceil(stat[3]*10/8)<<" bytes  "<<percentage<<"%"<<endl;  //  
229   percentage=(Double_t)(stat[4]*125)/fileDimension;
230   fStat<<"     Border Samples............"<<(ULong_t)TMath::Ceil(stat[4]*10/8)<<" bytes  "<<percentage<<"%"<<endl;  //  
231   fStat.close();
232   return 0;
233 }
234 ////////////////////////////////////////////////////////////////////////////////////////
235 Int_t AliTPCCompression::StoreTables(AliTPCHTable* table[],const Int_t NumTable){
236   //This method stores the tables in a sequence of binary file
237   char filename[15];
238   ofstream fTable;
239   for(Int_t k=0;k<NumTable;k++){
240     sprintf(filename,"Table%d.dat",k);
241 #ifndef __DECCXX 
242     fTable.open(filename,ios::binary);
243 #else
244     fTable.open(filename);
245 #endif
246     Int_t dim=table[k]->Size();
247     //Table dimension is written into a file
248     fTable.write((char*)(&dim),sizeof(Int_t));
249     //One table is written into a file
250     for(Int_t i=0;i<dim;i++){
251       UChar_t codeLen=table[k]->CodeLen()[i];
252       //      ULong_t code=(ULong_t)table[k]->Code()[i];
253       Double_t code=table[k]->Code()[i];
254       fTable.write((char*)(&codeLen),sizeof(UChar_t));
255       //fTable.write((char*)(&code),sizeof(ULong_t));
256       fTable.write((char*)(&code),sizeof(Double_t));
257     } //end for
258     fTable.close();
259   }//end for
260   return 0;
261 }
262 ////////////////////////////////////////////////////////////////////////////////////////
263 Int_t AliTPCCompression::CreateTableFormula(Double_t beta,ULong_t  M,Int_t dim,Int_t Type){
264   // Type = 0 for Bunch length
265   // Type = 1 for Time Gap
266   ULong_t freq;
267   Double_t sum=0;
268   Double_t min=10;
269   Double_t alpha=0;
270   Double_t A=0;
271   AliTPCHTable *Table=new AliTPCHTable(dim);
272   
273   freq=1;
274   Double_t FreqArray[1024];
275   for(Int_t i=0;i<1024;i++){
276     FreqArray[i]=0;
277   }
278   alpha=M*0.000000602+0.0104;
279   if (fVerbose)
280     cout<<"alpha "<<alpha<<endl;
281   for(Int_t x=0;x<dim;x++){
282     if (Type==1)
283       FreqArray[x]=TMath::Power((x+1),-beta)*TMath::Exp(-alpha*(x+1));
284     else
285       FreqArray[x]=TMath::Power((x+1),-beta);
286     sum+=FreqArray[x];
287     if (FreqArray[x]<min)min=FreqArray[x];
288   }//end for
289   if (fVerbose)
290     cout<<"Minimun Value "<<min<<endl;
291   A=1/sum;
292   if (fVerbose)
293     cout<<"A Value: "<<A<<endl;
294   for(Int_t x=0;x<dim;x++){
295     if (Type==0)//Bunch length
296       if (x>=3)//minimum bunch length
297         Table->SetValFrequency(x,A*FreqArray[x]*1000);
298       else
299         Table->SetValFrequency(x,0);
300     else //Time table
301       Table->SetValFrequency(x,A*FreqArray[x]);
302   }
303   Table->BuildHTable();
304   ofstream fTable;
305   char filename[15];
306   sprintf(filename,"Table%d.dat",Type); 
307 #ifndef __DECCXX 
308   fTable.open(filename,ios::binary);
309 #else
310   fTable.open(filename);
311 #endif
312   Int_t dimTable=Table->Size();
313   //Table dimension is written into a file
314   fTable.write((char*)(&dimTable),sizeof(Int_t));
315   //One table is written into a file
316   for(Int_t i=0;i<dimTable;i++){
317     UChar_t CodeLen=Table->CodeLen()[i];
318     Double_t Code=Table->Code()[i];
319     fTable.write((char*)(&CodeLen),sizeof(UChar_t));
320     fTable.write((char*)(&Code),sizeof(Double_t));
321   } //end for
322   fTable.close();
323   delete Table;
324   return 0;
325 }
326 ////////////////////////////////////////////////////////////////////////////////////////
327 Int_t AliTPCCompression::CreateTables(const char* fSource,const Int_t NumTables){
328   //Tables manager
329   /*
330     Table index:
331     0==> Bunch length values     
332     1==> Time Bin values 
333     2==> 1-samples bunch
334     3==> Central samples
335     4==> Border samples
336   */
337   Int_t n=10;// 10 bits per symbol 
338   AliTPCHTable ** table = new AliTPCHTable*[NumTables];
339   //The table is inizialized with the rigth number of rows 
340   for(Int_t i=0;i<NumTables;i++){
341     table[i]=new  AliTPCHTable((Int_t)(TMath::Power(2,n)));
342     table[i]->SetVerbose(fVerbose);
343   }
344   //The frequencies are calculated and the tables are filled
345   if (fVerbose)
346     cout<<"Filling tables...\n";
347   //The get the frequencies 
348   FillTables(fSource,table,NumTables);
349
350   //This part will be used in the table optimization phase
351   
352   for(Int_t i=0;i<NumTables;i++){
353     table[i]->CompleteTable(i);
354   }
355   
356   if(fVerbose){
357     cout<<"Entropy of Bunch length table........."<<table[0]->GetEntropy()<<endl;
358     cout<<"Entropy of Time bin table............."<<table[1]->GetEntropy()<<endl;
359     cout<<"Entropy of one Sample bunch table....."<<table[2]->GetEntropy()<<endl;
360     cout<<"Entropy of Central Sample table......."<<table[3]->GetEntropy()<<endl;
361     cout<<"Entropy Border Samples table.........."<<table[4]->GetEntropy()<<endl;
362   }
363   fStat.open("Statistics",ios::app);
364   fStat<<endl;
365   fStat<<"----------------- ENTROPY for castomized tables --------------------------"<<endl;
366   fStat<<"Entropy of Bunch length table......."<<table[0]->GetEntropy()<<endl;
367   fStat<<"Entropy of Time bin table..........."<<table[1]->GetEntropy()<<endl;
368   fStat<<"Entropy of one Sample bunch table..."<<table[2]->GetEntropy()<<endl;
369   fStat<<"Entropy of Central Sample table....."<<table[3]->GetEntropy()<<endl;
370   fStat<<"Entropy Border Samples table........"<<table[4]->GetEntropy()<<endl;
371   fStat.close();
372  
373   if (fVerbose)
374     cout<<"Tables filled \n";
375   
376   //Frequencies normalization
377   table[0]->NormalizeFrequencies();
378   table[1]->NormalizeFrequencies();
379   table[2]->NormalizeFrequencies();
380   table[3]->NormalizeFrequencies();
381   table[4]->NormalizeFrequencies();
382   
383   //Tables are saved in a sequence of text file and using the macro Histo.C is it possible to get
384   //a series of histograms rappresenting the frequency distribution
385   table[0]->StoreFrequencies("BunchLenFreq.txt");
386   table[1]->StoreFrequencies("TimeFreq.txt");
387   table[2]->StoreFrequencies("Sample1Freq.txt");
388   table[3]->StoreFrequencies("SCentralFreq.txt");
389   table[4]->StoreFrequencies("SBorderFreq.txt");
390   if (fVerbose)
391     cout<<"Creating Tables..\n";
392   //One Huffman tree is created for each table starting from the frequencies of the symbols
393   for(Int_t i=0;i<NumTables;i++){
394     table[i]->BuildHTable();
395     if (fVerbose==2){
396       cout<<"Number of elements inside the table:"<<table[i]->GetWordsNumber();
397       switch(i){
398       case 0:{
399         cout<<" (Bunch Length)"<<endl;
400         break;
401       }
402       case 1:{
403         cout<<" (Time Bin)"<<endl;
404         break;
405       }
406       case 2:{
407         cout<<" (1 Samples Bunch)"<<endl;
408         break;
409       }
410       case 3:{
411         cout<<" (Central Samples)"<<endl;
412         break;
413       }
414       case 4:{
415         cout<<" (Border Samples)"<<endl;
416         break;
417       }
418       }//end switch
419       table[i]->PrintTable();
420     }
421   }
422   //The tables are saved ad binary files
423   StoreTables(table,NumTables);
424   //The tables stored in memory are deleted; 
425   for(Int_t i=0;i<NumTables;i++)delete table[i];
426   delete [] table;
427   return 0;
428 }
429 /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
430 Int_t AliTPCCompression::RetrieveTables(AliTPCHTable* table[],Int_t NumTable){
431   //This method retrieve the Huffman tables from a sequence of binary files
432   if (fVerbose)
433     cout<<"Retrieving tables from files \n";
434   //  ULong_t code;
435   Double_t code;
436   UChar_t codeLen;
437   ifstream fTable;  
438   char filename[15];
439   //The following for loop is used to generate the Huffman trees acording to the tables
440   for(Int_t k=0;k<NumTable;k++){
441     Int_t dim;//this variable contains the table dimension
442     sprintf(filename,"Table%d.dat",k); 
443 #ifndef __DECCXX 
444     fTable.open(filename,ios::binary);
445 #else
446     fTable.open(filename);
447 #endif
448     if(!fTable){cout<<"File doesn't exist:"<<filename<<endl; exit(1);}
449     fTable.read((char*)(&dim),sizeof(Int_t));
450     if (fVerbose)
451       cout<<"Table dimension: "<<dim<<endl;
452     table[k]=new AliTPCHTable(dim);
453     for(Int_t i=0;i<dim;i++){
454       fTable.read((char*)(&codeLen),sizeof(UChar_t));
455       table[k]->SetCodeLen(codeLen,i);
456       //      fTable.read((char*)(&code),sizeof(ULong_t));
457       fTable.read((char*)(&code),sizeof(Double_t));
458       table[k]->SetCode(Mirror((ULong_t)code,codeLen),i);
459     }//end for 
460     fTable.close();
461   }//end for 
462   if (fVerbose)
463     cout<<"Trees generated \n";
464   //At this point the trees are been built
465   return 0;
466 }
467
468 Int_t AliTPCCompression::CreateTablesFromTxtFiles(Int_t NumTable){
469   //This method creates a set of binary tables, needed by the Huffman
470   //algorith, starting from a set of frequencies tables stored in form of
471   //txt files
472   if (fVerbose)
473     cout<<"Retrieving frequencies from txt files \n";
474   ifstream fTable;  
475   char filename[15];
476   //Tables are read from the files (Each codeword has been "Mirrored")
477   AliTPCHTable **table = new AliTPCHTable*[NumTable];
478   for(Int_t k=0;k<NumTable;k++){
479     sprintf(filename,"Table%d.txt",k); 
480     cout<<filename<<endl;
481     fTable.open(filename);
482     if(!fTable){cout<<"File doesn't exist: "<<filename<<endl; exit(1);}
483     Int_t symbol=0;
484     Double_t freq=0;
485     table[k]=new AliTPCHTable(1024);
486     while(!fTable.eof()){
487       fTable>>freq;
488       if (fTable.good()){
489         if (freq<0){
490           cout<<"Frequency cannot be negative !!!\n";
491           exit(1);
492         }
493         table[k]->SetValFrequency(symbol,freq);
494       }
495       symbol++;
496     }//end while
497     fTable.clear();
498     fTable.close();
499   }//end for
500   fStat.open("Statistics",ios::app);
501   fStat<<endl;
502   fStat<<"----------------- ENTROPY for external txt tables --------------------------"<<endl;
503   fStat<<"Entropy of Bunch length table......."<<table[0]->GetEntropy()<<endl;
504   fStat<<"Entropy of Time bin table..........."<<table[1]->GetEntropy()<<endl;
505   fStat<<"Entropy of one Sample bunch table..."<<table[2]->GetEntropy()<<endl;
506   fStat<<"Entropy of Central Sample table....."<<table[3]->GetEntropy()<<endl;
507   fStat<<"Entropy Border Samples table........"<<table[4]->GetEntropy()<<endl;
508   fStat.close();
509   for(Int_t k=0;k<NumTable;k++){
510     table[k]->BuildHTable();
511   }//end for
512   //The tables are saved ad binary files
513   StoreTables(table,NumTable);  
514   //The tables stored in memory are deleted; 
515   for(Int_t i=0;i<NumTable;i++)delete table[i];
516   delete [] table;
517   return 0;
518 }
519
520 ////////////////////////////////////////////////////////////////////////////////////////
521 /*                               COMPRESSION                                          */
522 ////////////////////////////////////////////////////////////////////////////////////////
523
524 void AliTPCCompression::StoreValue(ULong_t val,UChar_t len){
525   //This method stores the value "val" of "len" bits into the internal buffer "fBuffer"
526   if (len<=fFreeBitsBuffer){           // val is not splitted in two buffer
527     fFreeBitsBuffer-=len;
528     fBuffer=fBuffer<<len;
529     fBuffer=fBuffer|val;    
530     if(!fFreeBitsBuffer){              // if the buffer is full it is written into a file 
531       f.write((char*)(&fBuffer),sizeof(ULong_t));       
532       fFreeBitsBuffer=fDimBuffer;
533       fBuffer=0;
534     }
535   }//end if
536   else{                               //val has to be splitted in two buffers
537     fBuffer=fBuffer<<fFreeBitsBuffer;
538     ULong_t temp;
539     temp=val;
540     temp=temp>>(len-fFreeBitsBuffer);
541     fBuffer=fBuffer|temp;
542     f.write((char*)(&fBuffer),sizeof(ULong_t));
543     fFreeBitsBuffer=fDimBuffer-(len-fFreeBitsBuffer);
544     val=val<<fFreeBitsBuffer;
545     val=val>>fFreeBitsBuffer;
546     fBuffer=val;
547   }//end else
548   return;
549 }
550 //////////////////////////////////////////////////////////////////////////////////////////////////
551 void AliTPCCompression::Flush(){
552   //The last buffer cannot be completely full so to save it 
553   //into the output file it is first necessary to fill it with an hexadecimal pattern
554   if(fFreeBitsBuffer<fDimBuffer){
555     fBuffer=fBuffer<<fFreeBitsBuffer;
556     f.write((char*)(&fBuffer),sizeof(ULong_t));  
557   }//end if
558   return;
559 }
560 //////////////////////////////////////////////////////////////////////////////////////////////////
561 ULong_t AliTPCCompression::Mirror(ULong_t val,UChar_t len)const{
562   //This method inverts the digits of the number "val" and length "len"
563   //indicates the number of digits of the number considered in binary notation
564   ULong_t specular=0;
565   ULong_t mask=0x1;
566   ULong_t bit;
567   for(Int_t i=0;i<len;i++){
568     bit=val&mask;
569     bit=bit>>i;
570     specular=specular<<1;
571     specular=specular|bit;
572     mask=mask<<1;
573   }
574   return specular;
575 }
576
577 /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
578 Int_t AliTPCCompression::CompressDataOptTables(Int_t NumTable,const char* fSource,const char* fDest){
579   //This method compress an Altro format file using a general set of tables stored as binary files to be provided
580   if (fVerbose){
581     cout<<" BackWord COMPRESSION "<<endl;
582     cout<<"compression of the file "<<fSource<<" Output File: "<<fDest<<endl;
583   }
584   //Tables are read from the files (Each codeword has been "Mirrored")
585   AliTPCHTable **table = new AliTPCHTable*[NumTable];
586   RetrieveTables(table,NumTable);
587   //the output file is open
588 #ifndef __DECCXX 
589   f.open(fDest,ios::binary|ios::out);
590 #else
591   f.open(fDest,ios::out);
592 #endif
593   // Source file is open
594   AliTPCBuffer160 buff(fSource,0);
595   //coded words are written into a file
596   Int_t numWords,padNum,rowNum,secNum=0;
597   ULong_t  storedWords=0;
598   Int_t    value=0;
599   ULong_t  numPackets=0;
600   Double_t stat[5]={0.,0.,0.,0.,0.};
601   ULong_t  trailerNumbers=0;
602   Double_t numElem[5]={0,0,0,0,0};
603   Double_t fillWords=0.;
604   fStat.open("Statistics",ios::app);
605   fStat<<endl;
606   fStat<<"-------------------COMPRESSION STATISTICS----------"<<endl;
607   Int_t end=1;
608   while(buff.ReadTrailerBackward(numWords,padNum,rowNum,secNum) !=-1 ){
609     if(end){
610       fillWords=buff.GetFillWordsNum();
611       end=0;
612     }//endif
613
614     numPackets++;
615     if (numWords%4){
616       fillWords+=4-numWords%4;
617       for(Int_t j=0;j<(4-numWords%4);j++){
618         value=buff.GetNextBackWord();
619       }//end for
620     }//end if
621
622     Int_t packet[1024];
623     Int_t timePos[345];
624     Int_t tp=0;
625     for(Int_t i=0;i<345;i++)timePos[i]=0;
626     for(Int_t i=0;i<1024;i++)packet[i]=0;
627
628     Int_t nextTableType=0;
629     Int_t bunchLen=0;
630     Int_t count=0;
631     for(Int_t i=0;i<numWords;i++){
632       value=buff.GetNextBackWord();
633       packet[i]=value;
634       if(nextTableType==1){
635         timePos[tp]=i;
636         tp++;
637       }
638       NextTable(value,nextTableType,bunchLen,count);
639     }//end for
640     //computing the Time gap between two bunches
641     Int_t temp=0;
642     tp--;
643     Int_t previousTime=packet[timePos[tp]];
644     for(Int_t i=tp-1;i>=0;i--){
645       Int_t timPos=timePos[i];
646       Int_t bunchLen=packet[timPos-1]-2;
647       temp=packet[timPos];
648       packet[timPos]=packet[timPos]-previousTime-bunchLen;
649       previousTime=temp;
650     }//end for
651
652     nextTableType=0;
653     count=0;
654     bunchLen=0;
655     Int_t timeBin=0;
656     for(Int_t i=0;i<numWords;i++){
657       value=packet[i];
658       if(nextTableType==1)timeBin=value;
659       if(nextTableType>1){
660         //ULong_t val=(ULong_t)table[nextTableType]->Code()[value];     // val is the code
661         Double_t val=table[nextTableType]->Code()[value];     // val is the code
662         UChar_t len=table[nextTableType]->CodeLen()[value];  // len is the length (number of bits)of val
663         stat[nextTableType]+=len;
664         numElem[nextTableType]++;
665         StoreValue((ULong_t)val,len);
666         storedWords++;
667       }//end if
668       NextTable(value,nextTableType,bunchLen,count);
669       if(nextTableType==0){
670         //      ULong_t val=(ULong_t)table[1]->Code()[timeBin];     // val is the code
671         Double_t val=table[1]->Code()[timeBin];     // val is the code
672         UChar_t len=table[1]->CodeLen()[timeBin];  // len is the length (number of bits)of val
673         stat[1]+=len;
674         numElem[1]++;
675         StoreValue((ULong_t)val,len);
676         //      val=(ULong_t)table[nextTableType]->Code()[(bunchLen+2)];     // val is the code
677         val=table[nextTableType]->Code()[(bunchLen+2)];     // val is the code
678         len=table[nextTableType]->CodeLen()[(bunchLen+2)];  // len is the length (number of bits)of val
679         StoreValue((ULong_t)val,len);
680         stat[nextTableType]+=len;
681         numElem[nextTableType]++;
682         storedWords+=2;
683       }
684     }//end for
685     //Trailer
686     StoreValue(numWords,10);
687     StoreValue(padNum,10);
688     StoreValue(rowNum,10);
689     StoreValue(secNum,9);
690     StoreValue(1,1);
691     storedWords+=4;
692     trailerNumbers++;
693   }//end  while
694   StoreValue(numPackets,32);
695   if(fVerbose)
696     cout<<"Number of strored packets: "<<numPackets<<endl;
697   StoreValue(1,1);
698   //The last buffen cannot be completely full
699   Flush();
700   if(fVerbose)
701     cout<<"Number of stored words: "<<storedWords<<endl;
702   f.close();
703   //Tables are deleted
704   for(Int_t i=0;i<NumTable;i++){
705     delete table[i];
706   }//end for
707   delete [] table;
708   Double_t dimension=(ULong_t)TMath::Ceil((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])/8)+trailerNumbers*5;
709   fStat<<"Trailer Dimension in bytes......"<<trailerNumbers*5<<endl;
710   fStat<<"Data Dimension in bytes........."<<(ULong_t)TMath::Ceil((stat[0]+stat[1]+stat[2]+stat[3]+stat[4])/8)<<endl;
711   fStat<<"Compressed file dimension......."<<(ULong_t)dimension<<endl;
712   /*
713   fStat<<(ULong_t)trailerNumbers<<endl;
714   fStat<<(ULong_t)fillWords<<endl;
715   fStat<<(ULong_t)numElem[0]<<endl;
716   fStat<<(ULong_t)numElem[1]<<endl;
717   fStat<<(ULong_t)numElem[2]<<endl;
718   fStat<<(ULong_t)numElem[3]<<endl;
719   fStat<<(ULong_t)numElem[4]<<endl;
720   */
721   fillWords=(fillWords+numElem[0]+numElem[1]+numElem[2]+numElem[3]+numElem[4]+trailerNumbers*4)*10/8;
722   fStat<<"Original file dimension........."<<(ULong_t)fillWords<<endl;
723
724   Double_t ratio=(dimension/fillWords)*100;
725   fStat<<"Compression ratio (Compressed/Uncompressed)..."<<ratio<<"%"<<endl;
726   fStat<<endl;
727   if (numElem[0])
728     fStat<<"Bunch length size in bytes......"<<(ULong_t)TMath::Ceil(stat[0]/8)<<" Comppression.."<<(stat[0]/numElem[0])*10<<"%"<<endl;
729   if (numElem[1])  
730     fStat<<"Time gap size in bytes.........."<<(ULong_t)TMath::Ceil(stat[1]/8)<<" Comppression.."<<(stat[1]/numElem[1])*10<<"%"<<endl;
731   if (numElem[2]+numElem[3]+numElem[4])  
732     fStat<<"Amplitude values in bytes......."<<(ULong_t)TMath::Ceil((stat[2]+stat[3]+stat[4])/8)<<" Comppression.."<<
733       ((stat[2]+stat[3]+stat[4])/(numElem[2]+numElem[3]+numElem[4]))*10<<"%"<<endl;
734   if (numElem[2])
735   fStat<<"     One Samples in bytes............"<<(ULong_t)TMath::Ceil(stat[2]/8)<<" Comppression.."<<(stat[2]/numElem[2])*10<<"%"<<endl;
736   if (numElem[3])
737   fStat<<"     Central Samples size in bytes..."<<(ULong_t)TMath::Ceil(stat[3]/8)<<" Comppression.."<<(stat[3]/numElem[3])*10<<"%"<<endl;
738   if (numElem[4])
739   fStat<<"     Border Samples size in bytes...."<<(ULong_t)TMath::Ceil(stat[4]/8)<<" Comppression.."<<(stat[4]/numElem[4])*10<<"%"<<endl;
740   fStat<<endl;
741   fStat<<"Average number of bits per word"<<endl;
742   if (numElem[0])
743     fStat<<"Bunch length ......"<<stat[0]/numElem[0]<<endl;
744   if (numElem[1])
745     fStat<<"Time gap .........."<<stat[1]/numElem[1]<<endl;
746   if (numElem[2])
747     fStat<<"One Samples........"<<stat[2]/numElem[2]<<endl;
748   if (numElem[3])
749     fStat<<"Central Samples ..."<<stat[3]/numElem[3]<<endl;
750   if (numElem[4])
751     fStat<<"Border Samples....."<<stat[4]/numElem[4]<<endl;
752   fStat.close();
753   return 0;
754 }
755
756 ////////////////////////////////////////////////////////////////////////////////////////
757
758 ////////////////////////////////////////////////////////////////////////////////////////
759 /*                               DECOMPRESSION                                        */
760 ////////////////////////////////////////////////////////////////////////////////////////
761
762 void AliTPCCompression::CreateTreesFromFile(AliTPCHNode *RootNode[],const Int_t NumTables){
763   //For each table this method builds the associate Huffman tree starting from the codeword and 
764   //the codelength of each symbol 
765   if(fVerbose)
766     cout<<"Creating the Huffman trees \n";
767   AliTPCHNode *node=0;
768   // ULong_t code;
769   Double_t code;
770   UChar_t codeLen;
771   ifstream fTable;  
772   char filename[15];
773   //The following for loop is used to generate the Huffman trees acording to the tables
774   //loop over the tables
775   for(Int_t k=0;k<NumTables;k++){
776     RootNode[k]=new AliTPCHNode(); //RootNode is the root of the tree
777     Int_t dim=0;//this variable contains the table dimension
778     sprintf(filename,"Table%d.dat",k); 
779 #ifndef __DECCXX 
780     fTable.open(filename,ios::binary);
781 #else
782     fTable.open(filename);
783 #endif
784     if(!fTable){cout<<"Tables don't exist !!!"<<endl;exit(1);}
785     fTable.read((char*)(&dim),sizeof(Int_t));
786     if (fVerbose)
787       cout<<"Table dimension: "<<dim<<endl;
788     //loop over the words of one table
789     for(Int_t i=0;i<dim;i++){
790       fTable.read((char*)(&codeLen),sizeof(UChar_t));
791       //fTable.read((char*)(&code),sizeof(ULong_t));
792       fTable.read((char*)(&code),sizeof(Double_t));
793       node=RootNode[k];
794       for(Int_t j=1;j<=codeLen;j++){
795         ULong_t bit,val=0;
796         val=(ULong_t)TMath::Power(2,codeLen-j);
797         bit=(ULong_t)code&val; 
798         AliTPCHNode *temp=node;
799         if(bit){
800           node=node->GetRight();
801           if(!node){
802             node=new AliTPCHNode();
803             temp->SetRight(node);
804           }//end if 
805         }//end if
806         else{
807           node=node->GetLeft();
808           if(!node){
809             node=new AliTPCHNode();
810             temp->SetLeft(node);
811           }//end if
812         }//end else
813       }//end for
814       if(codeLen){
815         node->SetSymbol(i);
816         node->SetFrequency(codeLen);
817       }//end if
818     }//end for 
819     fTable.close();
820   }//end for 
821   if (fVerbose)
822     cout<<"Trees generated \n";
823   //At this point the trees are been built
824 }
825 //////////////////////////////////////////////////////////////////////////////////////////////////
826 void AliTPCCompression::DeleteHuffmanTree(AliTPCHNode* node){
827   //This function deletes all the nodes of an Huffman tree
828   //In an Huffman tree any internal node has always two children 
829   if (node){
830     DeleteHuffmanTree(node->GetLeft());
831     DeleteHuffmanTree(node->GetRight());
832     //    cout<<node->GetSymbol()<<"  "<<(Int_t)node->GetFrequency()<<endl;
833     delete node;
834   }
835 }
836 //////////////////////////////////////////////////////////////////////////////////////////////////
837 void AliTPCCompression::VisitHuffmanTree(AliTPCHNode* node){
838   //This function realizes an in order visit of a binary tree 
839   if (node){
840     cout<<node->GetSymbol()<<" "<<node->GetFrequency()<<endl;
841     VisitHuffmanTree(node->GetLeft());
842     VisitHuffmanTree(node->GetRight());
843   }
844 }
845 //////////////////////////////////////////////////////////////////////////////////////////////////
846 ULong_t AliTPCCompression::ReadWord(Int_t NumberOfBit){
847   //This method retrieves a word of a specific number of bits from the file through the internal buffer 
848   ULong_t result=0;
849   ULong_t bit=0;
850   for (Int_t i=0;i<NumberOfBit;i++){
851     if (fReadBits==32){
852       fPos-=sizeof(ULong_t);
853       f.seekg(fPos);
854       f.read((char*)(&fBuffer),sizeof(ULong_t));
855       fReadBits=0;
856     }//end if
857     ULong_t mask=0;
858     mask=(ULong_t)TMath::Power(2,fReadBits);
859     bit=fBuffer&mask;
860     bit=bit>>fReadBits;
861     fReadBits++;
862     bit=bit<<i;
863     result=result|bit;
864   }//end for
865   return result;
866 }
867 //////////////////////////////////////////////////////////////////////////////////////////////////
868 ULong_t AliTPCCompression::ReadWordBuffer(Int_t NumberOfBit){
869   //This method retrieves a word of a specific number of bits from the file through the buffer 
870   ULong_t result=0;
871   ULong_t bit=0;
872   for (Int_t i=0;i<NumberOfBit;i++){
873     if (fReadBits==32){
874       fPointBuffer-=8;
875       fBuffer=0;
876       for(Int_t i=0;i<4;i++){
877         ULong_t val=0;
878         val=*fPointBuffer;
879         val&=0xFF;
880         fPointBuffer++;
881         val<<=8*i;
882         fBuffer=fBuffer|val;
883       }//end for
884       fReadBits=0;
885     }//end if
886     ULong_t mask=0;
887     mask=(ULong_t)TMath::Power(2,fReadBits);
888     bit=fBuffer&mask;
889     bit=bit>>fReadBits;
890     fReadBits++;
891     bit=bit<<i;
892     result=result|bit;
893   }//end for
894   return result;
895 }
896
897 //////////////////////////////////////////////////////////////////////////////////////////////////
898 void AliTPCCompression::ReadTrailer(Int_t &WordsNumber,Int_t &PadNumber,Int_t &RowNumber,Int_t &SecNumber,Bool_t Memory){
899   //It retrieves a trailer 
900   if(Memory){
901     ReadWordBuffer(1);
902     SecNumber=ReadWordBuffer(9);
903     RowNumber=ReadWordBuffer(10);
904     PadNumber=ReadWordBuffer(10);
905     WordsNumber=ReadWordBuffer(10);
906   }
907   else{
908     ReadWord(1);
909     SecNumber=ReadWord(9);
910     RowNumber=ReadWord(10);
911     PadNumber=ReadWord(10);
912     WordsNumber=ReadWord(10);
913   }
914   return;
915 }
916 //////////////////////////////////////////////////////////////////////////////////////////////////
917 ULong_t AliTPCCompression::GetDecodedWord(AliTPCHNode* root,Bool_t Memory){
918   //This method retrieves a decoded word.
919   AliTPCHNode *node=root;
920   ULong_t symbol=0;
921   Bool_t decoded=0;
922   while(!decoded){
923     ULong_t bit=0;
924     if(Memory)
925       bit=ReadWordBuffer(1);
926     else
927       bit=ReadWord(1);
928     if(bit)
929       node=node->GetRight();
930     else
931       node=node->GetLeft();
932     if (!(node->GetLeft())){
933       symbol=node->GetSymbol();
934       decoded=1;
935     }
936   }//end while
937   return symbol;
938 }
939 //////////////////////////////////////////////////////////////////////////////////////////////////
940
941 Int_t AliTPCCompression::DecompressDataOptTables(Int_t NumTables,const char* fname,char* fDest){
942   //This method decompress a file using separate Huffman tables
943   if(fVerbose){
944     cout<<"   DECOMPRESSION:"<<endl;
945     cout<<"Source File "<<fname<<" Destination File "<<fDest<<endl; 
946   }
947   AliTPCHNode ** rootNode = new AliTPCHNode*[NumTables];
948   //Creation of the Huffman trees
949   CreateTreesFromFile(rootNode,NumTables);
950 #ifndef __DECCXX
951   f.open(fname,ios::binary|ios::in);
952 #else
953   f.open(fname,ios::in);
954 #endif
955   if(!f){cout<<"File doesn't exist:"<<fname<<endl;;return -1;}
956   //to go to the end of the file
957   f.seekg(0,ios::end);
958   //to get the file dimension in byte
959   fPos=f.tellg();
960   fPos-=sizeof(ULong_t);
961   f.seekg(fPos);
962   fReadBits=0;
963   fBuffer=0;
964   f.read((char*)(&fBuffer),sizeof(ULong_t));
965   Int_t bit=0;
966   ULong_t mask=0x1;
967   while(!bit){
968     bit=fBuffer&mask;
969     mask=mask<<1;
970     fReadBits++;
971   }
972   ULong_t packetNumber=ReadWord(sizeof(ULong_t)*8);
973   if(fVerbose){
974     cout<<"Number of Packect: "<<packetNumber<<endl;
975   }
976   AliTPCBuffer160 bufferFile(fDest,1);
977   ULong_t k=0;
978   ULong_t wordsRead=0; //number of read coded words 
979   while(k<packetNumber){
980     Int_t numWords,padNumber,rowNumber,secNumber=0;
981     ReadTrailer(numWords,padNumber,rowNumber,secNumber,kFALSE);
982     k++;
983     wordsRead+=4;
984     Int_t previousTime=-1;
985     Int_t time=0;
986     Int_t nextTableType=0;
987     Int_t bunchLen=0;
988     Int_t count=0;
989     for(Int_t i=0;i<numWords;i++){
990       ULong_t symbol=GetDecodedWord(rootNode[nextTableType],kFALSE);
991       wordsRead++;
992       //Time reconstruction
993       if (nextTableType==1){
994         if (previousTime!=-1){
995           previousTime=symbol+previousTime+bunchLen;
996         }
997         else previousTime=symbol;
998         time=previousTime;
999       }
1000       if(nextTableType>1)
1001         bufferFile.FillBuffer(symbol);
1002       NextTable(symbol,nextTableType,bunchLen,count); 
1003       if(nextTableType==0){
1004         bufferFile.FillBuffer(time);
1005         bufferFile.FillBuffer(bunchLen+2);
1006         bunchLen=0;
1007       }
1008     }//end for
1009     bufferFile.WriteTrailer(numWords,padNumber,rowNumber,secNumber);
1010   }//end while
1011   if(fVerbose){
1012     cout<<"Number of decoded words:"<<wordsRead<<endl;
1013   }
1014   f.close();
1015   //The trees are deleted 
1016   for(Int_t j=0;j<NumTables;j++){
1017       DeleteHuffmanTree(rootNode[j]);
1018   }//end for
1019   delete [] rootNode;
1020   return 0; 
1021 }
1022
1023 //////////////////////////////////////////////////////////////////////////////////////////////////
1024 Int_t AliTPCCompression::Decompress(AliTPCHNode *RootNode[],const Int_t NumTables,char* PointBuffer,ULong_t BufferSize,UShort_t out[],ULong_t &dim){
1025   //This method decompress a file using separate Huffman tables
1026
1027   fPointBuffer=PointBuffer+BufferSize-4;
1028   fReadBits=0;
1029   fBuffer=0;
1030   
1031   for(Int_t i=0;i<4;i++){
1032     ULong_t val=0;
1033     val=*fPointBuffer;
1034     val&=0xFF;
1035     fPointBuffer++;
1036     val<<=8*i;
1037     fBuffer=fBuffer|val;
1038   }//end for
1039   Int_t bit=0;
1040   ULong_t mask=0x1;
1041   while(!bit){
1042     bit=fBuffer&mask;
1043     mask=mask<<1;
1044     fReadBits++;
1045   }//end while
1046   ULong_t packetNumber=ReadWordBuffer(sizeof(ULong_t)*8); //32 bits
1047   if (fVerbose){
1048     cout<<"First one has been found "<<endl;
1049     cout<<"Number of packets:"<<packetNumber<<endl;
1050   }//end if
1051   ULong_t k=0;
1052   ULong_t wordsRead=0; //number of read coded words
1053   while(k<packetNumber){
1054     Int_t numWords,padNumber,rowNumber,secNumber=0;
1055     ReadTrailer(numWords,padNumber,rowNumber,secNumber,kTRUE);
1056     out[dim]=numWords;
1057     dim++;
1058     out[dim]=padNumber;
1059     dim++;
1060     out[dim]=rowNumber;
1061     dim++;
1062     out[dim]=secNumber;
1063     dim++;
1064     //ftxt<<"S:"<<secNumber<<" R:"<<rowNumber<<" P:"<<padNumber<<" W:"<<numWords<<endl;
1065     //    padDigits->SetPadID(padNumber,rowNumber,secNumber,DDL);
1066     k++;
1067     wordsRead+=4;
1068     Int_t previousTime=-1;
1069     Int_t time=0;
1070     Int_t nextTableType=0;
1071     Int_t bunchLen=0;
1072     Int_t count=0;
1073     Int_t timeDigit=0;
1074     for(Int_t i=0;i<numWords;i++){
1075       ULong_t symbol=GetDecodedWord(RootNode[nextTableType],kTRUE);
1076       wordsRead++;
1077       //Time reconstruction
1078       if (nextTableType==1){
1079         if (previousTime!=-1){
1080           previousTime=symbol+previousTime+bunchLen;
1081         }
1082         else previousTime=symbol;
1083         time=previousTime;
1084         out[dim]=bunchLen+2;
1085         dim++;
1086         out[dim]=time;
1087         dim++;
1088         timeDigit=time-bunchLen;
1089       }
1090       if(nextTableType>1){
1091         //
1092         //ftxt<<symbol<<endl;
1093         out[dim]=symbol;
1094         dim++;
1095         timeDigit++;
1096         //padDigits->SetDigits(symbol,timeDigit);
1097       }
1098       NextTable(symbol,nextTableType,bunchLen,count); 
1099       if(nextTableType==0){
1100         //
1101         //ftxt<<time<<endl;
1102         //  ftxt<<(bunchLen+2)<<endl;
1103         bunchLen=0;
1104       }
1105     }//end for
1106   }//end while
1107   return 0; 
1108 }
1109
1110 //////////////////////////////////////////////////////////////////////////////////////////////////
1111 Int_t AliTPCCompression::DestroyTables(AliTPCHNode *RootNode[],const Int_t NumTables){
1112   //The trees are deleted
1113   for(Int_t j=0;j<NumTables;j++){
1114     DeleteHuffmanTree(RootNode[j]);
1115   }//end for
1116   if(fVerbose)
1117     cout<<"Huffman trees destroyed"<<endl;
1118   return 0;
1119 }
1120 //////////////////////////////////////////////////////////////////////////////////////////////////
1121
1122 void AliTPCCompression::ReadAltroFormat(char* fileOut,char* fileIn)const{
1123   //This method creates a text file containing the same information stored in 
1124   //an Altro file. The information in the text file is organized pad by pad and 
1125   //and for each pad it consists in a sequence of bunches (Bunch length +2,
1126   //Time bin of the last amplitude sample in the bunch, amplitude values)
1127   //It is used mainly for debugging
1128   ofstream ftxt(fileOut);
1129   AliTPCBuffer160 buff(fileIn,0);
1130   Int_t numWords,padNum,rowNum,secNum=0;
1131   Int_t value=0;
1132   if (fVerbose) cout<<"Creating a txt file from an Altro Format file"<<endl;
1133   while(buff.ReadTrailerBackward(numWords,padNum,rowNum,secNum) !=-1 ){
1134     ftxt<<"S:"<<secNum<<" R:"<<rowNum<<" P:"<<padNum<<" W:"<<numWords<<endl;
1135     if (numWords%4){
1136       for(Int_t j=0;j<(4-numWords%4);j++){
1137         value=buff.GetNextBackWord();
1138       }//end for
1139     }//end if
1140     for(Int_t i=0;i<numWords;i++){
1141       value=buff.GetNextBackWord();
1142       ftxt<<value<<endl;
1143     }//end for
1144   }//end while
1145   ftxt.close();
1146   return;
1147 }
1148
1149 //////////////////////////////////////////////////////////////////////////////////////////