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