]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDDLRawData.cxx
Fix for transient fSDigits, AliITSRawStream classes adapted to changed AliRawReader...
[u/mrichter/AliRoot.git] / TPC / AliTPCDDLRawData.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 /* $Id$ */
16
17
18 //This class conteins all the methods to create raw data 
19 //as par a given DDL.
20 //It produces DDL with both compressed and uncompressed format.
21 //For compression we use the optimized table wich needs 
22 //to be provided.
23
24 #include <TObjArray.h>
25 #include <Riostream.h>
26 #include <stdio.h>
27 #include <stdlib.h>
28 #include "AliTPCCompression.h"
29 #include "AliTPCBuffer160.h"
30 #include "AliTPCDDLRawData.h"
31
32 ClassImp(AliTPCDDLRawData)
33 ////////////////////////////////////////////////////////////////////////////////////////
34
35 AliTPCDDLRawData::AliTPCDDLRawData(const AliTPCDDLRawData &source){
36   // Copy Constructor
37   fVerbose=source.fVerbose;
38   return;
39 }
40
41 AliTPCDDLRawData& AliTPCDDLRawData::operator=(const AliTPCDDLRawData &source){
42   //Assigment operator
43   fVerbose=source.fVerbose;
44   return *this;
45 }
46
47
48 ////////////////////////////////////////////////////////////////////////////
49 void AliTPCDDLRawData::RawData(Int_t LDCsNumber,Int_t EventNumber){
50   //Raw data slides generation
51   //Number of DDL=2*36+4*36=216
52   //2 DDL for each inner sector
53   //4 DDL for each outer sector
54   Int_t ddlPerFile=216/LDCsNumber;
55   Int_t offset=1;
56   if (216%LDCsNumber) ddlPerFile++;
57   cout<<"Number of DDL per slide: "<<ddlPerFile<<endl;
58   ifstream f;
59 #ifndef __DECCXX
60   f.open("AliTPCDDL.dat",ios::binary);
61 #else
62   f.open("AliTPCDDL.dat");
63 #endif
64   if(!f){cout<<"File doesn't exist !!"<<endl;return;}
65   struct DataPad{
66     Int_t Sec;
67     Int_t SubSec;
68     Int_t Row;
69     Int_t Pad;
70     Int_t Dig;
71     Int_t Time;
72   };
73   DataPad data;
74
75   //AliTPCBuffer160 is used in write mode to generate AltroFormat.dat file
76   Int_t sliceNumber=1;
77   char  filename[15];
78   sprintf(filename,"Ev%dTPCslice%d",EventNumber,sliceNumber); 
79   cout<<"   Creating "<<filename<<endl;
80   AliTPCBuffer160 *buffer=new AliTPCBuffer160(filename,1);
81
82   ULong_t count=0;
83   Int_t pSecNumber=-1;  //Previous Sector number
84   Int_t pRowNumber=-1;  //Previous Row number  
85   Int_t pPadNumber=-1;  //Previous Pad number
86   Int_t pTimeBin=-1;    //Previous Time-Bin
87   Int_t pSubSector=-1;  //Previous Sub Sector
88   Int_t bunchLength=0;
89   Int_t countDDL=0;
90   Int_t nwords=0;
91   ULong_t numPackets=0;
92   while (f.read((char*)(&data),sizeof(data))){
93     count++;
94     if (pPadNumber==-1){
95       pSecNumber=data.Sec;
96       pRowNumber=data.Row;
97       pPadNumber=data.Pad;
98       pTimeBin=data.Time;
99       pSubSector=data.SubSec;
100       //size magic word sector number sub-sector number 0 for TPC 0 for uncompressed
101       buffer->WriteMiniHeader(0,pSecNumber,pSubSector,0,0);//Dummy;
102       bunchLength=1;
103       buffer->FillBuffer(data.Dig-offset);
104       nwords++;
105     }//end if
106     else{
107       if ( (data.Time==(pTimeBin+1)) &&
108            (pPadNumber==data.Pad) &&
109            (pRowNumber==data.Row) &&
110            (pSecNumber==data.Sec)){
111         bunchLength++;
112       }//end if
113       else{
114         buffer->FillBuffer(pTimeBin);
115         buffer->FillBuffer(bunchLength+2);
116         nwords+=2;
117         if ((pPadNumber!=data.Pad)||(pRowNumber!=data.Row)||(pSecNumber!=data.Sec)){
118           //Trailer is formatted and inserted!!
119           buffer->WriteTrailer(nwords,pPadNumber,pRowNumber,pSecNumber);
120           numPackets++;
121           nwords=0;
122
123           if(pSubSector!=data.SubSec){
124             countDDL++;
125             if(countDDL==ddlPerFile){
126               //size magic word sector number sub-sector number 0 for TPC 0 for uncompressed
127               buffer->Flush();
128               buffer->WriteMiniHeader(1,pSecNumber,pSubSector,0,0);
129               //cout<<"Mini header for DDL:"<<PSecNumber<<" Sub-sec:"<<PSubSector<<endl;
130               delete buffer;
131               sliceNumber++;
132               sprintf(filename,"Ev%dTPCslice%d",EventNumber,sliceNumber);
133               cout<<"   Creating "<<filename<<endl;
134               buffer=new AliTPCBuffer160(filename,1);
135               buffer->WriteMiniHeader(0,data.Sec,data.SubSec,0,0);//Dummy;
136               countDDL=0;
137             }//end if
138             else{
139               buffer->Flush();
140               buffer->WriteMiniHeader(1,pSecNumber,pSubSector,0,0);
141               buffer->WriteMiniHeader(0,data.Sec,data.SubSec,0,0);//Dummy;
142             }
143             pSubSector=data.SubSec;
144           }//end if
145         }//end if
146         
147         bunchLength=1;
148         pPadNumber=data.Pad;
149         pRowNumber=data.Row;
150         pSecNumber=data.Sec;
151       }//end else
152       pTimeBin=data.Time;
153       buffer->FillBuffer(data.Dig-offset);
154       nwords++;
155     }//end else
156   }//end while
157   buffer->FillBuffer(pTimeBin);
158   buffer->FillBuffer(bunchLength+2);
159   nwords+=2;
160   buffer->WriteTrailer(nwords,pPadNumber,pRowNumber,pSecNumber);
161   //write the  M.H.
162   buffer->Flush();
163   buffer->WriteMiniHeader(1,pSecNumber,pSubSector,0,0);
164   //cout<<"Mini header for D D L:"<<pSecNumber<<" Sub-sec:"<<pSubSector<<endl;
165   delete buffer;
166   cout<<"Number of digits: "<<count<<endl;
167   f.close();
168   return;
169 }
170 ////////////////////////////////////////////////////////////////////////////
171 ////////////////////////////////////////////////////////////////////////////
172
173
174 Int_t AliTPCDDLRawData::RawDataCompDecompress(Int_t LDCsNumber,Int_t EventNumber,Int_t Comp){
175   //This method is used to compress and decompress the slides
176   static const Int_t kNumTables=5;
177   char filename[20];
178   char dest[20];
179   fstream f;
180   ULong_t size=0;
181   //Int_t MagicWord,DDLNumber,SecNumber,SubSector,Detector;
182   Int_t flag=0;
183   for(Int_t i=1;i<=LDCsNumber;i++){
184     if(!Comp){
185       sprintf(filename,"Ev%dTPCslice%d",EventNumber,i);
186       sprintf(dest,"Ev%dTPCslice%d.comp",EventNumber,i);
187     }
188     else{
189       sprintf(filename,"Ev%dTPCslice%d.comp",EventNumber,i);
190       sprintf(dest,"Ev%dTPCslice%d.decomp",EventNumber,i);
191     }
192 #ifndef __DECCXX
193     f.open(filename,ios::binary|ios::in);
194 #else
195     f.open(filename,ios::in);
196 #endif
197     if(!f){cout<<"BE CAREFUL!! There isn't enough data to generate "<<LDCsNumber<<" slices"<<endl;break;}
198     if (fVerbose)
199       cout<<filename<<"  "<<dest<<endl;
200     ofstream fdest;
201 #ifndef __DECCXX
202     fdest.open(dest,ios::binary);
203 #else
204     fdest.open(dest);
205 #endif
206     //loop over the DDL block 
207     //Each block contains a Mini Header followed by raw data (ALTRO FORMAT)
208     //The number of block is ceil(216/LDCsNumber)
209     ULong_t miniHeader[3];
210     //here the Mini Header is read
211     while( (f.read((char*)(miniHeader),sizeof(ULong_t)*3)) ){
212       size=miniHeader[0];
213       // cout<<"Data size:"<<size<<endl;
214       //Int_t dim=sizeof(ULong_t)+sizeof(Int_t)*5;
215       //cout<<" Sec "<<SecNumber<<" SubSector "<<SubSector<<" size "<<size<<endl;
216       //open the temporay File
217       ofstream fo;
218       char temp[15]="TempFile";
219 #ifndef __DECCXX
220       fo.open(temp,ios::binary);
221 #else
222       fo.open(temp);
223 #endif
224       Int_t car=0;
225       for(ULong_t j=0;j<size;j++){
226         f.read((char*)(&car),1);
227         fo.write((char*)(&car),1);
228       }//end for
229       fo.close();
230       //The temp file is compressed or decompressed
231       AliTPCCompression *util = new AliTPCCompression();
232       util->SetVerbose(0);
233       if(!Comp){
234         util->CompressDataOptTables(kNumTables,temp,"TempCompDecomp");
235       }
236       else
237         util->DecompressDataOptTables(kNumTables,temp,"TempCompDecomp");
238       delete util;
239       //the temp compressed file is open and copied to the final file fdest
240       ifstream fi;
241 #ifndef __DECCXX
242       fi.open("TempCompDecomp",ios::binary);
243 #else
244       fi.open("TempCompDecomp");
245 #endif
246       fi.seekg(0,ios::end);
247       size=fi.tellg();
248       fi.seekg(0);
249       //The Mini Header is updated (size and Compressed flag) 
250       //and written into the output file
251       miniHeader[0]=size;
252       if(!Comp)
253         flag=1;
254       else
255         flag=0;
256       ULong_t aux=0x0;
257       flag<<=8;
258       aux|=flag;
259       miniHeader[2]=miniHeader[2]|aux;
260       fdest.write((char*)(miniHeader),sizeof(ULong_t)*3);
261       //The compressem temp file is copied into the output file fdest
262       for(ULong_t j=0;j<size;j++){
263         fi.read((char*)(&car),1);
264         fdest.write((char*)(&car),1);
265       }//end for
266       fi.close();
267     }//end while
268     f.clear();
269     f.close();
270     fdest.close();
271     remove("TempFile");
272     remove("TempCompDecomp");
273   }//end for
274   return 0;
275 }
276
277 /////////////////////////////////////////////////////////////////////////////////
278 void AliTPCDDLRawData::RawDataAltro()const{
279   //This method is used to build the Altro format from AliTPCDDL.dat
280   //It is used to debug the code and creates the tables used in the compresseion phase
281   Int_t offset=1;
282   ifstream f;
283 #ifndef __DECCXX
284   f.open("AliTPCDDL.dat",ios::binary);
285 #else
286   f.open("AliTPCDDL.dat");
287 #endif
288   if(!f){cout<<"File doesn't exist !!"<<endl;return;}
289   struct DataPad{
290     Int_t Sec;
291     Int_t SubSec;
292     Int_t Row;
293     Int_t Pad;
294     Int_t Dig;
295     Int_t Time;
296   };
297   DataPad data;
298
299   //AliTPCBuffer160 is used in write mode to generate AltroFormat.dat file
300   char  filename[30]="AltroFormatDDL.dat";
301   cout<<"   Creating "<<filename<<endl;
302   AliTPCBuffer160 *buffer=new AliTPCBuffer160(filename,1);
303
304   ULong_t count=0;
305   Int_t pSecNumber=-1;  //Previous Sector number
306   Int_t pRowNumber=-1;  //Previous Row number  
307   Int_t pPadNumber=-1;  //Previous Pad number
308   Int_t pTimeBin=-1;    //Previous Time-Bin
309   Int_t bunchLength=0;
310   Int_t nwords=0;
311   ULong_t numPackets=0;
312   while (f.read((char*)(&data),sizeof(data))){
313     count++;
314     if (pPadNumber==-1){
315       pSecNumber=data.Sec;
316       pRowNumber=data.Row;
317       pPadNumber=data.Pad;
318       pTimeBin=data.Time;
319       bunchLength=1;
320       buffer->FillBuffer(data.Dig-offset);
321       nwords++;
322     }//end if
323     else{
324       if ( (data.Time==(pTimeBin+1)) &&
325            (pPadNumber==data.Pad) &&
326            (pRowNumber==data.Row) &&
327            (pSecNumber==data.Sec)){
328         bunchLength++;
329       }//end if
330       else{
331         buffer->FillBuffer(pTimeBin);
332         buffer->FillBuffer(bunchLength+2);
333         nwords+=2;
334         if ((pPadNumber!=data.Pad)||(pRowNumber!=data.Row)||(pSecNumber!=data.Sec)){
335           //Trailer is formatted and inserted!!
336           buffer->WriteTrailer(nwords,pPadNumber,pRowNumber,pSecNumber);
337           numPackets++;
338           nwords=0;
339         }//end if
340         
341         bunchLength=1;
342         pPadNumber=data.Pad;
343         pRowNumber=data.Row;
344         pSecNumber=data.Sec;
345       }//end else
346       pTimeBin=data.Time;
347       buffer->FillBuffer(data.Dig-offset);
348       nwords++;
349     }//end else
350   }//end while
351   buffer->FillBuffer(pTimeBin);
352   buffer->FillBuffer(bunchLength+2);
353   nwords+=2;
354   buffer->WriteTrailer(nwords,pPadNumber,pRowNumber,pSecNumber);
355   delete buffer;
356   cout<<"Number of digits: "<<count<<endl;
357   f.close(); 
358   return;
359 }
360
361 /////////////////////////////////////////////////////////////////////////
362 void AliTPCDDLRawData::RawDataAltroDecode(Int_t LDCsNumber,Int_t EventNumber,Int_t Comp){
363   //This method merges the slides in only one file removing at the same 
364   //time all the mini headers. The file so obtained must be Altro format
365   //complaiant.
366   //It is used mainly in the debugging phase 
367   char filename[15];
368   char dest[30];
369   fstream f;
370   if(!Comp)
371     sprintf(dest,"AltroDDLRecomposed.dat");
372   else
373     sprintf(dest,"AltroDDLRecomposedDec.dat");
374   ofstream fdest;
375
376 #ifndef __DECCXX
377   fdest.open(dest,ios::binary);
378 #else
379   fdest.open(dest);
380 #endif
381   ULong_t size=0;
382   //Int_t MagicWord,DDLNumber,SecNumber,SubSector,Detector,flag=0;
383   for(Int_t i=1;i<=LDCsNumber;i++){
384     if(!Comp){
385       sprintf(filename,"Ev%dTPCslice%d",EventNumber,i);  
386     }
387     else{
388       sprintf(filename,"Ev%dTPCslice%d.decomp",EventNumber,i);  
389     }
390 #ifndef __DECCXX
391     f.open(filename,ios::binary|ios::in);
392 #else
393     f.open(filename,ios::in);
394 #endif
395     if(!f){cout<<"BE CAREFUL!! There isn't enough data to generate "<<LDCsNumber<<" slices"<<endl;break;}
396     //loop over the DDL block 
397     //Each block contains a Mini Header followed by raw data (ALTRO FORMAT)
398     //The number of block is ceil(216/LDCsNumber)
399     ULong_t miniHeader[3];
400     //here the Mini Header is read
401     //cout<<filename<<endl;
402     while( (f.read((char*)(miniHeader),sizeof(ULong_t)*3)) ){
403       Int_t car=0;
404       size=miniHeader[0];
405       for(ULong_t j=0;j<size;j++){
406         f.read((char*)(&car),1);
407         fdest.write((char*)(&car),1);
408       }//end for
409     }//end while
410     f.clear();
411     f.close();
412   }//end for
413   fdest.close();
414   return;
415 }
416