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