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