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