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