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