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