1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 //_________________________________________________________________________
18 // Class to calculate calibration parameters from beam tests etc.
19 // First pass - one should calcuate pedestals, in the second pass - gains.
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'.
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.
32 // Finally fills database.
35 // AliPHOSCalibrator * c = new AliPHOSCalibrator("path/galice.root") ;
36 // c->AddRun("path2/galice.root") ;
37 // c->ScanPedestals();
38 // c->CalculatePedestals();
39 // c->WritePedestals();
41 // c->CalculateGains() ;
44 //*-- Author : D.Peressounko (RRC KI)
45 //////////////////////////////////////////////////////////////////////////////
47 // --- ROOT system ---
50 #include "TObjString.h"
51 // --- Standard library ---
53 // --- AliRoot header files ---
54 #include "AliPHOSGetter.h"
55 #include "AliPHOSCalibrator.h"
56 #include "AliPHOSConTableDB.h"
57 #include "AliPHOSCalibrManager.h"
58 #include "AliPHOSCalibrationData.h"
60 ClassImp(AliPHOSCalibrator)
63 //____________________________________________________________________________
64 AliPHOSCalibrator::AliPHOSCalibrator():TTask("AliPHOSCalibrator","Default")
66 //Default constuctor for root. Normally should not be used
74 fConTableDB = "Beamtest2002" ;
75 fConTableDBFile = "ConTableDB.root" ;
77 //____________________________________________________________________________
78 AliPHOSCalibrator::AliPHOSCalibrator(const char* file, const char* title):
79 TTask("AliPHOSCalibrator",title)
81 //Constructor which should normally be used.
82 //file: path/galice.root - header file
83 //title: branch name of PHOS reconstruction (e.g. "Default")
86 fRunList = new TList() ;
87 fRunList->SetOwner() ;
88 fRunList->Add(new TObjString(file)) ;
90 fPedPat = 257 ; //Patterns for different kind of events
99 fAcceptCorr = 10 ; //Maximal deviation from mean, considered as normal
101 fGainAcceptCorr = 5 ; //Factor for gain deviation
107 fConTableDB = "Beamtest2002" ;
108 fConTableDBFile = "ConTableDB.root" ;
111 //____________________________________________________________________________
112 AliPHOSCalibrator::~AliPHOSCalibrator()
122 delete fhPedestalsWid ;
128 //____________________________________________________________________________
129 void AliPHOSCalibrator::AddRun(const char * filename)
131 //Adds one more run to list of runs, which will be scanned in ScanXXX methods
133 TObjString * fn = new TObjString(filename) ;
135 fRunList=new TList() ;
136 fRunList->SetOwner() ;
141 TIter next(fRunList) ;
143 while((r=(TObjString *)(next()))){
144 if(fn->String().CompareTo(r->String())==0){
145 Error("Run already in list: ",filename) ;
153 //____________________________________________________________________________
154 void AliPHOSCalibrator::Exec(Option_t * option)
156 // reads parameters and does the calibration
157 ScanPedestals(option);
158 CalculatePedestals();
164 //____________________________________________________________________________
165 void AliPHOSCalibrator::Init(void)
167 // intializes everything
169 //check if ConTableDB already read
171 TFile * v = gROOT->GetFile(fConTableDBFile) ;
173 v = TFile::Open(fConTableDBFile) ;
175 Fatal("Can not open file with Connection Table DB:",fConTableDBFile) ;
178 fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
181 fNch = fctdb->GetNchanels() ;
182 fhPedestals = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
183 fhPedestalsWid= new TH1F("hPedestalsWid","Pedestals width",fNch,0.,fNch) ;
184 fhGains = new TH1F("hGains","Gains ",fNch,0.,fNch) ;
185 fhGainsWid = new TH1F("hGainsWid","Gains width",fNch,0.,fNch) ;
187 //____________________________________________________________________________
188 void AliPHOSCalibrator::SetConTableDB(const char * file,const char * name)
190 //Reads Connection Table database with name "name" from file "file"
192 if(file==0 || name == 0){
193 Error("Please, specify file with database"," and its title") ;
196 if(fctdb && strcmp(fctdb->GetTitle(),name)==0) //already read
205 TFile * v = gROOT->GetFile(fConTableDBFile) ;
207 v = TFile::Open(fConTableDBFile) ;
209 Error("Can not open file with Connection Table DB:",fConTableDBFile) ;
212 fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
215 //____________________________________________________________________________
216 void AliPHOSCalibrator::PlotPedestal(Int_t chanel)
218 //Plot histogram for a given channel, filled in Scan method
219 if(fPedHistos && fPedHistos->GetEntriesFast()){
220 static_cast<TH1F*>(fPedHistos->At(chanel))->Draw() ;
223 printf("Histograms not created yet! \n") ;
226 //____________________________________________________________________________
227 void AliPHOSCalibrator::PlotPedestals(void)
229 // draws pedestals distribution
230 fhPedestals->Draw() ;
232 //____________________________________________________________________________
233 void AliPHOSCalibrator::PlotGain(Int_t chanel)
235 //Plot histogram for a given channel, filled in Scan method
236 if(fGainHistos && fGainHistos->GetEntriesFast()){
237 static_cast<TH1F*>(fGainHistos->At(chanel))->Draw() ;
240 printf("Histograms not created yet! \n") ;
243 //____________________________________________________________________________
244 void AliPHOSCalibrator::PlotGains(void)
246 // draws gains distribution
249 //____________________________________________________________________________
250 void AliPHOSCalibrator::ScanPedestals(Option_t * option )
252 //scan all files in list fRunList and fill pedestal hisgrams
253 //option: "clear" - clear pedestal histograms filled up to now
254 // "deb" - plot file name currently processed
259 if(fPedHistos && strstr(option,"clear"))
260 fPedHistos->Delete() ;
262 fPedHistos = new TObjArray(fNch) ;
264 //Create histos for each channel, fills them and extracts mean values.
265 //First - prepare histos
267 for(ich=0;ich<fNch ;ich++){
268 TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
272 TString name("Pedestal for channel ") ;
274 fPedHistos->AddAt(new TH1F(n,name,fNChan,0,fNChan),ich) ;
278 TIter next(fRunList) ;
280 while((file = static_cast<TObjString *>(next()))){
281 if(strstr(option,"deb"))
282 printf("Processing file %s \n ",file->String().Data()) ;
283 AliPHOSGetter * gime = AliPHOSGetter::Instance(file->String().Data(),GetTitle()) ;
285 for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
286 gime->Event(ievent,"D") ;
287 if(gime->EventPattern() == fPedPat ){
289 TClonesArray * digits = gime->Digits() ;
290 for(idigit = 0; idigit<digits->GetEntriesFast() ; idigit++){
291 AliPHOSDigit * digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
292 ich = fctdb->AbsId2Raw(digit->GetId());
295 Float_t amp = digit->GetAmp() ;
296 TH1F * hh = dynamic_cast<TH1F*>(fPedHistos->At(ich)) ;
305 //____________________________________________________________________________
306 void AliPHOSCalibrator::CalculatePedestals()
308 //Fit histograms, filled in ScanPedestals method with Gaussian
309 //find mean and width, check deviation from mean for each channel.
311 if(!fPedHistos || !fPedHistos->At(0)){
312 Error("CalculatePedestals","You should run ScanPedestals first!") ;
316 //Now fit results with Gauss
317 TF1 * gs = new TF1("gs","gaus",0.,10000.) ;
319 for(ich=0;ich<fNch ;ich++){
320 TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
321 Int_t max = h->GetMaximumBin() ;
322 Axis_t xmin = max/2. ;
323 Axis_t xmax = max*3/2 ;
324 gs->SetRange(xmin,xmax) ;
326 par[0] = h->GetBinContent(max) ;
329 gs->SetParameters(par[0],par[1],par[2]) ;
331 gs->GetParameters(par) ;
332 fhPedestals->SetBinContent(ich,par[1]) ;
333 fhPedestals->SetBinError(ich,par[2]) ;
334 fhPedestalsWid->Fill(ich,par[2]) ;
338 //now check reasonability of results
339 TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
340 fhPedestals->Fit("p0","Q") ;
342 p0->GetParameters(&meanPed);
343 for(ich=0;ich<fNch ;ich++){
344 Float_t ped = fhPedestals->GetBinContent(ich) ;
345 if(ped < 0 || ped > meanPed+fAcceptCorr){
346 TString out("Pedestal of channel ") ;
350 out+= "it is too far from mean " ;
352 Error("PHOSCalibrator",out) ;
358 //____________________________________________________________________________
359 void AliPHOSCalibrator::ScanGains(Option_t * option)
361 //Scan all runs, listed in fRunList and fill histograms for all channels
362 //options: "clear" - clean histograms, filled up to now
363 // "deb" - print current file name
364 // "narrow" - scan only narrow beam events
368 if(fGainHistos && strstr(option,"clear"))
369 fGainHistos->Delete() ;
371 if(strstr(option,"deball"))
372 printf("creating array for %d channels \n",fNch) ;
373 fGainHistos = new TObjArray(fNch) ;
376 //Create histos for each channel, fills them and extracts mean values.
377 //First - prepare histos
379 if(!fGainHistos->GetEntriesFast()){
381 for(ich=0;ich<fNch ;ich++){
384 TString name("Gains for channel ") ;
386 fGainHistos->AddAt(new TH1F(n,name,fNGainBins,0,fGainMax),ich) ;
387 // static_cast<TH1F*>(fGainHistos->At(ich))->Sumw2() ;
391 Bool_t all =!(Bool_t)strstr(option,"narrow") ;
394 TIter next(fRunList) ;
396 while((file = static_cast<TObjString *>(next()))){
398 AliPHOSGetter * gime = AliPHOSGetter::Instance(file->String().Data(),GetTitle()) ;
401 for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
402 gime->Event(ievent,"D") ;
403 if(gime->EventPattern() == fNBPat ||
404 (all && gime->EventPattern() == fWBPat)){
407 TClonesArray * digits = gime->Digits() ;
408 AliPHOSDigit * digit ;
411 for(idigit = 0; idigit<digits->GetEntriesFast() ; idigit++){
412 digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
413 if(digit->GetAmp() > max){
415 max = digit->GetAmp() ;
418 digit = static_cast<AliPHOSDigit *>(digits->At(imax) ) ;
419 Int_t ich = fctdb->AbsId2Raw(digit->GetId());
421 Float_t pedestal = fhPedestals->GetBinContent(ich) ;
422 const Float_t kshowerInCrystall = 0.9 ;
423 Float_t beamEnergy = gime->BeamEnergy() ;
424 Float_t gain = beamEnergy*kshowerInCrystall/
425 (digit->GetAmp() - pedestal) ;
426 static_cast<TH1F*>(fGainHistos->At(ich))->Fill(gain) ;
430 if(strstr(option,"deb"))
431 printf("Hadled %d events \n",handled) ;
434 //____________________________________________________________________________
435 void AliPHOSCalibrator::CalculateGains(void)
439 if(!fGainHistos || !fGainHistos->GetEntriesFast()){
440 Error("CalculateGains","You should run ScanGains first!") ;
444 //Fit results with Landau
445 TF1 * gs = new TF1("gs","landau",0.,10000.) ;
447 for(ich=0;ich<fNch ;ich++){
448 TH1F * h = static_cast<TH1F *>(fGainHistos->At(ich)) ;
449 Int_t bmax = h->GetMaximumBin() ;
450 Axis_t center = h->GetBinCenter(bmax) ;
451 Axis_t xmin = center - 0.01 ;
452 Axis_t xmax = center + 0.02 ;
453 gs->SetRange(xmin,xmax) ;
455 par[0] = h->GetBinContent(bmax) ;
458 gs->SetParameters(par[0],par[1],par[2]) ;
460 gs->GetParameters(par) ;
461 fhGains->SetBinContent(ich,par[1]) ;
462 fhGains->SetBinError(ich,par[2]) ;
463 fhGainsWid->Fill(ich,par[2]) ;
467 //now check reasonability of results
468 TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
469 fhGains->Fit("p0","Q") ;
471 p0->GetParameters(&meanGain);
472 for(ich=0;ich<fNch ;ich++){
473 Float_t gain = fhGains->GetBinContent(ich) ;
474 if(gain < meanGain/fGainAcceptCorr || gain > meanGain*fGainAcceptCorr){
475 TString out("Gain of channel ") ;
479 out+= "it is too far from mean " ;
481 Error("PHOSCalibrator",out) ;
488 //_____________________________________________________________________________
489 void AliPHOSCalibrator::WritePedestals(const char * version,
490 Int_t begin,Int_t end)
492 //Write calculated data to file using AliPHOSCalibrManager
493 //version and validitirange (begin-end) will be used to identify data
496 Error("WritePedestals","\n Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
500 AliPHOSCalibrationData ped("Pedestals",version);
501 for(Int_t i=0; i<fNch;i++){
502 Int_t absid=fctdb->Raw2AbsId(i) ;
503 ped.SetData(absid,fhPedestals->GetBinContent(i)) ;
504 ped.SetDataCheck(absid,fhPedestalsWid->GetBinContent(i)) ;
507 //evaluate validity range
509 TIter next(fRunList) ;
513 while((file=((TObjString*)next()))){
514 TString s = file->GetString() ;
515 TString ss = s(s.Last('_'),s.Last('.'));
517 if(sscanf(ss.Data(),"%d",&tmp)){
524 ped.SetValidityRange(ibegin,iend) ;
527 ped.SetValidityRange(begin,end) ;
529 //check, may be Manager instance already configured?
530 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
532 Warning("Write Pedestals","Using database file 'PHOSBTCalibration.root'") ;
533 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
535 cmngr->WriteData(&ped) ;
537 //_____________________________________________________________________________
538 void AliPHOSCalibrator::ReadPedestals(const char * version,
541 //Read data from file using AliPHOSCalibrManager
542 //version and range will be used to choose proper data
544 AliPHOSCalibrationData ped("Pedestals",version);
545 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
547 Warning("ReadPedestals","Using database file 'PHOSBTCalibration.root'") ;
548 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
550 cmngr->ReadFromRoot(ped,range) ;
551 Int_t npeds=ped.NChannels() ;
552 fNch = fctdb->GetNchanels() ;
555 fhPedestals = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
556 for(Int_t i=0;i<npeds;i++){
557 Int_t raw =fctdb->AbsId2Raw(i) ;
559 fhPedestals->SetBinContent(raw-1,ped.Data(i)) ;
560 fhPedestals->SetBinError(raw-1,ped.DataCheck(i)) ;
564 //_____________________________________________________________________________
565 void AliPHOSCalibrator::ReadGains(const char * version,
568 //Read data from file using AliPHOSCalibrManager
569 //version and range will be used to choose proper data
571 AliPHOSCalibrationData gains("Gains",version);
572 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
574 Warning("ReadGainss","Using database file 'PHOSBTCalibration.root'") ;
575 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
577 cmngr->ReadFromRoot(gains,range) ;
578 Int_t npeds=gains.NChannels() ;
579 fNch = fctdb->GetNchanels() ;
582 fhGains = new TH1F("hGainss","Gains mean",fNch,0.,fNch) ;
583 for(Int_t i=0;i<npeds;i++){
584 Int_t raw =fctdb->AbsId2Raw(i) ;
586 fhGains->SetBinContent(raw-1,gains.Data(i)) ;
587 fhGains->SetBinError(raw-1,gains.DataCheck(i)) ;
591 //_____________________________________________________________________________
592 void AliPHOSCalibrator::WriteGains(const char * version,
593 Int_t begin,Int_t end)
595 //Write gains through AliPHOSCalibrManager
596 //version and validity range(begin-end) are used to identify data
599 Error("WriteGains","\n Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
603 AliPHOSCalibrationData gains("Gains",version);
604 for(Int_t i=0; i<fNch;i++){
605 Int_t absid=fctdb->Raw2AbsId(i) ;
606 gains.SetData(absid,fhGains->GetBinContent(i)) ;
607 gains.SetDataCheck(absid,fhGainsWid->GetBinContent(i)) ;
610 TIter next(fRunList) ;
614 while((file=((TObjString*)next()))){
615 TString s = file->GetString() ;
616 TSubString ss = s(s.Last('_'),s.Last('.'));
618 if(sscanf(ss.Data(),"%d",&tmp)){
625 gains.SetValidityRange(ibegin,iend) ;
628 gains.SetValidityRange(begin,end) ;
629 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
631 Warning("WriteGains","Using database file 'PHOSBTCalibration.root'") ;
632 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
634 cmngr->WriteData(&gains) ;
636 //_____________________________________________________________________________
637 void AliPHOSCalibrator::Print(const Option_t * option)const
640 printf("--------------AliPHOSCalibrator-----------------\n") ;
641 printf("Files to handle:\n") ;
642 TIter next(fRunList) ;
644 while((r=(TObjString *)(next())))
645 printf(" %s\n",r->GetName()) ;
647 printf("Name of ConTableDB:.....................%s\n",fConTableDB.Data()) ;
648 printf("File of ConTableDB:.....................%s\n",fConTableDBFile.Data() ) ;
649 printf("Maximal deviation from mean Gain (factor):.%f\n",fGainAcceptCorr) ;
650 printf("Maximal deviation of Pedestal from mean:...%f\n",fAcceptCorr) ;
651 printf("Range used in Gain histos:..............%f\n",fGainMax) ;
652 printf("Number of bins in Gain histos:..........%d\n",fNGainBins) ;
653 printf("Number of channels to calibrate:........%d\n",fNch) ;
654 printf("Number of bins in pedestal histos:......%d\n",fNChan) ;
655 printf("trigger pattern for PEDESTAL events:....%d\n",fPedPat) ;
656 printf("trigger pattern for PULSER events:......%d\n",fPulPat) ;
657 printf("trigger pattern for LED events:.........%d\n",fLEDPat) ;
658 printf("trigger pattern for WIDE BEAM events:...%d\n",fWBPat) ;
659 printf("trigger pattern for NARROW BEAM events:.%d\n",fNBPat) ;
660 printf("--------------------------------------------------\n") ;