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