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