]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/comp/AliL3Compress.cxx
Moving to static AliL3Transform.
[u/mrichter/AliRoot.git] / HLT / comp / AliL3Compress.cxx
1 //$Id$
2
3 // Author: Anders Vestbo <mailto:vestbo$fi.uib.no>
4 //*-- Copyright &copy ASV
5
6 #include <stdio.h>
7 #include <stream.h>
8 #include <stdlib.h>
9 #include <TH1.h>
10 #include <TH2.h>
11 #include <TRandom.h>
12
13 #include "AliL3Compress.h"
14 #include "AliL3TrackArray.h"
15 #include "AliL3ModelTrack.h"
16 #include "AliL3Transform.h"
17 #include "AliL3MemHandler.h"
18 #ifdef use_aliroot
19 #include "AliL3FileHandler.h"
20 #endif
21 #include "bitio.h"
22
23 //_____________________________________________________________
24 //
25 //  AliL3Compress
26 //
27 // Class for compressing and uncompressing data.
28
29 ClassImp(AliL3Compress)
30
31 AliL3Compress::AliL3Compress()
32 {
33   fTracks=0;
34   SetBitNumbers(0,0,0,0);
35   fSlice =0;
36   fPatch=0;
37   fDigits=0;
38   fDPt=0;
39 }
40
41 AliL3Compress::AliL3Compress(Int_t slice,Int_t patch,Char_t *path)
42 {
43   fSlice=slice;
44   fPatch=patch;
45   SetBitNumbers(0,0,0,0);
46   fTracks=0;
47   fDigits=0;
48   fDPt=0;
49   sprintf(fPath,"%s",path);
50 }
51
52 AliL3Compress::~AliL3Compress()
53 {
54   if(fTracks)
55     delete fTracks;
56   if(fDigits)
57     delete [] fDigits;
58   if(fDPt)
59     delete [] fDPt;
60 }
61
62 void AliL3Compress::SetBitNumbers(Int_t pad,Int_t time,Int_t charge,Int_t shape)
63 {
64   fNumPadBits=pad;
65   fNumTimeBits=time;
66   fNumChargeBits=charge;
67   fNumShapeBits=shape;
68 }
69
70 void AliL3Compress::WriteFile(AliL3TrackArray *tracks)
71 {
72   Char_t fname[100];
73   sprintf(fname,"%s/tracks_m_%d_%d.raw",fPath,fSlice,fPatch);
74   FILE *file = fopen(fname,"w");
75   Short_t ntracks = tracks->GetNTracks();
76   //cout<<"Writing "<<ntracks<<" tracks to file"<<endl;
77     
78   Int_t count=0;
79   AliL3ClusterModel *clusters=0;
80   AliL3TrackModel *model=0;
81   for(Int_t i=0; i<ntracks; i++)
82     {
83       AliL3ModelTrack *track = (AliL3ModelTrack*)tracks->GetCheckedTrack(i);
84       if(!track) continue;
85       
86       //Do not save useless tracks or clusters:
87       //if(track->CheckClustersQuality(2) == 0) 
88       //continue;
89       
90       track->FillModel();
91       model = track->GetModel();
92       if(model->fNClusters==0) continue;
93       clusters = track->GetClusters();
94       //cout<<"Writing "<<(int)model->fNClusters<<" clusters"<<endl;
95       if(fwrite(model,sizeof(AliL3TrackModel),1,file)!=1) break;
96       //cout<<"Writing "<<(int)model->fNClusters<<" clusters to file"<<endl;
97       if(fwrite(clusters,model->fNClusters*sizeof(AliL3ClusterModel),1,file)!=1) break;
98       //track->Print();
99       count++;
100       
101     }
102   cout<<"Wrote "<<count<<" tracks "<<endl;
103   fclose(file);
104 }
105
106 void AliL3Compress::ReadFile(Char_t which)
107 {
108   //Read the trackfile.
109
110   Char_t fname[100];
111   if(which == 'm')
112     sprintf(fname,"%s/tracks_m_%d_%d.raw",fPath,fSlice,fPatch);
113   else if(which == 'u')
114     sprintf(fname,"%s/tracks_u_%d_%d.raw",fPath,fSlice,fPatch);
115   else
116     {
117       cerr<<"AliL3Compress::ReadFile() : Wrong option"<<endl;
118       return;
119     }
120
121   FILE *file = fopen(fname,"r");
122   if(!file)
123     {
124       cerr<<"Cannot open file "<<fname<<endl;
125       return;
126     }
127
128   if(fTracks)
129     delete fTracks;
130   fTracks = new AliL3TrackArray("AliL3ModelTrack");
131   
132   cout<<"Reading file "<<fname<<endl;
133   while(!feof(file))
134     {
135       AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->NextTrack();
136       track->Init(fSlice,fPatch);
137       AliL3TrackModel *model = track->GetModel();
138       AliL3ClusterModel *clusters = track->GetClusters();
139       //cout<<"Reading model "<<(int)model<<endl;
140       if(fread(model,sizeof(AliL3TrackModel),1,file)!=1) break;
141       //cout<<"Reading clusters "<<(int)clusters<<endl;
142       if(fread(clusters,(model->fNClusters)*sizeof(AliL3ClusterModel),1,file)!=1) break;
143       //cout<<"Filling track"<<endl;
144       track->FillTrack();
145       //track->Print();
146     }
147
148   fTracks->RemoveLast();
149   cout<<"Read "<<fTracks->GetNTracks()<<" tracks from file"<<endl;
150   fclose(file);
151 }
152
153 void AliL3Compress::CompressFile()
154 {
155   if(fNumTimeBits==0)
156     {
157       cerr<<"AliL3Compress::CompressFile() : Bitnumbers not set"<<endl;
158       return;
159     }
160   
161   Char_t fname[100];
162   sprintf(fname,"%s/tracks_c_%d_%d.raw",fPath,fSlice,fPatch);
163   BIT_FILE *output = OpenOutputBitFile(fname);
164   
165   sprintf(fname,"%s/tracks_m_%d_%d.raw",fPath,fSlice,fPatch);
166   FILE *input = fopen(fname,"r");
167   if(!input)
168     {
169       cerr<<"AliL3Compress::CompressFile() : Error opening file: "<<fname<<endl;
170       return;
171     }
172
173   AliL3TrackModel track;
174   AliL3ClusterModel cluster;
175   Int_t temp;
176   Short_t power;
177   
178   Int_t timeo,pado,chargeo,shapeo;
179   timeo=pado=chargeo=shapeo=0;
180   while(!feof(input))
181     {
182       if(fread(&track,sizeof(AliL3TrackModel),1,input)!=1) break;
183       
184       if(output->mask != 0x80) //Write the current byte to file.
185         {
186           //cerr<<"\nAliL3Compress::CompressFile() : Writing overhead bits!!!"<<endl;
187           if(putc(output->rack,output->file )!=output->rack)
188             cerr<<"AliL3Compress::ComressFile : Error writing to bitfile"<<endl;
189           output->mask=0x80;
190           output->rack=0;
191         }
192       
193       //Write track parameters:
194       fwrite(&track,sizeof(AliL3TrackModel),1,output->file);
195       for(Int_t i=0; i<track.fNClusters; i++)
196         {
197           if(fread(&cluster,sizeof(AliL3ClusterModel),1,input)!=1) break;
198           
199           //Write empty flag:
200           temp = (Int_t)cluster.fPresent;
201           OutputBit(output,temp);
202           if(!temp) continue;
203           
204           //Write time information:
205           temp = (Int_t)cluster.fDTime;
206           if(temp<0)
207             OutputBit(output,0);
208           else
209             OutputBit(output,1);
210           power = 1<<(fNumTimeBits-1);
211           if(abs(temp)>=power)
212             {
213               timeo++;
214               temp=power - 1;
215             }
216           OutputBits(output,abs(temp),(fNumTimeBits-1));
217           
218           //Write pad information:
219           temp = (Int_t)cluster.fDPad;
220           if(temp<0)
221             OutputBit(output,0);
222           else
223             OutputBit(output,1);
224           power = 1<<(fNumPadBits-1);
225           if(abs(temp)>=power)
226             {
227               pado++;
228               temp=power - 1;
229             }
230           OutputBits(output,abs(temp),(fNumPadBits-1));
231           
232           //Write charge information:
233           temp = (Int_t)cluster.fDCharge;
234           if(temp<0)
235             OutputBit(output,0);
236           else
237             OutputBit(output,1);
238           power = 1<<(fNumChargeBits-1);
239           if(abs(temp)>=power)
240             {
241               chargeo++;
242               temp=power - 1;
243             }
244           OutputBits(output,abs(temp),(fNumChargeBits-1));
245                   
246           //Write shape information:
247           temp = (Int_t)cluster.fDSigmaY2;
248           if(temp<0)
249             OutputBit(output,0);
250           else
251             OutputBit(output,1);
252           power = 1<<(fNumShapeBits-1);
253           if(abs(temp) >= power)
254             {
255               shapeo++;
256               temp = power - 1;
257             }
258           OutputBits(output,abs(temp),(fNumShapeBits-1));
259           
260           temp = (Int_t)cluster.fDSigmaZ2;
261           if(temp<0)
262             OutputBit(output,0);
263           else
264             OutputBit(output,1);
265           power = 1<<(fNumShapeBits-1);
266           if(abs(temp) >= power)
267             {
268               shapeo++;
269               temp=power - 1;
270             }
271           OutputBits(output,abs(temp),(fNumShapeBits-1));
272         }
273     }
274   
275   fclose(input);
276   CloseOutputBitFile(output);
277
278   cout<<endl<<"Saturations: "<<endl
279       <<"Pad "<<pado<<endl
280       <<"Time "<<timeo<<endl
281       <<"Charge "<<chargeo<<endl
282       <<"Shape "<<shapeo<<endl;
283 }
284
285 void AliL3Compress::ExpandFile()
286 {
287   if(fNumTimeBits==0)
288     {
289       cerr<<"AliL3Compress::ExpandFile() : Bitnumbers not set"<<endl;
290       return;
291     }
292   
293   Char_t fname[100];
294   sprintf(fname,"%s/tracks_c_%d_%d.raw",fPath,fSlice,fPatch);
295   BIT_FILE *input = OpenInputBitFile(fname);
296   
297   sprintf(fname,"%s/tracks_u_%d_%d.raw",fPath,fSlice,fPatch);
298   FILE *output = fopen(fname,"w");
299   if(!output)
300     {
301       cerr<<"AliL3Compress::ExpandFile() : Error opening file: "<<fname<<endl;
302       return;
303     }
304
305   AliL3TrackModel trackmodel;
306   AliL3ClusterModel *clusters=0;
307   Int_t count=0;
308   
309   clusters = new AliL3ClusterModel[(AliL3Transform::GetNRows(fPatch))];
310   while(!feof(input->file))
311     {
312       input->mask=0x80;//make sure we read a new byte from file.
313       
314       //Read and write track:
315       if(fread(&trackmodel,sizeof(AliL3TrackModel),1,input->file)!=1) break;
316       fwrite(&trackmodel,sizeof(AliL3TrackModel),1,output);
317       
318       for(Int_t i=0; i<AliL3Transform::GetNRows(fPatch); i++)
319         {
320           Int_t temp,sign;
321           
322           //Read empty flag:
323           temp = InputBit(input);
324           if(!temp) 
325             {
326               clusters[i].fPresent=kFALSE;
327               continue;
328             }
329           clusters[i].fPresent=kTRUE;
330           
331           //Read time information:
332           sign=InputBit(input);
333           temp = InputBits(input,(fNumTimeBits-1));
334           if(!sign)
335             temp*=-1;
336           clusters[i].fDTime = temp;
337           
338           //Read pad information:
339           sign=InputBit(input);
340           temp = InputBits(input,(fNumPadBits-1));
341           if(!sign)
342             temp*=-1;
343           clusters[i].fDPad = temp;
344           
345           //Read charge information:
346           sign = InputBit(input);
347           temp=InputBits(input,(fNumChargeBits-1));
348           if(!sign)
349             temp*=-1;
350           clusters[i].fDCharge = temp;
351           
352           //Read shape information:
353           sign = InputBit(input);
354           temp = InputBits(input,(fNumShapeBits-1));
355           if(!sign)
356             temp*=-1;
357           clusters[i].fDSigmaY2 = temp;
358           
359           sign = InputBit(input);
360           temp = InputBits(input,(fNumShapeBits-1));
361           if(!sign)
362             temp*=-1;
363           clusters[i].fDSigmaZ2 = temp;
364         }
365       
366
367       count++;
368       fwrite(clusters,(trackmodel.fNClusters)*sizeof(AliL3ClusterModel),1,output);
369       
370     }
371   
372   delete [] clusters;
373   fclose(output);
374   CloseInputBitFile(input);
375 }
376
377 void AliL3Compress::CreateDigitArray(Int_t maxnumber)
378 {
379   fNUsed=0;
380   fNDigits = 0;
381   fMaxDigits=maxnumber;
382   if(fDigits) delete [] fDigits;
383   fDigits = new AliL3RandomDigitData[maxnumber];
384   if(fDPt) delete [] fDPt;
385   fDPt = new AliL3RandomDigitData*[maxnumber];
386 }
387
388 void AliL3Compress::RestoreData(Char_t which='u')
389 {
390   //Restore the data.
391   //which == u : restore compressed data
392   //which == m : restore uncompressed data
393   
394
395   ReadFile(which);
396   
397   CreateDigitArray(10000000);
398   
399   Float_t pad,time,sigmaY2,sigmaZ2;
400   Int_t charge;
401   for(Int_t j=AliL3Transform::GetFirstRow(fPatch); j<=AliL3Transform::GetLastRow(fPatch); j++)
402     {
403       cout<<"Building clusters on row "<<j<<endl;
404       for(Int_t i=0; i<fTracks->GetNTracks(); i++)
405         {
406           AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i);
407           if(!track) continue;
408           if(!track->GetPad(j,pad) || 
409              !track->GetTime(j,time) || 
410              !track->GetClusterCharge(j,charge) ||
411              !track->GetXYWidth(j,sigmaY2) || 
412              !track->GetZWidth(j,sigmaZ2))
413             continue;
414           
415           CreateDigits(j,pad,time,charge,sigmaY2,sigmaZ2);
416         }
417     }
418   
419   QSort(fDPt,0,fNDigits);
420 }
421
422 void AliL3Compress::PrintDigits(Int_t padrow=-1)
423 {
424   Int_t pad,time,charge,row;
425   for(Int_t i=0; i<fNDigits; i++)
426     {
427       row = fDPt[i]->fRow;
428       if(padrow > 0)
429         if(row != padrow) continue;
430       pad = fDPt[i]->fPad;
431       time = fDPt[i]->fTime;
432       charge = fDPt[i]->fCharge;
433       if(i>0 && row != fDPt[i-1]->fRow)
434         cout<<"---Padrow "<<row<<"---"<<endl;
435       cout<<"Padrow "<<row<<" Pad "<<pad<<" time "<<time<<" charge "<<charge<<endl;
436     }
437 }
438
439 void AliL3Compress::WriteRestoredData()
440 {
441   Char_t fname[100];
442   
443   //Get the remaining raw data array:
444   AliL3MemHandler *mem = new AliL3MemHandler();
445   sprintf(fname,"%s/remains_%d_%d.raw",fPath,fSlice,fPatch);
446   mem->SetBinaryInput(fname);
447   UInt_t numdigits;
448   AliL3DigitRowData *origRow = mem->CompBinary2Memory(numdigits);
449   mem->CloseBinaryInput();
450   
451   //Allocate memory for the merged data:
452   UInt_t size = mem->GetAllocatedSize() + fNDigits*sizeof(AliL3DigitData);
453   cout<<"Allocating "<<size<<" bytes for merged data array "<<endl;
454   Byte_t *data = new Byte_t[size];
455   memset(data,0,size);
456   AliL3DigitRowData *tempRow = (AliL3DigitRowData*)data;
457
458   Int_t ndigits,action,charge;
459   UShort_t pad,time;
460       
461   UInt_t digit_counter;
462   Int_t row_counter=0;
463   for(Int_t i=AliL3Transform::GetFirstRow(fPatch); i<=AliL3Transform::GetLastRow(fPatch); i++)
464     {
465       tempRow->fRow = i;
466       ndigits=0;
467       AliL3DigitData *origDig = origRow->fDigitData;
468       AliL3DigitData *tempDig = tempRow->fDigitData;
469       if((Int_t)origRow->fRow != i)
470         cerr<<"AliL3Compress::WriteRestoredData() : Mismatching row numbering "<<(Int_t)origRow->fRow<<" "<<i<<endl;
471
472       //cout<<"Writing row "<<i<<" with "<<(Int_t)origRow->fNDigit<<" old digits"<<endl;
473       digit_counter=0;
474       
475       while(1)
476         {
477           while(digit_counter < origRow->fNDigit)
478             {
479               pad = origDig[digit_counter].fPad;
480               time = origDig[digit_counter].fTime;
481               charge = origDig[digit_counter].fCharge;
482               digit_counter++;
483               while( (action=ComparePoints(i,pad,time)) == 1)
484                 {
485                   tempDig[ndigits].fPad = fDPt[fNUsed]->fPad;
486                   tempDig[ndigits].fTime = fDPt[fNUsed]->fTime;
487                   tempDig[ndigits].fCharge = fDPt[fNUsed]->fCharge;
488                   ndigits++;
489                   fNUsed++;
490
491                 }
492               if(action == 0)
493                 {
494                   tempDig[ndigits].fPad = pad;
495                   tempDig[ndigits].fTime = time;
496                   tempDig[ndigits].fCharge = charge;
497                   ndigits++;
498                 }
499             }
500           
501           if(fNUsed >= fNDigits) 
502             {
503               //cerr<<"AliL3Compress::WriteRestoredData() : Array out of range : "<<fNUsed<<" "<<fNDigits<<endl;
504               break;
505             }
506           if(fDPt[fNUsed]->fRow != i) //we are on a new row
507             break;
508           tempDig[ndigits].fPad = fDPt[fNUsed]->fPad;
509           tempDig[ndigits].fTime = fDPt[fNUsed]->fTime;
510           tempDig[ndigits].fCharge = fDPt[fNUsed]->fCharge;
511           ndigits++;
512           fNUsed++;
513         }
514       //cout<<"Writing "<<ndigits<<" digits on row "<<i<<endl;
515       if(ndigits > 4)
516         {
517           row_counter++;
518         }
519       tempRow->fNDigit = ndigits;
520       Int_t size = sizeof(AliL3DigitData)*tempRow->fNDigit + sizeof(AliL3DigitRowData);
521       Byte_t *byte_pt = (Byte_t*)tempRow;
522       byte_pt += size;
523       tempRow = (AliL3DigitRowData*)byte_pt;
524       mem->UpdateRowPointer(origRow);
525     }
526   
527   if(row_counter != AliL3Transform::GetNRows(fPatch))
528     cerr<<"AliL3Compress::WriteRestoredData() : Written rows: "<<row_counter<<" total rows "<<AliL3Transform::GetNRows(fPatch)<<endl;
529   
530   mem->Free();  
531   sprintf(fname,"%s/restored_%d_%d.raw",fPath,fSlice,fPatch);
532   mem->SetBinaryOutput(fname);
533   mem->Memory2CompBinary((UInt_t)AliL3Transform::GetNRows(fPatch),(AliL3DigitRowData*)data);
534   mem->CloseBinaryOutput();
535   
536   delete [] data;
537   delete mem;
538   
539 }
540
541 void AliL3Compress::CreateDigits(Int_t row,Float_t pad,Float_t time,Int_t charge,Float_t sigmaY2,Float_t sigmaZ2)
542 {
543   //Create raw data out of the cluster.
544   
545
546     
547   //Convert back the sigmas:
548   /*
549     Float_t padw,timew;
550     if(fPatch < 2)
551     padw = AliL3Transform::GetPadPitchWidthLow();
552     else
553     padw = AliL3Transform::GetPadPitchWidthUp();
554     timew = AliL3Transform::GetZWidth();
555   */
556   
557   /*
558     if(fPatch < 3)
559     sigmaY2 = sigmaY2/2.07;
560     sigmaY2 = sigmaY2/0.108;
561   */
562   
563   //sigmaY2 = sigmaY2/(padw*padw);
564   
565   //sigmaY2 = sigmaY2;// - 1./12;
566   
567   /*
568     if(fPatch < 3)
569     sigmaZ2 = sigmaZ2/1.77;
570     sigmaZ2 = sigmaZ2/0.169;
571   */
572   
573   //sigmaZ2 = sigmaZ2/(timew*timew);
574   //sigmaZ2 = sigmaZ2;// - 1./12;
575   
576   if(sigmaY2 <= 0 || sigmaZ2 <= 0)
577     {
578       cerr<<"AliL3Compress::CreateDigits() : Wrong sigmas : "<<sigmaY2<<" "<<sigmaZ2;
579       cerr<<" on row "<<row<<" pad "<<pad<<" time "<<time<<endl;
580       return;
581     }
582   
583   TRandom *random = new TRandom();
584   
585   Int_t entries=1000;
586   TH1F *hist1 = new TH1F("hist1","",AliL3Transform::GetNPads(row),0,AliL3Transform::GetNPads(row)-1);
587   TH1F *hist2 = new TH1F("hist2","",AliL3Transform::GetNTimeBins(),0,AliL3Transform::GetNTimeBins()-1);
588   TH2F *hist3 = new TH2F("hist3","",AliL3Transform::GetNPads(row),0,AliL3Transform::GetNPads(row)-1,AliL3Transform::GetNTimeBins(),0,AliL3Transform::GetNTimeBins()-1);
589
590   
591   
592   //Create the distributions in pad and time:
593   for(Int_t i=0; i<entries; i++)
594     {
595       hist1->Fill(random->Gaus(pad,sqrt(sigmaY2)));
596       hist2->Fill(random->Gaus(time,sqrt(sigmaZ2)));
597     }
598     
599   //Create the cluster:
600   Int_t bin1,bin2;
601   Double_t content1,content2,dpad,dtime;
602   for(Int_t i=0; i<hist1->GetEntries(); i++)
603     {
604       bin1 = hist1->GetBin(i);
605       content1 = hist1->GetBinContent(bin1);
606       if((Int_t)content1==0) continue;
607       content1 = charge*content1/entries;
608       dpad = hist1->GetBinCenter(bin1);
609       for(Int_t j=0; j<hist2->GetEntries(); j++)
610         {
611           bin2 = hist2->GetBin(j);
612           content2 = hist2->GetBinContent(bin2);
613           if((Int_t)content2==0) continue;
614           content2 = content1*content2/entries;
615           dtime = hist2->GetBinCenter(bin2);
616           hist3->Fill(dpad,dtime,content2);
617         }
618     }
619   
620   Int_t local_dig=0;
621   //Fill it into the digit array:
622   for(Int_t i=0; i<hist3->GetNbinsX(); i++)
623     {
624       for(Int_t j=0; j<hist3->GetNbinsY(); j++)
625         {
626           bin1 = hist3->GetBin(i,j);
627           content1 = hist3->GetBinContent(bin1);
628           if((Int_t)content1 < 3) continue;
629           if(content1 >= 1024)
630             content1 = 1023;
631           if(fNDigits >= fMaxDigits)
632             {
633               cerr<<"AliL3Compress::CreateDigits() : Array index out of range : "<<fNDigits<<endl;
634               return;
635             }
636           fDigits[fNDigits].fCharge=(Int_t)content1;
637           fDigits[fNDigits].fRow = row;
638           fDigits[fNDigits].fPad = (Int_t)hist3->GetXaxis()->GetBinCenter(i);
639           fDigits[fNDigits].fTime = (Int_t)hist3->GetYaxis()->GetBinCenter(j);
640           fDPt[fNDigits] = &fDigits[fNDigits];
641           fNDigits++;
642           local_dig++;
643         }
644     }
645   //if(local_dig < 5)
646   //  cout<<"Small cluster "<<local_dig<<" pad "<<(Int_t)fDigits[fNDigits-1].fPad<<" time "<<(Int_t)fDigits[fNDigits-1].fTime<<endl;
647   
648   delete random;
649   delete hist1;
650   delete hist2;
651   delete hist3;
652 }
653
654 void AliL3Compress::PrintCompRatio()
655 {
656   Char_t fname[100];
657   sprintf(fname,"%s/remains_%d_%d.raw",fPath,fSlice,fPatch);
658   AliL3MemHandler *mem = new AliL3MemHandler();
659   if(!mem->SetBinaryInput(fname))
660     {
661       cerr<<"AliL3Compress::PrintCompRatio(): Error opening file: "<<fname<<endl;
662       return;
663     } 
664   UInt_t ndigits;
665   AliL3DigitRowData *rowPt = (AliL3DigitRowData*)mem->CompBinary2Memory(ndigits);
666   mem->CloseBinaryInput();
667
668   Int_t digit_counter=0;
669   for(Int_t i=AliL3Transform::GetFirstRow(fPatch); i<=AliL3Transform::GetLastRow(fPatch); i++)
670     {
671       digit_counter += rowPt->fNDigit;
672       mem->UpdateRowPointer(rowPt);
673     }
674   delete mem;
675   
676   sprintf(fname,"%s/tracks_c_%d_%d.raw",fPath,fSlice,fPatch);
677   FILE *file1 = fopen(fname,"r");
678   if(!file1)
679     {
680       cerr<<"AliL3Compress::PrintCompRatio(): Error opening file: "<<fname<<endl;
681       return;
682     }
683   fseek(file1,0,SEEK_END);
684   UInt_t filesize1 = (UInt_t)ftell(file1);
685   fclose(file1);
686   
687   sprintf(fname,"%s/digits_%d_%d.raw",fPath,fSlice,fPatch);
688   FILE *file2 = fopen(fname,"r");
689   if(!file2)
690     {
691       cerr<<"AliL3Compress::PrintCompRatio(): Error opening file: "<<fname<<endl;
692       return;
693     }
694   fseek(file2,0,SEEK_END);
695   UInt_t filesize2 = (UInt_t)ftell(file2);
696   fclose(file2);
697   
698   cout<<"----------------------"<<endl;
699   cout<<"Original file size   : "<<filesize2<<endl;
700   cout<<"Compressed file size : "<<filesize1<<endl;
701   cout<<"Remaining digits     : "<<digit_counter<<endl;
702   cout<<"Compression ratio    : "<<(Float_t)(filesize1 + (10*digit_counter)/8)/(Float_t)(filesize2)<<endl;
703   
704 }
705
706 void AliL3Compress::QSort(AliL3RandomDigitData **a, Int_t first, Int_t last)
707 {
708   
709   // Sort array of AliL3RandomDigitData pointers using a quicksort algorithm.
710   // Uses CompareDigits() to compare objects.
711   // Thanks to Root!
712   
713   static AliL3RandomDigitData *tmp;
714   static int i;           // "static" to save stack space
715   int j;
716   
717   while (last - first > 1) {
718     i = first;
719     j = last;
720     for (;;) {
721       while (++i < last && CompareDigits(a[i], a[first]) < 0)
722         ;
723       while (--j > first && CompareDigits(a[j], a[first]) > 0)
724         ;
725       if (i >= j)
726         break;
727       
728       tmp  = a[i];
729       a[i] = a[j];
730       a[j] = tmp;
731     }
732     if (j == first) {
733       ++first;
734       continue;
735     }
736     tmp = a[first];
737     a[first] = a[j];
738     a[j] = tmp;
739     if (j - first < last - (j + 1)) {
740       QSort(a, first, j);
741       first = j + 1;   // QSort(j + 1, last);
742     } else {
743       QSort(a, j + 1, last);
744       last = j;        // QSort(first, j);
745     }
746   }
747 }
748
749 void AliL3Compress::WriteRootFile(Char_t *newrootfile)
750 {
751 #ifdef use_aliroot
752   Char_t fname[100];
753   AliL3MemHandler *mem = new AliL3MemHandler();
754   sprintf(fname,"%s/restored_%d_%d.raw",fPath,fSlice,fPatch);
755   mem->SetBinaryInput(fname);
756   UInt_t ndigits;
757   AliL3DigitRowData *rowPt = (AliL3DigitRowData*)mem->CompBinary2Memory(ndigits);
758   mem->CloseBinaryInput();
759   
760   sprintf(fname,"%s/digitfile",fPath);
761   
762   AliL3FileHandler *file = new AliL3FileHandler();
763   if(!file->SetAliInput(fname))
764     {
765       cerr<<"AliL3Compress::WriteRootFile() : Error opening file: "<<fname<<endl;
766       return;
767     }
768   const Int_t rows[2] = {AliL3Transform::GetFirstRow(fPatch),AliL3Transform::GetLastRow(fPatch)};
769   file->Init(fSlice,fPatch,rows);
770   file->AliDigits2RootFile(rowPt,newrootfile);
771   file->CloseAliInput();
772
773   delete mem;
774   delete file;
775 #endif
776
777   return;
778 }