]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSCalibrator.cxx
pyquen added.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSCalibrator.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
17 //_________________________________________________________________________
18 // Class to calculate calibration parameters  from beam tests etc.
19 //  First pass - one should calcuate pedestals, in the second pass - gains.
20 //
21 //     To calculate pedestals we scan pedestals events and fill histos for each
22 //     channel. Then in each histo we find maximum and fit it with Gaussian in the 
23 //     visinity of the maximum. Extracted mean and width of distribution are put into 
24 //     resulting histogram which cheks results for 'reasonability'.
25 //
26 //     To evaluate gains: scans beam events and calculate gain in the cristall with 
27 //     maximal energy deposition, assuming, that 90% of energy deposited in it. 
28 //     For each channel separate histogramm is filled. When scan is finished, 
29 //     these histograms are fitted with Gaussian and finds mean and width, which 
30 //     are put in final histogram. Finaly gains are checked for deviation from mean. 
31
32 //     Finally fills database.
33 //
34 // Use Case:
35 //       AliPHOSCalibrator * c = new AliPHOSCalibrator("path/galice.root") ;
36 //       c->AddRun("path2/galice.root") ; 
37 //       c->ScanPedestals();
38 //       c->CalculatePedestals();             
39 //       c->WritePedestals();
40 //       c->ScanGains() ;
41 //       c->CalculateGains() ;
42 //       c->WriteGains() ;
43 //
44 //*-- Author : D.Peressounko  (RRC KI) 
45 //////////////////////////////////////////////////////////////////////////////
46
47 // --- ROOT system ---
48 #include "TF1.h"
49 #include "TFile.h"
50 #include "TObjString.h"
51 #include "TROOT.h"
52
53 // --- Standard library ---
54
55 // --- AliRoot header files ---
56 #include "AliPHOSCalibrManager.h"
57 #include "AliPHOSCalibrationData.h"
58 #include "AliPHOSCalibrator.h"
59 #include "AliPHOSConTableDB.h"
60 #include "AliPHOSGetter.h"
61
62 ClassImp(AliPHOSCalibrator)
63
64
65 //____________________________________________________________________________ 
66   AliPHOSCalibrator::AliPHOSCalibrator():TTask("AliPHOSCalibrator","Default") 
67 {
68   //Default constuctor for root. Normally should not be used
69   fRunList=0 ;
70   fNch = 0 ;
71   fPedHistos = 0 ;
72   fGainHistos = 0 ;
73   fhPedestals = 0 ;
74   fhPedestalsWid = 0 ;
75   fctdb = 0 ;
76   fConTableDB = "Beamtest2002" ;
77   fConTableDBFile = "ConTableDB.root" ;
78 }
79 //____________________________________________________________________________ 
80 AliPHOSCalibrator::AliPHOSCalibrator(const char* file, const char* title):
81   TTask("AliPHOSCalibrator",title) 
82
83   //Constructor which should normally be used.
84   //file: path/galice.root  - header file
85   //title: branch name of PHOS reconstruction (e.g. "Default")
86  
87
88   fRunList = new TList() ;
89   fRunList->SetOwner() ;
90   fRunList->Add(new TObjString(file)) ;
91   fNch = 0 ;
92   fPedPat = 257 ;  //Patterns for different kind of events
93   fPulPat = 33 ;
94   fLEDPat = 129 ;
95   fWBPat  = 1027 ;
96   fNBPat  = 1029 ;
97
98   fNChan  = 100 ;  
99   fGainMax = 0.1 ;
100   fNGainBins= 100 ;
101   fAcceptCorr = 10 ;     //Maximal deviation from mean, considered as normal 
102
103   fGainAcceptCorr = 5 ;  //Factor for gain deviation
104   fPedHistos = 0 ;
105   fGainHistos = 0 ;
106   fhPedestals = 0 ;
107   fhPedestalsWid = 0 ;
108   fctdb = 0 ;
109   fConTableDB = "Beamtest2002" ;
110   fConTableDBFile = "ConTableDB.root" ;
111 }
112
113 //____________________________________________________________________________ 
114   AliPHOSCalibrator::~AliPHOSCalibrator()
115 {
116   // dtor
117   if(fPedHistos)
118     delete fPedHistos ;
119   if(fGainHistos)
120     delete fGainHistos ;
121   if(fhPedestals)
122     delete fhPedestals ;
123   if(fhPedestalsWid)
124     delete fhPedestalsWid ;
125   if(fctdb)
126     delete fctdb ;
127   if(fRunList)
128     delete fRunList ;
129 }
130 //____________________________________________________________________________ 
131 void AliPHOSCalibrator::AddRun(const char * filename)
132 {
133   //Adds one more run to list of runs, which will be scanned in ScanXXX methods
134   
135   TObjString * fn = new TObjString(filename) ;
136   if(!fRunList){
137     fRunList=new TList() ;
138     fRunList->SetOwner() ;
139     fRunList->Add(fn) ;
140     return ;
141   }
142   else{
143     TIter next(fRunList) ;
144     TObjString * r ;
145     while((r=(TObjString *)(next()))){
146       if(fn->String().CompareTo(r->String())==0){
147         Error("Run already in list: ",filename) ;
148         return ;
149       }
150     }
151     fRunList->Add(fn) ;
152   }
153   
154 }
155 //____________________________________________________________________________ 
156 void AliPHOSCalibrator::Exec(Option_t * option)
157 {
158   // reads parameters and does the calibration
159   ScanPedestals(option);
160   CalculatePedestals();             
161   WritePedestals();
162   ScanGains(option) ;
163   CalculateGains() ;
164   WriteGains() ;
165 }
166 //____________________________________________________________________________ 
167 void AliPHOSCalibrator::Init(void)
168 {
169   // intializes everything
170
171   //check if ConTableDB already read 
172   if(!fctdb){     
173     TFile * v = gROOT->GetFile(fConTableDBFile) ;
174     if(!v)
175       v = TFile::Open(fConTableDBFile) ;
176     if(!v){
177       Fatal("Can not open file with Connection Table DB:",fConTableDBFile) ;
178       return ;
179     }  
180     fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
181   }
182
183   fNch = fctdb->GetNchanels() ;
184   fhPedestals   = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
185   fhPedestalsWid= new TH1F("hPedestalsWid","Pedestals width",fNch,0.,fNch) ;
186   fhGains       = new TH1F("hGains","Gains ",fNch,0.,fNch) ; 
187   fhGainsWid    = new TH1F("hGainsWid","Gains width",fNch,0.,fNch) ; 
188 }
189 //____________________________________________________________________________ 
190 void AliPHOSCalibrator::SetConTableDB(const char * file,const char * name)
191 {
192   //Reads Connection Table database with name "name" from file "file" 
193
194   if(file==0 || name == 0){
195     Error("Please, specify file with database"," and its title") ;
196     return ;
197   }
198   if(fctdb && strcmp(fctdb->GetTitle(),name)==0) //already read
199     return ;
200
201   //else read new one
202   if(fctdb){
203     delete fctdb ;
204     fctdb = 0;
205   }
206
207   TFile * v = gROOT->GetFile(fConTableDBFile) ;
208   if(!v)
209     v = TFile::Open(fConTableDBFile) ;
210   if(!v){
211     Error("Can not open file with Connection Table DB:",fConTableDBFile) ;
212     return ;
213   }  
214   fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
215   
216 }
217 //____________________________________________________________________________ 
218 void AliPHOSCalibrator::PlotPedestal(Int_t chanel)
219 {
220   //Plot histogram for a given channel, filled in Scan method
221   if(fPedHistos && fPedHistos->GetEntriesFast()){
222     static_cast<TH1F*>(fPedHistos->At(chanel))->Draw() ;
223   }
224   else{
225     printf("Histograms not created yet! \n") ;
226   } 
227 }
228 //____________________________________________________________________________ 
229 void AliPHOSCalibrator::PlotPedestals(void)
230 {
231   // draws pedestals distribution
232   fhPedestals->Draw() ;
233 }
234 //____________________________________________________________________________ 
235 void AliPHOSCalibrator::PlotGain(Int_t chanel)
236 {
237   //Plot histogram for a given channel, filled in Scan method
238   if(fGainHistos && fGainHistos->GetEntriesFast()){
239     static_cast<TH1F*>(fGainHistos->At(chanel))->Draw() ;
240   }
241   else{
242     printf("Histograms not created yet! \n") ;
243   } 
244 }
245 //____________________________________________________________________________ 
246 void AliPHOSCalibrator::PlotGains(void)
247 {
248   // draws gains distribution
249   fhGains->Draw() ;
250 }
251 //____________________________________________________________________________ 
252 void AliPHOSCalibrator::ScanPedestals(Option_t * option )
253 {
254   //scan all files in list fRunList and fill pedestal hisgrams
255   //option: "clear" - clear pedestal histograms filled up to now
256   //        "deb" - plot file name currently processed
257
258   if(!fctdb)
259     Init() ;
260
261   if(fPedHistos && strstr(option,"clear"))
262     fPedHistos->Delete() ;
263   if(!fPedHistos)
264     fPedHistos = new TObjArray(fNch) ;
265
266   //Create histos for each channel, fills them and extracts mean values.
267   //First - prepare histos
268   Int_t ich ;
269   for(ich=0;ich<fNch ;ich++){
270     TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
271     if(!h ){
272       TString n("hPed");
273       n+=ich ;
274       TString name("Pedestal for channel ") ;
275       name += ich ;
276       fPedHistos->AddAt(new TH1F(n,name,fNChan,0,fNChan),ich) ; 
277     }
278   }
279
280   TIter next(fRunList) ;
281   TObjString * file ;
282   while((file = static_cast<TObjString *>(next()))){
283     if(strstr(option,"deb"))
284       printf("Processing file %s \n ",file->String().Data()) ;
285     AliPHOSGetter * gime = AliPHOSGetter::Instance(file->String().Data(),GetTitle()) ;
286     Int_t ievent ;
287     for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
288       gime->Event(ievent,"D") ;
289       if(gime->EventPattern() == fPedPat ){
290         Int_t idigit ;
291         TClonesArray * digits = gime->Digits() ;
292         for(idigit = 0; idigit<digits->GetEntriesFast() ;  idigit++){
293           AliPHOSDigit * digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
294           ich = fctdb->AbsId2Raw(digit->GetId());
295
296           if(ich>=0){
297             Float_t amp = digit->GetAmp() ;
298             TH1F * hh = dynamic_cast<TH1F*>(fPedHistos->At(ich)) ;
299             hh->Fill(amp) ;
300           }
301         }
302       }
303     }
304   }
305    
306 }
307 //____________________________________________________________________________ 
308 void AliPHOSCalibrator::CalculatePedestals()
309 {
310   //Fit histograms, filled in ScanPedestals method with Gaussian
311   //find mean and width, check deviation from mean for each channel.
312
313   if(!fPedHistos || !fPedHistos->At(0)){
314     Error("CalculatePedestals","You should run ScanPedestals first!") ;
315     return ;
316   }
317
318   //Now fit results with Gauss
319   TF1 * gs = new TF1("gs","gaus",0.,10000.) ;
320   Int_t ich ;
321   for(ich=0;ich<fNch ;ich++){
322     TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
323     Int_t max = h->GetMaximumBin() ;
324     Axis_t xmin = max/2. ;
325     Axis_t xmax = max*3/2 ;
326     gs->SetRange(xmin,xmax) ;
327     Double_t par[3] ;
328     par[0] = h->GetBinContent(max) ;
329     par[1] = max ;
330     par[2] = max/3 ;
331     gs->SetParameters(par[0],par[1],par[2]) ;
332     h->Fit("gs","QR") ;
333     gs->GetParameters(par) ;
334     fhPedestals->SetBinContent(ich,par[1]) ;
335     fhPedestals->SetBinError(ich,par[2]) ;
336     fhPedestalsWid->Fill(ich,par[2]) ;
337   }
338   delete gs ;
339
340   //now check reasonability of results
341   TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
342   fhPedestals->Fit("p0","Q") ;
343   Double_t meanPed ;
344   p0->GetParameters(&meanPed);
345   for(ich=0;ich<fNch ;ich++){
346     Float_t ped =  fhPedestals->GetBinContent(ich) ;
347     if(ped < 0  || ped > meanPed+fAcceptCorr){
348       TString out("Pedestal of channel ") ;
349       out+=ich ;     
350       out+=" is ";
351       out+= ped ;
352       out+= "it is too far from mean " ;
353       out+= meanPed ;
354       Error("PHOSCalibrator",out) ;
355     }
356   }
357   delete p0 ;
358     
359 }
360 //____________________________________________________________________________ 
361 void AliPHOSCalibrator::ScanGains(Option_t * option)
362 {
363   //Scan all runs, listed in fRunList and fill histograms for all channels
364   //options: "clear" - clean histograms, filled up to now
365   //         "deb" - print current file name
366   //         "narrow" - scan only narrow beam events
367
368   if(!fctdb)
369     Init() ;
370   if(fGainHistos && strstr(option,"clear"))
371     fGainHistos->Delete() ;
372   if(!fGainHistos){
373     if(strstr(option,"deball"))
374       printf("creating array for %d channels \n",fNch) ;            
375     fGainHistos   = new TObjArray(fNch) ;
376   }
377
378   //Create histos for each channel, fills them and extracts mean values.
379   //First - prepare  histos
380
381   if(!fGainHistos->GetEntriesFast()){
382     Int_t ich ;
383     for(ich=0;ich<fNch ;ich++){   
384       TString n("hGain");
385       n+=ich ;
386       TString name("Gains for channel ") ;
387       name += ich ;
388       fGainHistos->AddAt(new TH1F(n,name,fNGainBins,0,fGainMax),ich) ; 
389       //      static_cast<TH1F*>(fGainHistos->At(ich))->Sumw2() ;
390     }
391   }
392
393   Bool_t all =!(Bool_t)strstr(option,"narrow") ;
394
395
396   TIter next(fRunList) ;
397   TObjString * file ;
398   while((file = static_cast<TObjString *>(next()))){
399     
400     AliPHOSGetter * gime = AliPHOSGetter::Instance(file->String().Data(),GetTitle()) ;
401     Int_t handled = 0;
402     Int_t ievent ;
403     for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
404       gime->Event(ievent,"D") ;
405       if(gime->EventPattern() == fNBPat || 
406          (all && gime->EventPattern() == fWBPat)){
407         handled ++ ;
408         Int_t idigit ;
409         TClonesArray * digits = gime->Digits() ;
410         AliPHOSDigit * digit ;
411         Int_t max = 0 ;
412         Int_t imax = 0;
413         for(idigit = 0; idigit<digits->GetEntriesFast() ;  idigit++){
414           digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
415           if(digit->GetAmp() > max){
416             imax = idigit ;
417             max = digit->GetAmp() ;
418           }
419         }
420         digit = static_cast<AliPHOSDigit *>(digits->At(imax) ) ;
421         Int_t ich = fctdb->AbsId2Raw(digit->GetId());
422         if(ich>=0){
423           Float_t pedestal = fhPedestals->GetBinContent(ich) ;
424           const Float_t kshowerInCrystall = 0.9 ;
425           Float_t beamEnergy = gime->BeamEnergy() ;
426           Float_t gain = beamEnergy*kshowerInCrystall/
427                          (digit->GetAmp() - pedestal) ;
428           static_cast<TH1F*>(fGainHistos->At(ich))->Fill(gain) ;
429         } 
430       }
431     }
432     if(strstr(option,"deb"))
433       printf("Hadled %d events \n",handled) ;
434   }
435 }   
436 //____________________________________________________________________________ 
437 void AliPHOSCalibrator::CalculateGains(void)
438 {
439   //calculates gain
440
441   if(!fGainHistos || !fGainHistos->GetEntriesFast()){
442     Error("CalculateGains","You should run ScanGains first!") ; 
443     return ;
444   }
445
446   //Fit results with Landau
447   TF1 * gs = new TF1("gs","landau",0.,10000.) ;
448   Int_t ich ;
449   for(ich=0;ich<fNch ;ich++){
450     TH1F * h = static_cast<TH1F *>(fGainHistos->At(ich)) ;
451     Int_t bmax = h->GetMaximumBin() ;
452     Axis_t center = h->GetBinCenter(bmax) ;
453     Axis_t xmin = center - 0.01 ;
454     Axis_t xmax = center + 0.02 ;
455     gs->SetRange(xmin,xmax) ;
456     Double_t par[3] ;
457     par[0] = h->GetBinContent(bmax) ;
458     par[1] = center ;
459     par[2] = 0.001 ;
460     gs->SetParameters(par[0],par[1],par[2]) ;
461     h->Fit("gs","QR") ;
462     gs->GetParameters(par) ;
463     fhGains->SetBinContent(ich,par[1]) ;
464     fhGains->SetBinError(ich,par[2]) ;
465     fhGainsWid->Fill(ich,par[2]) ;
466   }
467   delete gs ;
468
469   //now check reasonability of results
470   TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
471   fhGains->Fit("p0","Q") ;
472   Double_t meanGain ;
473   p0->GetParameters(&meanGain);
474   for(ich=0;ich<fNch ;ich++){
475     Float_t gain =  fhGains->GetBinContent(ich) ;
476     if(gain < meanGain/fGainAcceptCorr  || gain > meanGain*fGainAcceptCorr){
477       TString out("Gain of channel ") ;
478       out+=ich ;     
479       out+=" is ";
480       out+= gain ;
481       out+= "it is too far from mean " ;
482       out+= meanGain ;
483       Error("PHOSCalibrator",out) ;
484     }
485   }
486   delete p0 ;
487     
488 }
489
490 //_____________________________________________________________________________
491 void AliPHOSCalibrator::WritePedestals(const char * version,
492                                          Int_t begin,Int_t end)
493 {
494   //Write calculated data to file using AliPHOSCalibrManager
495   //version and validitirange (begin-end) will be used to identify data 
496
497   if(!fctdb){
498     Error("WritePedestals","\n           Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
499     return ;
500   }
501   //fill data   
502   AliPHOSCalibrationData ped("Pedestals",version);
503   for(Int_t i=0; i<fNch;i++){
504     Int_t absid=fctdb->Raw2AbsId(i) ;
505     ped.SetData(absid,fhPedestals->GetBinContent(i)) ;
506     ped.SetDataCheck(absid,fhPedestalsWid->GetBinContent(i)) ;
507   }
508
509   //evaluate validity range
510   if(begin==0){
511     TIter next(fRunList) ;
512     Int_t ibegin=99999;
513     Int_t iend=0 ;
514     TObjString * file ;
515     while((file=((TObjString*)next()))){
516        TString s = file->GetString() ;
517        TString ss = s(s.Last('_'),s.Last('.'));
518        Int_t tmp ;
519        if(sscanf(ss.Data(),"%d",&tmp)){
520          if(ibegin<tmp)
521            ibegin=tmp ;  
522          if(iend>tmp)
523            iend=tmp ;
524        }
525     } 
526     ped.SetValidityRange(ibegin,iend) ;
527   }
528   else    
529     ped.SetValidityRange(begin,end) ;
530
531   //check, may be Manager instance already configured?
532   AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
533   if(!cmngr){
534     Warning("Write Pedestals","Using database file 'PHOSBTCalibration.root'") ;
535     cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
536   }
537   cmngr->WriteData(&ped) ;
538 }       
539 //_____________________________________________________________________________
540 void AliPHOSCalibrator::ReadPedestals(const char * version,
541                                          Int_t range)
542
543   //Read data from file using AliPHOSCalibrManager 
544   //version and range will be used to choose proper data
545
546   AliPHOSCalibrationData ped("Pedestals",version);
547   AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
548   if(!cmngr){
549    Warning("ReadPedestals","Using database file 'PHOSBTCalibration.root'") ;
550    cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
551   }
552   cmngr->ReadFromRoot(ped,range) ;
553   Int_t npeds=ped.NChannels() ;
554   fNch = fctdb->GetNchanels() ;
555   if(fhPedestals)
556     delete fhPedestals ;
557   fhPedestals   = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
558   for(Int_t i=0;i<npeds;i++){
559     Int_t raw =fctdb->AbsId2Raw(i) ;
560     if(raw){
561       fhPedestals->SetBinContent(raw-1,ped.Data(i)) ;
562       fhPedestals->SetBinError(raw-1,ped.DataCheck(i)) ;
563     }
564   }
565 }       
566 //_____________________________________________________________________________
567 void AliPHOSCalibrator::ReadGains(const char * version,
568                                          Int_t range)
569
570   //Read data from file using AliPHOSCalibrManager 
571   //version and range will be used to choose proper data
572
573   AliPHOSCalibrationData gains("Gains",version);
574   AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
575   if(!cmngr){
576     Warning("ReadGainss","Using database file 'PHOSBTCalibration.root'") ;
577     cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
578   }
579   cmngr->ReadFromRoot(gains,range) ;
580   Int_t npeds=gains.NChannels() ;
581   fNch = fctdb->GetNchanels() ;
582   if(fhGains)
583     delete fhGains ;
584   fhGains   = new TH1F("hGainss","Gains mean",fNch,0.,fNch) ;
585   for(Int_t i=0;i<npeds;i++){
586     Int_t raw =fctdb->AbsId2Raw(i) ;
587     if(raw){
588       fhGains->SetBinContent(raw-1,gains.Data(i)) ;
589       fhGains->SetBinError(raw-1,gains.DataCheck(i)) ;
590     }
591   }
592 }       
593 //_____________________________________________________________________________
594 void AliPHOSCalibrator::WriteGains(const char * version,
595                                      Int_t begin,Int_t end)
596
597   //Write gains through AliPHOSCalibrManager
598   //version and validity range(begin-end) are used to identify data
599
600   if(!fctdb){
601     Error("WriteGains","\n        Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
602     return ;
603   }
604
605   AliPHOSCalibrationData gains("Gains",version);
606   for(Int_t i=0; i<fNch;i++){
607     Int_t absid=fctdb->Raw2AbsId(i) ;
608     gains.SetData(absid,fhGains->GetBinContent(i)) ;
609     gains.SetDataCheck(absid,fhGainsWid->GetBinContent(i)) ;
610   }
611   if(begin==0){
612     TIter next(fRunList) ;
613     Int_t ibegin=99999;
614     Int_t iend=0 ;
615     TObjString * file ;
616     while((file=((TObjString*)next()))){
617        TString s = file->GetString() ;
618        TSubString ss = s(s.Last('_'),s.Last('.'));
619        Int_t tmp ;
620        if(sscanf(ss.Data(),"%d",&tmp)){
621          if(ibegin<tmp)
622            ibegin=tmp ;  
623          if(iend>tmp)
624            iend=tmp ;
625        }
626     } 
627     gains.SetValidityRange(ibegin,iend) ;
628   }
629   else    
630     gains.SetValidityRange(begin,end) ;
631   AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
632   if(!cmngr){
633     Warning("WriteGains","Using database file 'PHOSBTCalibration.root'") ;
634     cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
635   }
636   cmngr->WriteData(&gains) ;
637 }       
638 //_____________________________________________________________________________
639 void AliPHOSCalibrator::Print()const 
640 {
641   // prints everything
642   printf("--------------AliPHOSCalibrator-----------------\n") ;
643   printf("Files to handle:\n") ;
644   TIter next(fRunList) ;
645   TObjString * r ;
646   while((r=(TObjString *)(next())))
647       printf("                %s\n",r->GetName()) ;
648
649   printf("Name of ConTableDB:.....................%s\n",fConTableDB.Data()) ;
650   printf("File of ConTableDB:.....................%s\n",fConTableDBFile.Data() ) ;
651   printf("Maximal deviation from mean Gain (factor):.%f\n",fGainAcceptCorr) ;
652   printf("Maximal deviation of Pedestal from mean:...%f\n",fAcceptCorr) ; 
653   printf("Range used in Gain histos:..............%f\n",fGainMax) ;
654   printf("Number of bins in Gain histos:..........%d\n",fNGainBins) ;
655   printf("Number of channels to calibrate:........%d\n",fNch) ;
656   printf("Number of bins in pedestal histos:......%d\n",fNChan) ;
657   printf("trigger pattern for PEDESTAL events:....%d\n",fPedPat) ;
658   printf("trigger pattern for PULSER events:......%d\n",fPulPat) ;
659   printf("trigger pattern for LED events:.........%d\n",fLEDPat) ;
660   printf("trigger pattern for WIDE BEAM events:...%d\n",fWBPat) ;
661   printf("trigger pattern for NARROW BEAM events:.%d\n",fNBPat) ;
662   printf("--------------------------------------------------\n") ;
663 }