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