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