]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSRaw2Digits.cxx
Polishing to fulfill coding conventions
[u/mrichter/AliRoot.git] / PHOS / AliPHOSRaw2Digits.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 // Class designed to convert raw data to aliroot format. We assume, that
20 // prototype is situated in the center of 3 PHOS module and convert prototype
21 // outpur to AliPHOSDigits. In addition, we fill branch of TreeE with 
22 // AliPHOSBeamTestEvent, contaning description of event(triggers etc).
23 // Note, that one byte per channel in raw data is transvormed to class 
24 // AliPHOSDigit, so finale zise increase ~ 100 times. So output can be split 
25 // into peases of reasonable size: each file can not contain more than 
26 // fMaxPerFile: if there are still events in raw file, then new directory 
27 // is created and header+digits files are written to it.
28 // 
29 // Use Case:
30 //   AliPHOSRaw2Digits * r = new AliPHOSRaw2Digits("path/input.file") ;
31 //                   //note, that it can be gzipped file!
32 //   //Set position of the target in the given run.
33 //   //Z axis along beam direction, from target to prototype (0-surface of prototype)
34 //   //X axis along columns of prototype (0-center of prototype)
35 //   //Y axis along raws of prototype    (0-center of prototype)
36 //   Double_t pos[3]={0,0,-452.} ;
37 //   r->SetTargetPosition(pos) ;
38 //   //Read/create connection Table:
39 //   TFile f("ConTableDB.root") ;
40 //   AliPHOSConTableDB * cdb = f.Get("AliPHOSConTableDB") ;
41 //   f.Close() ;
42 //   r->SetConTableDB(cdb) ;
43 //   r->ExecuteTask() ;
44 //
45 // As a result files galice.root and PHOS.Digits.root should be produced in 
46 // current dir, and, possibly, dirs 1,2,3... each with galice.root and PHOS.Digits.root,
47 // where the rest of data are written. 
48 //
49 /*-- Author: Maxim Volkov (RRC KI)
50              Dmitri Peressounko (RRC KI & SUBATECH)
51              Yuri Kharlov (IHEP & SUBATECH)     */
52
53 //////////////////////////////////////////////////////////////////////////////
54
55 // --- ROOT system ---
56 #include "TClonesArray.h"
57 #include "TFile.h"
58 #include "TTree.h"
59 #include "TSystem.h"
60 #include "Bytes.h"
61
62 // --- Standard library ---
63 #include <sys/stat.h>
64 #include <fcntl.h>
65 #include <unistd.h>
66 #include <ctype.h>
67 #include <stdio.h>
68 #include <stdlib.h>
69 #include <unistd.h>
70 #include <string.h>
71 #include <netinet/in.h>
72
73 // --- AliRoot header files ---
74 #include "AliPHOSDigit.h"
75 #include "AliPHOSConTableDB.h"
76 #include "AliPHOSBeamTestEvent.h"
77 #include "AliPHOSRaw2Digits.h"
78 #include "AliPHOSv1.h"
79 #include "../EVGEN/AliGenBox.h"
80 #include "AliRun.h"
81
82 ClassImp(AliPHOSRaw2Digits)
83   
84   
85 //____________________________________________________________________________ 
86   AliPHOSRaw2Digits::AliPHOSRaw2Digits():TTask() 
87 {
88 //As one can easily see, this is constructor.
89   fInName="";  
90   fMK1 = 0x0123CDEF ;
91   fMK2 = 0x80708070 ;
92   fMK3 = 0x4321ABCD ;
93   fMK4 = 0x80618061 ;
94   fCKW = 0x4640E400 ;
95   fDebug = kFALSE;             //  Debug flag
96   fToSplit = kFALSE ;
97   fIsInitialized = kFALSE ;
98   fTarget[0] = 0 ;
99   fTarget[1] = 0 ;
100   fTarget[2] = 0 ;
101   fDigits = 0 ;
102   fHeaderFile = 0 ;
103   fDigitsFile = 0 ;
104   fPHOSHeader =0 ;
105   fMaxPerFile = 20000 ;  //Maximal number of events in root file.
106   fEvent = 0 ;
107   fctdb = 0;
108 }
109 //____________________________________________________________________________ 
110   AliPHOSRaw2Digits::AliPHOSRaw2Digits(const char * filename,Bool_t toSplit):TTask("Default","") 
111 {
112 //this constructor should be normally used. Parameters: imput file and should we produce output in split mode.
113   fInName=filename;
114   fToSplit = toSplit ;
115   TString outname("") ;
116   if(fToSplit)
117     outname = "galice.root" ;
118   else{
119     outname =fInName ;
120     outname.ToLower() ;
121     outname.ReplaceAll(".fz",".root") ;
122     outname.ReplaceAll(".gz","") ;
123   }
124   SetTitle(outname) ;
125
126   fMK1 = 0x0123CDEF ;
127   fMK2 = 0x80708070 ;
128   fMK3 = 0x4321ABCD ;
129   fMK4 = 0x80618061 ;
130   fCKW = 0x4640E400 ;
131   fDebug = kFALSE;             //  Debug flag
132   fIsInitialized = kFALSE ;
133   fTarget[0] = 0 ;
134   fTarget[1] = 0 ;
135   fTarget[2] = 0 ;
136   fDigits = 0 ;
137   fPHOSHeader =0 ;
138   fHeaderFile = 0 ;
139   fDigitsFile = 0 ;
140   fMaxPerFile = 20000 ;
141   fEvent = 0 ;
142   fctdb = 0;
143 }
144 //____________________________________________________________________________ 
145 AliPHOSRaw2Digits::~AliPHOSRaw2Digits()
146 {
147 //destructor
148   if(fPHOSHeader)
149     fPHOSHeader->Delete() ;
150   if(fDigits){
151     fDigits->Delete() ;
152     delete fDigits ;
153   }
154   
155 }
156 //____________________________________________________________________________ 
157 void AliPHOSRaw2Digits::Exec(Option_t * option){
158   //This is steering method performing all the conversion
159
160   if(!fIsInitialized) //need initialization
161     if(!Init())       //failed to initialize
162       return ;
163
164   ProcessRawFile() ;
165
166
167 //____________________________________________________________________________ 
168 Bool_t AliPHOSRaw2Digits::Init(void){
169   //Makes initialization of contaniers
170
171   if(fIsInitialized)
172     return kTRUE;
173
174   //Make container for digits
175   fDigits = new TClonesArray("AliPHOSDigit",1000) ;
176   fPHOSHeader = new  AliPHOSBeamTestEvent() ;
177   fIsInitialized = kTRUE ;
178   return StartRootFiles() ;
179
180 }
181 //____________________________________________________________________________ 
182 Bool_t AliPHOSRaw2Digits::StartRootFiles(void ){
183   //Create PHOS geometry, sets magnetic field to zero, 
184   //create Generator - to store target position, 
185   //opens out file, creates TreeE 
186
187   //create gAlice if nececcary
188   if(!gAlice)
189     new AliRun("gAlice","The ALICE Off-line Simulation Framework") ;
190
191   //Create PHOS
192   if(!gAlice->GetModule("PHOS"))
193     new AliPHOSv1("PHOS","GPS2") ;
194
195   //Set Magnetic field
196   gAlice->SetField(0,2);  
197
198   //Set positin of the virtex
199   AliGenerator * gener = gAlice->Generator() ; 
200   if(!gener)    
201     gener = new AliGenBox(1);
202   Float_t ox = fTarget[1]; 
203   Float_t oy = fTarget[2]+460.; 
204   Float_t oz = fTarget[0];
205   gener->SetOrigin(ox, oy, oz);
206
207   //make directory 
208   Int_t nRootFile = (fEvent+1)/fMaxPerFile ;    
209   if(nRootFile){
210     char dname[20];
211     sprintf(dname,"%d",nRootFile) ;
212     if(gSystem->AccessPathName(dname)) //strange return: 0 if exists
213       if(gSystem->MakeDirectory(dname)!=0)
214         Fatal("StartRootFiles","Can not make directory %s \n",dname) ;
215     
216     if(!gSystem->ChangeDirectory(dname))
217       Fatal("StartRootFiles","Can not cd to %s\n",dname) ;
218   }
219
220   //  Create the output file
221   TString outname("") ;
222   if(strstr(GetTitle(),"root")){
223     outname=GetTitle();
224   }
225   else{
226     outname = fInName ;
227     outname.ToLower() ;
228     outname.ReplaceAll(".fz",".root") ;
229   }
230
231   fHeaderFile = new TFile(outname,"recreate");
232   fHeaderFile->SetCompressionLevel(2);
233   
234   // Create the Root Trees
235   gAlice->MakeTree("E") ;
236   
237   //Fill now TreeE
238   Int_t splitlevel = 0 ;
239   Int_t bufferSize = 32000 ;    
240   TBranch * headerBranch = gAlice->TreeE()->Branch("AliPHOSBeamTestEvent", 
241                                                    "AliPHOSBeamTestEvent", 
242                                                    &fPHOSHeader,bufferSize,splitlevel);
243   headerBranch->SetName("AliPHOSBeamTestEvent") ;
244
245   if(fToSplit){
246     fDigitsFile = new TFile("PHOS.Digits.root","recreate") ;
247     fDigitsFile->SetCompressionLevel(2) ;
248   }
249   return kTRUE ;
250 }
251 //____________________________________________________________________________ 
252 Bool_t AliPHOSRaw2Digits::CloseRootFiles(void ){
253   //cleans everething to start next root file
254   if(fHeaderFile){
255     printf("writing gAlice \n") ;
256     fHeaderFile->cd() ;
257     gAlice->Write(0,TObject::kOverwrite);
258     gAlice->TreeE()->Write(0,TObject::kOverwrite);
259   }
260
261   delete gAlice ;
262   
263   if(fHeaderFile){
264     fHeaderFile->Close() ;
265     delete fHeaderFile ;
266     fHeaderFile = 0;
267   }   
268         
269   if(fDigitsFile){
270     fDigitsFile->Close() ;
271     delete fDigitsFile ;
272     fDigitsFile = 0 ;
273   }
274   
275   Int_t nRootFile = (fEvent-1)/fMaxPerFile ;    
276   if(nRootFile){
277     if(!gSystem->ChangeDirectory("../")){
278      Fatal("CloseRootFile","Can not return to initial dir \n") ;      
279      return kFALSE ;
280     }
281   }
282   return kTRUE ;
283 }
284 //____________________________________________________________________________ 
285 Bool_t AliPHOSRaw2Digits::ProcessRawFile(){
286   
287   //Method opens zebra file and reads successively bytes from it,
288   //filling corresponding fields in header and digits.
289   
290   
291   fStatus= -3 ;
292   //First of all, open file and check if it is a zebra file
293   
294   char command[256];
295   sprintf(command,"zcat %s",fInName.Data());
296   FILE *dataFile = popen(command, "r");
297   if (dataFile == NULL) {
298     Warning("ProcessRawFile", " Cannot open file %s\n", fInName.Data() ) ;
299     perror(fInName.Data()) ;
300     fStatus = -1 ;
301     return kFALSE ;
302   }
303
304   // Check if byte ordering is little-endian 
305   UInt_t w = 0x12345678;
306   Int_t swapo = memcmp(&w, "\x78\x56\x34\x12", sizeof(UInt_t)) == 0;
307   if(fDebug)
308     Info("ProcessRawFile", "swapo=%f\n", swapo ) ;
309   
310   
311   UInt_t recBuf[300] ;
312
313   // Read physical record control words 
314   UInt_t nb = 8*sizeof(UInt_t);
315   Int_t n = fread(recBuf, nb, 1, dataFile);
316   if(static_cast<UInt_t>(n) != 1) {
317     if (n < 0 )
318       perror(fInName.Data());
319     else
320       Error("ProcessRawFile", "Could not read physical record control words" ) ;
321     fStatus = -2 ;
322     return kFALSE;
323   }
324   
325   if(fDebug)
326     Info("ProcessRawFile", "recbuf[0] = %d\n", recBuf[0] );
327   
328   // Check if it is a ZEBRA file and if the words are swapped 
329   UInt_t swapi = 0 ;
330   if (recBuf[0] != fMK1) {
331     Swab4(recBuf, &w, 1);
332     if (w != fMK1) {
333       Error("ProcessRawFile", "Does not look like a ZEBRA file\n" ) ;
334       pclose(dataFile) ;
335       fStatus = -2 ;
336       return kFALSE;
337     }
338     swapi=1 ;
339   }
340   
341   if(fDebug){
342     TString message ; 
343     message  = "        w = %f\n" ; 
344     message += "    swapi = %f\n" ; 
345     Info("ProcessRawFile", message.Data(), w, swapi ) ; 
346   }
347   
348   // Get number of words in physical record 
349   UInt_t  nwphr ;
350   if (swapi)
351     Swab4(&recBuf[4],&nwphr,1);
352   else 
353     nwphr = recBuf[4];
354   nwphr*=2; // 1998 -- Now we have 2 records 150 words each 
355
356
357   //start loop over data  
358   // Read rest of record 
359   nb = (nwphr-8)*sizeof(UInt_t);
360   n = fread(&recBuf[8], nb, 1, dataFile) ;
361   if (static_cast<UInt_t>(n) != 1) {
362     if (n < 0 ){
363       perror(fInName.Data());
364       fStatus = -2 ;
365       return kFALSE;
366     }
367   }
368   nb = nwphr *sizeof(UInt_t);
369
370   UInt_t userVector[16] ;
371   UInt_t zheader[12];    
372   UShort_t pattern ;
373   UShort_t scanning[32] ;
374   UShort_t charge[12];
375   UInt_t scaler[12]; 
376   UShort_t tdc2228[32];
377   
378   //read untill the end of file
379   fEvent=0 ;
380   while(1){
381
382     //    StartNewEvent() ;
383     fDigits->Delete() ;
384     if((fEvent%fMaxPerFile == 0) && fEvent ){
385       CloseRootFiles() ;
386       StartRootFiles() ;
387     }
388     gAlice->SetEvent(fEvent%fMaxPerFile) ;
389
390     //Set Beam Energy
391     fPHOSHeader->SetBeamEnergy(fBeamEnergy) ;
392           
393     Int_t i ;
394     for(i=0;i<16;i++)
395       userVector[i]=*(recBuf+21+i);
396     if(!swapi)
397       Swab4(userVector, userVector, 16);     
398     fPHOSHeader->SetUserVector(userVector) ;
399     
400     
401     // ZEBRA event header
402     for(i=0;i<12;i++)
403       zheader[i]=*(recBuf+47+i);
404     if(swapi)
405       Swab4(zheader, zheader, 12);
406     fPHOSHeader->SetHeader(zheader) ;
407     
408     // Swap input 
409     if (swapi)
410       Swab4(recBuf, recBuf, nwphr);
411     
412     /* Physical record control words */
413     UInt_t * recptr = recBuf;  //Pointer to current position
414
415     if(recptr[7] != 1) {
416       Error("ProcessRawFile", "Cannot handle fast blocks" ) ; 
417       fStatus = -2 ;
418       return kFALSE;
419     }    
420     recptr += 8;
421     
422     // Logical record control words   
423     UInt_t lrtyp = recptr[1];
424     if (lrtyp != 3) {
425       Error("ProcessRawFile", "Can not handle logical record type %d", lrtyp ) ;
426       fStatus = -2 ;
427       return kFALSE;
428     }
429     
430     recptr += 2;
431     if (recptr[0] != fCKW) {
432       Error("ProcessRawFile", "Bad check word" ) ;
433       fStatus = -2 ;
434       return kFALSE;
435     }
436     
437     UInt_t  nwuh = recptr[9];
438     recptr += 10+nwuh;
439     
440     // Bank system words 
441     UInt_t nd = recptr[8];            /* Number of data words */
442     recptr += 10;
443     
444     // Data words 
445     UInt_t evtno = recptr[2];                   /* Event number */
446     
447     if(fDebug)
448        Info("ProcessRawFile", "evtno= %d", evtno);
449     
450     UInt_t nh = recptr[4];                   /* Number of header words in data bank */
451     recptr += nh;
452     
453     // Unswap data from VME 
454     if (swapi)
455       Swab4(recptr, recptr, nd-nh-3);
456     
457     // Give buffer to monitor program 
458     //  UInt_t esize = nd-nh-3;
459     //    if (swapo)
460     //       Swab2(recptr, recptr, esize); 
461     // Two byte data are correct after this operation. 
462     //But we're in trouble if the data array contains 4 byte data!
463     
464     // From now on deal with VME data (MSB first, or network byte order).
465     
466     
467     // Fill the event with data from ADCs
468     UChar_t *byteptr=(UChar_t*)recptr;
469     
470     // Trigger bit register  
471     pattern=net2host(*(UShort_t*)byteptr);
472     fPHOSHeader->SetPattern(pattern) ;
473     byteptr+=sizeof(UShort_t);
474     
475     // Either peak ADCs, 10 modulesX8=80 channels, 
476     //or Kurchatov 64+2 channel ADC 
477     //(the rest of the channels padded with 0xffff) 
478     for(i=0;i<80;i++){
479       Int_t peak = static_cast<Int_t>(net2host(*(UShort_t*)byteptr));
480       //make digit
481       Int_t absID = fctdb->Raw2AbsId(i) ;
482       if(absID > 0)
483         new((*fDigits)[i])AliPHOSDigit(-1,absID,peak,0.,i) ;
484       if(fDebug){
485         if(peak>(UShort_t)1000)
486           Info("ProcessRawFile", "event= %d peak[%d] = %f", fEvent, i, peak);
487       }
488       byteptr+=sizeof(UShort_t);
489     }
490     
491     // Scanning ADCs, 4 modulesX8=32 channels
492     for(i=0;i<32;i++){
493       scanning[i]=net2host(*(UShort_t*)byteptr);
494       byteptr+=sizeof(UShort_t);
495     }
496     fPHOSHeader->SetScanning(scanning) ;
497     
498     // Charge ADCs, 1 moduleX12=12 channels
499     for(i=0;i<12;i++){
500       charge[i]=net2host(*(UShort_t*)byteptr);
501       byteptr+=sizeof(UShort_t);
502     }
503     fPHOSHeader->SetCharge(charge) ;
504     
505     // Scalers, 1 moduleX12=12 (4 byte) channels
506     for(i=0;i<12;i++){
507       scaler[i]=net2host(*(UInt_t*)byteptr);
508       byteptr+=sizeof(UInt_t);
509     }
510     fPHOSHeader->SetScaler(scaler) ;
511     
512     // LeCroy TDC 2228A, 4 moduleX8=32 channels
513     for(i=0;i<8;i++){
514       tdc2228[i]=net2host(*(UShort_t*)byteptr);
515       byteptr+=sizeof(UShort_t);
516     }
517     fPHOSHeader->SetTDC(tdc2228) ;
518
519     WriteDigits() ;
520     if(fDebug)
521       Info("ProcessRawFile", "event= %d written", fEvent) ;
522  
523     // Read next record 
524     UInt_t nb = nwphr *sizeof(UInt_t);
525     n = fread( recBuf, nb,1,dataFile);
526     if (n < 0 ){
527       perror(fInName);
528       fStatus = -2 ;
529       return kFALSE;
530     }
531     if (static_cast<UInt_t>(n) != 1) {
532       pclose(dataFile) ;
533       fStatus = 1 ;
534       return kTRUE ; //all read
535     }
536     fEvent++ ;
537   }
538   CloseRootFiles() ;
539   
540   fStatus = 1 ;  
541   return kTRUE ;  
542 }
543 //____________________________________________________________________________ 
544 void AliPHOSRaw2Digits::Swab4(void *from, void *to, size_t nwords)const{
545   // The function swaps 4 bytes: byte#3<-->byte#0, byte#2<-->byte#1 
546   register char *pf=static_cast<char*>(from) ;
547   register char *pt=static_cast<char*>(to) ;
548   register char c;
549   while (nwords-- > 0 ) {
550     c = pf[0];
551     pt[0] = pf[3];
552     pt[3] = c;
553     c = pf[1];
554     pt[1] = pf[2];
555     pt[2] = c;
556     pf += 4;
557     pt += 4;
558   }
559 }
560
561 //____________________________________________________________________________ 
562 void AliPHOSRaw2Digits::Swab2(void *from, void *to, size_t nwords)const
563
564 //The function swaps 2x2 bytes: byte#0<-->byte#1, byte#2<-->byte#3 
565   register char *pf=static_cast<char*>(from) ;
566   register char *pt=static_cast<char*>(to);
567   register char c;   
568   while (nwords-- > 0 ) {
569     c = pf[0];
570     pt[0] = pf[1];
571     pt[1] = c;
572     c = pf[2];
573     pt[2] = pf[3];
574     pt[3] = c;
575     pf += 4;
576     pt += 4;
577   }
578 }
579
580 //____________________________________________________________________________ 
581 void AliPHOSRaw2Digits::WriteDigits(void){
582   //In this method we create TreeD, write digits and Raw2Digits to it
583   // and write Header to TreeE. Finally we write TreeD to root file 
584   
585   //Start from Digits
586   fDigits->Sort() ;
587   fDigits->Expand(fDigits->GetEntriesFast()) ;
588   for(Int_t i=0;i<fDigits->GetEntriesFast(); i++)
589     static_cast<AliPHOSDigit*>(fDigits->At(i))->SetIndexInList(i) ;
590
591   char hname[30];
592   sprintf(hname,"TreeD%d",fEvent%fMaxPerFile);
593   TTree * treeD = new TTree(hname,"Digits");
594   //treeD->Write(0,TObject::kOverwrite);
595   
596   // -- create Digits branch
597   Int_t bufferSize = 32000 ;    
598   TBranch * digitsBranch = treeD->Branch("PHOS",&fDigits,bufferSize);
599   digitsBranch->SetTitle("Default");
600   
601   // -- Create Digitizer branch
602   Int_t splitlevel = 0 ;
603   const AliPHOSRaw2Digits * d = this ;
604   TBranch * digitizerBranch = treeD->Branch("AliPHOSRaw2Digits", 
605                                             "AliPHOSRaw2Digits", &d,bufferSize,splitlevel); 
606   digitizerBranch->SetTitle("Default");
607
608   if(fDigitsFile)
609     fDigitsFile->cd() ;
610   digitsBranch->Fill() ;
611   digitizerBranch->Fill() ; 
612   treeD->Write(0,TObject::kOverwrite);
613  
614   delete treeD ;
615
616   //Write header
617   fHeaderFile->cd() ;
618   gAlice->TreeE()->Fill();
619 }
620 //____________________________________________________________________________ 
621 void AliPHOSRaw2Digits::Print(Option_t * option)const{
622   //prints current configuration and status.
623
624   printf("----------AliPHOSRaw2Digits---------- \n") ;
625   printf("Current input  File: %s\n",fInName.Data()) ; 
626   printf("Current output File: %s\n", GetTitle()); 
627   printf("Events processes in the last file %d\n",fEvent) ; 
628   printf("Input file status\n") ;
629   switch (fStatus){
630   case 0: printf("`Have not processed yet'\n") ;
631     break ;
632   case 1: printf("`Processed normally'\n") ;
633     break ;
634   case -1: printf("`File not found'\n") ;
635     break ;
636   case -2: printf("`Error in reading'\n") ;
637     break ;
638   case -3: printf("'Interupted'\n") ;
639   default: ;
640   }
641   printf("Connection table: " );
642   if(fctdb)
643     printf("%s %s \n",fctdb->GetName(), fctdb->GetTitle() ) ; 
644   else
645     printf(" no DB \n" );
646
647 }