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