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