1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
21 /*-- Author: Maxim Volkov (RRC KI)
22 Dmitri Peressounko (RRC KI & SUBATECH)
23 Yuri Kharlov (IHEP & SUBATECH) */
25 //////////////////////////////////////////////////////////////////////////////
27 // --- ROOT system ---
28 #include "TClonesArray.h"
32 // --- Standard library ---
33 //#include <sys/types.h>
42 #include <netinet/in.h>
44 // --- AliRoot header files ---
45 #include "AliPHOSDigit.h"
46 #include "AliPHOSConTableDB.h"
47 #include "AliPHOSBeamTestEvent.h"
48 #include "AliPHOSRaw2Digits.h"
49 #include "AliPHOSv1.h"
50 #include "../EVGEN/AliGenBox.h"
53 ClassImp(AliPHOSRaw2Digits)
56 //____________________________________________________________________________
57 AliPHOSRaw2Digits::AliPHOSRaw2Digits():TTask()
65 fDebug = kFALSE; // Debug flag
66 fIsInitialized = kFALSE ;
75 //____________________________________________________________________________
76 AliPHOSRaw2Digits::AliPHOSRaw2Digits(const char * filename):TTask("Default","")
79 TString outname(fInName) ;
81 outname.ReplaceAll(".fz",".root") ;
82 outname.ReplaceAll(".gz","") ;
90 fDebug = kFALSE; // Debug flag
91 fIsInitialized = kFALSE ;
101 //____________________________________________________________________________
102 AliPHOSRaw2Digits::~AliPHOSRaw2Digits()
105 fPHOSHeader->Delete() ;
112 //____________________________________________________________________________
113 void AliPHOSRaw2Digits::Exec(Option_t * option){
114 //This is steering method performing all the conversion
116 if(!fIsInitialized) //need initialization
117 if(!Init()) //failed to initialize
124 //____________________________________________________________________________
125 Bool_t AliPHOSRaw2Digits::Init(void){
126 //Create PHOS geometry, sets magnetic field to zero,
127 //create Generator - to store target position,
128 //opens out file, creates TreeE and make initialization of contaniers
135 new AliPHOSv1("PHOS","GPS2") ;
138 gAlice->SetField(0,2);
140 //Set positin of the virtex
141 AliGenBox *gener = new AliGenBox(1);
142 Float_t ox = fTarget[1];
143 Float_t oy = fTarget[2]-460.;
144 Float_t oz = fTarget[0];
145 gener->SetOrigin(ox, oy, oz);
147 // Create the output file
148 TString outname("") ;
149 if(strstr(GetTitle(),"root")){
155 outname.ReplaceAll(".fz",".root") ;
157 TFile *rootfile = new TFile(outname,"recreate");
158 rootfile->SetCompressionLevel(2);
160 // Create the Root Trees
161 gAlice->MakeTree("E") ;
163 //Make container for digits
164 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
167 fPHOSHeader = new AliPHOSBeamTestEvent() ;
168 Int_t splitlevel = 0 ;
169 Int_t bufferSize = 32000 ;
170 TBranch * headerBranch = gAlice->TreeE()->Branch("AliPHOSBeamTestEvent",
171 "AliPHOSBeamTestEvent",
172 &fPHOSHeader,bufferSize,splitlevel);
173 headerBranch->SetName("AliPHOSBeamTestEvent") ;
175 fIsInitialized = kTRUE ;
178 //____________________________________________________________________________
179 Bool_t AliPHOSRaw2Digits::ProcessRawFile(){
181 //Method opens zebra file and reads successively bytes from it,
182 //filling corresponding fields in header and digits.
186 //First of all, open file and check if it is a zebra file
189 sprintf(command,"zcat %s",fInName.Data());
190 FILE *dataFile = popen(command, "r");
191 if (dataFile == NULL) {
192 Warning("ProcessRawFile", " Cannot open file %s\n", fInName.Data() ) ;
193 perror(fInName.Data()) ;
197 printf("Open pipe: %s\n",command);
199 // Check if byte ordering is little-endian
200 UInt_t w = 0x12345678;
201 Int_t swapo = memcmp(&w, "\x78\x56\x34\x12", sizeof(UInt_t)) == 0;
203 Info("ProcessRawFile", "swapo=%f\n", swapo ) ;
208 // Read physical record control words
209 UInt_t nb = 8*sizeof(UInt_t);
210 Int_t n = fread(recBuf, nb, 1, dataFile);
211 if(static_cast<UInt_t>(n) != 1) {
213 perror(fInName.Data());
215 Error("ProcessRawFile", "Could not read physical record control words" ) ;
221 Info("ProcessRawFile", "recbuf[0] = %d\n", recBuf[0] );
223 // Check if it is a ZEBRA file and if the words are swapped
225 if (recBuf[0] != fMK1) {
226 Swab4(recBuf, &w, 1);
228 Error("ProcessRawFile", "Does not look like a ZEBRA file\n" ) ;
238 message = " w = %f\n" ;
239 message += " swapi = %f\n" ;
240 Info("ProcessRawFile", message.Data(), w, swapi ) ;
243 // Get number of words in physical record
246 Swab4(&recBuf[4],&nwphr,1);
249 nwphr*=2; // 1998 -- Now we have 2 records 150 words each
252 //start loop over data
253 // Read rest of record
254 nb = (nwphr-8)*sizeof(UInt_t);
255 n = fread(&recBuf[8], nb, 1, dataFile) ;
256 if (static_cast<UInt_t>(n) != 1) {
258 perror(fInName.Data());
263 nb = nwphr *sizeof(UInt_t);
265 UInt_t userVector[16] ;
268 UShort_t scanning[32] ;
271 UShort_t tdc2228[32];
273 //read untill the end of file
279 gAlice->SetEvent(fEvent) ;
283 userVector[i]=*(recBuf+21+i);
285 Swab4(userVector, userVector, 16);
286 fPHOSHeader->SetUserVector(userVector) ;
289 // ZEBRA event header
291 zheader[i]=*(recBuf+47+i);
293 Swab4(zheader, zheader, 12);
294 fPHOSHeader->SetHeader(zheader) ;
298 Swab4(recBuf, recBuf, nwphr);
300 /* Physical record control words */
301 UInt_t * recptr = recBuf; //Pointer to current position
304 Error("ProcessRawFile", "Cannot handle fast blocks" ) ;
310 // Logical record control words
311 UInt_t lrtyp = recptr[1];
313 Error("ProcessRawFile", "Can not handle logical record type %d", lrtyp ) ;
319 if (recptr[0] != fCKW) {
320 Error("ProcessRawFile", "Bad check word" ) ;
325 UInt_t nwuh = recptr[9];
329 UInt_t nd = recptr[8]; /* Number of data words */
333 UInt_t evtno = recptr[2]; /* Event number */
336 Info("ProcessRawFile", "evtno= %d", evtno);
338 UInt_t nh = recptr[4]; /* Number of header words in data bank */
341 // Unswap data from VME
343 Swab4(recptr, recptr, nd-nh-3);
345 // Give buffer to monitor program
346 // UInt_t esize = nd-nh-3;
348 // Swab2(recptr, recptr, esize);
349 // Two byte data are correct after this operation.
350 //But we're in trouble if the data array contains 4 byte data!
352 // From now on deal with VME data (MSB first, or network byte order).
355 // Fill the event with data from ADCs
356 UChar_t *byteptr=(UChar_t*)recptr;
358 // Trigger bit register
359 pattern=ntohs(*(UShort_t*)byteptr);
360 fPHOSHeader->SetPattern(pattern) ;
361 byteptr+=sizeof(UShort_t);
363 // Either peak ADCs, 10 modulesX8=80 channels,
364 //or Kurchatov 64+2 channel ADC
365 //(the rest of the channels padded with 0xffff)
367 Int_t peak = static_cast<Int_t>(ntohs(*(UShort_t*)byteptr));
369 Int_t absID = fctdb->Raw2AbsId(i) ;
371 new((*fDigits)[i])AliPHOSDigit(-1,absID,peak,0.,i) ;
373 if(peak>(UShort_t)1000)
374 Info("ProcessRawFile", "event= %d peak[%d] = %f", fEvent, i, peak);
376 byteptr+=sizeof(UShort_t);
379 // Scanning ADCs, 4 modulesX8=32 channels
381 scanning[i]=ntohs(*(UShort_t*)byteptr);
382 byteptr+=sizeof(UShort_t);
384 fPHOSHeader->SetScanning(scanning) ;
386 // Charge ADCs, 1 moduleX12=12 channels
388 charge[i]=ntohs(*(UShort_t*)byteptr);
389 byteptr+=sizeof(UShort_t);
391 fPHOSHeader->SetCharge(charge) ;
393 // Scalers, 1 moduleX12=12 (4 byte) channels
395 scaler[i]=ntohl(*(UInt_t*)byteptr);
396 byteptr+=sizeof(UInt_t);
398 fPHOSHeader->SetScaler(scaler) ;
400 // LeCroy TDC 2228A, 4 moduleX8=32 channels
402 tdc2228[i]=ntohs(*(UShort_t*)byteptr);
403 byteptr+=sizeof(UShort_t);
405 fPHOSHeader->SetTDC(tdc2228) ;
409 Info("ProcessRawFile", "event= %d written", fEvent) ;
412 UInt_t nb = nwphr *sizeof(UInt_t);
413 n = fread( recBuf, nb,1,dataFile);
419 if (static_cast<UInt_t>(n) != 1) {
422 return kTRUE ; //all read
430 //____________________________________________________________________________
431 void AliPHOSRaw2Digits::Swab4(void *from, void *to, size_t nwords){
432 // The function swaps 4 bytes: byte#3<-->byte#0, byte#2<-->byte#1
433 register char *pf=static_cast<char*>(from) ;
434 register char *pt=static_cast<char*>(to) ;
436 while (nwords-- > 0 ) {
448 //____________________________________________________________________________
449 void AliPHOSRaw2Digits::Swab2(void *from, void *to, size_t nwords)
450 { //The function swaps 2x2 bytes: byte#0<-->byte#1, byte#2<-->byte#3
451 register char *pf=static_cast<char*>(from) ;
452 register char *pt=static_cast<char*>(to);
454 while (nwords-- > 0 ) {
466 //____________________________________________________________________________
467 void AliPHOSRaw2Digits::FinishRun(){
468 //Write geometry and header tree
469 gAlice->Write(0,TObject::kOverwrite);
470 gAlice->TreeE()->Write(0,TObject::kOverwrite);
473 //____________________________________________________________________________
474 void AliPHOSRaw2Digits::WriteDigits(void){
475 //In this method we create TreeD, write digits and Raw2Digits to it
476 // and write Header to TreeE. Finally we write TreeD to root file
480 fDigits->Expand(fDigits->GetEntriesFast()) ;
481 for(Int_t i=0;i<fDigits->GetEntriesFast(); i++)
482 static_cast<AliPHOSDigit*>(fDigits->At(i))->SetIndexInList(i) ;
485 sprintf(hname,"TreeD%d",fEvent);
486 TTree * treeD = new TTree(hname,"Digits");
487 //treeD->Write(0,TObject::kOverwrite);
489 // -- create Digits branch
490 Int_t bufferSize = 32000 ;
491 TBranch * digitsBranch = treeD->Branch("PHOS",&fDigits,bufferSize);
492 digitsBranch->SetTitle("Default");
494 // -- Create Digitizer branch
495 Int_t splitlevel = 0 ;
496 const AliPHOSRaw2Digits * d = this ;
497 TBranch * digitizerBranch = treeD->Branch("AliPHOSRaw2Digits",
498 "AliPHOSRaw2Digits", &d,bufferSize,splitlevel);
499 digitizerBranch->SetTitle("Default");
501 digitsBranch->Fill() ;
502 digitizerBranch->Fill() ;
503 treeD->Write(0,TObject::kOverwrite);
508 gAlice->TreeE()->Fill();
510 //____________________________________________________________________________
511 void AliPHOSRaw2Digits::Print(Option_t * option)const{
514 message = "----------AliPHOSRaw2Digits---------- \n" ;
515 message += "Input stream \n" ;
516 message += "Current input File: %s\n" ;
517 message += "Current output File: %s\n" ;
518 message += "Events processes in the last file %d\n" ;
519 message += "Input file status\n" ;
521 case 0: message += "`Have not processed yet'\n" ;
523 case 1: message += "`Processed normally'\n" ;
525 case -1: message += "`File not found'\n" ;
527 case -2: message += "`Error in reading'\n" ;
529 case -3: message += "'Interupted'\n" ;
532 message += "Connection table: " ;
534 message += "%s %s \n" ;
536 message += " no DB \n" ;
538 Info("Print", message.Data(),
542 fctdb->GetName(), fctdb->GetTitle() ) ;