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