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 ---
52 #include "TObjString.h"
53 // --- Standard library ---
55 // --- AliRoot header files ---
56 #include "AliPHOSGetter.h"
57 #include "AliPHOSCalibrator.h"
58 #include "AliPHOSConTableDB.h"
59 #include "AliPHOSCalibrManager.h"
60 #include "AliPHOSCalibrationData.h"
62 ClassImp(AliPHOSCalibrator)
65 //____________________________________________________________________________
66 AliPHOSCalibrator::AliPHOSCalibrator():TTask("AliPHOSCalibrator","Default")
68 //Default constuctor for root. Normally should not be used
76 fConTableDB = "Beamtest2002" ;
77 fConTableDBFile = "ConTableDB.root" ;
79 //____________________________________________________________________________
80 AliPHOSCalibrator::AliPHOSCalibrator(const char* file, const char* title,Bool_t toSplit):
81 TTask("AliPHOSCalibrator",title)
82 { //Constructor which should normally be used.
83 //file: path/galice.root - header file
84 //title: branch name of PHOS reconstruction (e.g. "Default")
85 //toSplit: wether we work in Split mode?
87 fRunList = new TList() ;
88 fRunList->SetOwner() ;
89 fRunList->Add(new TObjString(file)) ;
92 fPedPat = 257 ; //Patterns for different kind of events
101 fAcceptCorr = 10 ; //Maximal deviation from mean, considered as normal
103 fGainAcceptCorr = 5 ; //Factor for gain deviation
109 fConTableDB = "Beamtest2002" ;
110 fConTableDBFile = "ConTableDB.root" ;
113 //____________________________________________________________________________
114 AliPHOSCalibrator::~AliPHOSCalibrator()
123 delete fhPedestalsWid ;
129 //____________________________________________________________________________
130 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 ScanPedestals(option);
157 CalculatePedestals();
163 //____________________________________________________________________________
164 void AliPHOSCalibrator::Init(void){
166 //check if ConTableDB already read
168 TFile * v = gROOT->GetFile(fConTableDBFile) ;
170 v = TFile::Open(fConTableDBFile) ;
172 Fatal("Can not open file with Connection Table DB:",fConTableDBFile) ;
175 fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
178 fNch = fctdb->GetNchanels() ;
179 fhPedestals = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
180 fhPedestalsWid= new TH1F("hPedestalsWid","Pedestals width",fNch,0.,fNch) ;
181 fhGains = new TH1F("hGains","Gains ",fNch,0.,fNch) ;
182 fhGainsWid = new TH1F("hGainsWid","Gains width",fNch,0.,fNch) ;
184 //____________________________________________________________________________
185 void AliPHOSCalibrator::SetConTableDB(const char * file,const char * name){
186 //Reads Connection Table database with name "name" from file "file"
188 if(file==0 || name == 0){
189 Error("Please, specify file with database"," and its title") ;
192 if(fctdb && strcmp(fctdb->GetTitle(),name)==0) //already read
201 TFile * v = gROOT->GetFile(fConTableDBFile) ;
203 v = TFile::Open(fConTableDBFile) ;
205 Error("Can not open file with Connection Table DB:",fConTableDBFile) ;
208 fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
211 //____________________________________________________________________________
212 void AliPHOSCalibrator::PlotPedestal(Int_t chanel){
213 //Plot histogram for a given channel, filled in Scan method
214 if(fPedHistos && fPedHistos->GetEntriesFast()){
215 static_cast<TH1F*>(fPedHistos->At(chanel))->Draw() ;
218 printf("Histograms not created yet! \n") ;
221 //____________________________________________________________________________
222 void AliPHOSCalibrator::PlotPedestals(void){
223 fhPedestals->Draw() ;
225 //____________________________________________________________________________
226 void AliPHOSCalibrator::PlotGain(Int_t chanel){
227 //Plot histogram for a given channel, filled in Scan method
228 if(fGainHistos && fGainHistos->GetEntriesFast()){
229 static_cast<TH1F*>(fGainHistos->At(chanel))->Draw() ;
232 printf("Histograms not created yet! \n") ;
235 //____________________________________________________________________________
236 void AliPHOSCalibrator::PlotGains(void){
239 //____________________________________________________________________________
240 void AliPHOSCalibrator::ScanPedestals(Option_t * option ){
241 //scan all files in list fRunList and fill pedestal hisgrams
242 //option: "clear" - clear pedestal histograms filled up to now
243 // "deb" - plot file name currently processed
248 if(fPedHistos && strstr(option,"clear"))
249 fPedHistos->Delete() ;
251 fPedHistos = new TObjArray(fNch) ;
253 //Create histos for each channel, fills them and extracts mean values.
254 //First - prepare histos
256 for(ich=0;ich<fNch ;ich++){
257 TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
261 TString name("Pedestal for channel ") ;
263 fPedHistos->AddAt(new TH1F(n,name,fNChan,0,fNChan),ich) ;
267 TIter next(fRunList) ;
269 while((file = static_cast<TObjString *>(next()))){
270 if(strstr(option,"deb"))
271 printf("Processing file %s \n ",file->String().Data()) ;
272 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(file->String().Data(),GetTitle(),fToSplit) ;
274 for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
275 gime->Event(ievent,"D") ;
276 if(gime->EventPattern() == fPedPat ){
278 TClonesArray * digits = gime->Digits() ;
279 for(idigit = 0; idigit<digits->GetEntriesFast() ; idigit++){
280 AliPHOSDigit * digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
281 ich = fctdb->AbsId2Raw(digit->GetId());
284 Float_t amp = digit->GetAmp() ;
285 TH1F * hh = dynamic_cast<TH1F*>(fPedHistos->At(ich)) ;
294 //____________________________________________________________________________
295 void AliPHOSCalibrator::CalculatePedestals(){
296 //Fit histograms, filled in ScanPedestals method with Gaussian
297 //find mean and width, check deviation from mean for each channel.
299 if(!fPedHistos || !fPedHistos->At(0)){
300 Error("CalculatePedestals","You should run ScanPedestals first!") ;
304 //Now fit results with Gauss
305 TF1 * gs = new TF1("gs","gaus",0.,10000.) ;
307 for(ich=0;ich<fNch ;ich++){
308 TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
309 Int_t max = h->GetMaximumBin() ;
310 Axis_t xmin = max/2. ;
311 Axis_t xmax = max*3/2 ;
312 gs->SetRange(xmin,xmax) ;
314 par[0] = h->GetBinContent(max) ;
317 gs->SetParameters(par[0],par[1],par[2]) ;
319 gs->GetParameters(par) ;
320 fhPedestals->SetBinContent(ich,par[1]) ;
321 fhPedestals->SetBinError(ich,par[2]) ;
322 fhPedestalsWid->Fill(ich,par[2]) ;
326 //now check reasonability of results
327 TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
328 fhPedestals->Fit("p0","Q") ;
330 p0->GetParameters(&meanPed);
331 for(ich=0;ich<fNch ;ich++){
332 Float_t ped = fhPedestals->GetBinContent(ich) ;
333 if(ped < 0 || ped > meanPed+fAcceptCorr){
334 TString out("Pedestal of channel ") ;
338 out+= "it is too far from mean " ;
340 Error("PHOSCalibrator",out) ;
346 //____________________________________________________________________________
347 void AliPHOSCalibrator::ScanGains(Option_t * option){
348 //Scan all runs, listed in fRunList and fill histograms for all channels
349 //options: "clear" - clean histograms, filled up to now
350 // "deb" - print current file name
351 // "narrow" - scan only narrow beam events
355 if(fGainHistos && strstr(option,"clear"))
356 fGainHistos->Delete() ;
358 if(strstr(option,"deball"))
359 printf("creating array for %d channels \n",fNch) ;
360 fGainHistos = new TObjArray(fNch) ;
363 //Create histos for each channel, fills them and extracts mean values.
364 //First - prepare histos
366 if(!fGainHistos->GetEntriesFast()){
368 for(ich=0;ich<fNch ;ich++){
371 TString name("Gains for channel ") ;
373 fGainHistos->AddAt(new TH1F(n,name,fNGainBins,0,fGainMax),ich) ;
374 // static_cast<TH1F*>(fGainHistos->At(ich))->Sumw2() ;
378 Bool_t all =!(Bool_t)strstr(option,"narrow") ;
381 TIter next(fRunList) ;
383 while((file = static_cast<TObjString *>(next()))){
385 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(file->String().Data(),GetTitle(),fToSplit) ;
388 for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
389 gime->Event(ievent,"D") ;
390 if(gime->EventPattern() == fNBPat ||
391 (all && gime->EventPattern() == fWBPat)){
394 TClonesArray * digits = gime->Digits() ;
395 AliPHOSDigit * digit ;
398 for(idigit = 0; idigit<digits->GetEntriesFast() ; idigit++){
399 digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
400 if(digit->GetAmp() > max){
402 max = digit->GetAmp() ;
405 digit = static_cast<AliPHOSDigit *>(digits->At(imax) ) ;
406 Int_t ich = fctdb->AbsId2Raw(digit->GetId());
408 Float_t pedestal = fhPedestals->GetBinContent(ich) ;
409 const Float_t showerInCrystall = 0.9 ;
410 Float_t beamEnergy = gime->BeamEnergy() ;
411 Float_t gain = beamEnergy*showerInCrystall/
412 (digit->GetAmp() - pedestal) ;
413 static_cast<TH1F*>(fGainHistos->At(ich))->Fill(gain) ;
417 if(strstr(option,"deb"))
418 printf("Hadled %d events \n",handled) ;
421 //____________________________________________________________________________
422 void AliPHOSCalibrator::CalculateGains(void){
424 if(!fGainHistos || !fGainHistos->GetEntriesFast()){
425 Error("CalculateGains","You should run ScanGains first!") ;
429 //Fit results with Landau
430 TF1 * gs = new TF1("gs","landau",0.,10000.) ;
432 for(ich=0;ich<fNch ;ich++){
433 TH1F * h = static_cast<TH1F *>(fGainHistos->At(ich)) ;
434 Int_t bmax = h->GetMaximumBin() ;
435 Axis_t center = h->GetBinCenter(bmax) ;
436 Axis_t xmin = center - 0.01 ;
437 Axis_t xmax = center + 0.02 ;
438 gs->SetRange(xmin,xmax) ;
440 par[0] = h->GetBinContent(bmax) ;
443 gs->SetParameters(par[0],par[1],par[2]) ;
445 gs->GetParameters(par) ;
446 fhGains->SetBinContent(ich,par[1]) ;
447 fhGains->SetBinError(ich,par[2]) ;
448 fhGainsWid->Fill(ich,par[2]) ;
452 //now check reasonability of results
453 TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
454 fhGains->Fit("p0","Q") ;
456 p0->GetParameters(&meanGain);
457 for(ich=0;ich<fNch ;ich++){
458 Float_t gain = fhGains->GetBinContent(ich) ;
459 if(gain < meanGain/fGainAcceptCorr || gain > meanGain*fGainAcceptCorr){
460 TString out("Gain of channel ") ;
464 out+= "it is too far from mean " ;
466 Error("PHOSCalibrator",out) ;
473 //_____________________________________________________________________________
474 void AliPHOSCalibrator::WritePedestals(const char * version,
475 Int_t begin,Int_t end){
476 //Write calculated data to file using AliPHOSCalibrManager
477 //version and validitirange (begin-end) will be used to identify data
480 Error("WritePedestals","\n Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
484 AliPHOSCalibrationData ped("Pedestals",version);
485 for(Int_t i=0; i<fNch;i++){
486 Int_t absid=fctdb->Raw2AbsId(i) ;
487 ped.SetData(absid,fhPedestals->GetBinContent(i)) ;
488 ped.SetDataCheck(absid,fhPedestalsWid->GetBinContent(i)) ;
491 //evaluate validity range
493 TIter next(fRunList) ;
497 while((file=((TObjString*)next()))){
498 TString s = file->GetString() ;
499 TString ss = s(s.Last('_'),s.Last('.'));
501 if(sscanf(ss.Data(),"%d",&tmp)){
508 ped.SetValidityRange(ibegin,iend) ;
511 ped.SetValidityRange(begin,end) ;
513 //check, may be Manager instance already configured?
514 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
516 Warning("Write Pedestals","Using database file 'PHOSBTCalibration.root'") ;
517 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
519 cmngr->WriteData(&ped) ;
521 //_____________________________________________________________________________
522 void AliPHOSCalibrator::ReadPedestals(const char * version,
524 { //Read data from file using AliPHOSCalibrManager
525 //version and range will be used to choose proper data
527 AliPHOSCalibrationData ped("Pedestals",version);
528 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
530 Warning("ReadPedestals","Using database file 'PHOSBTCalibration.root'") ;
531 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
533 cmngr->ReadFromRoot(ped,range) ;
534 Int_t npeds=ped.NChannels() ;
535 fNch = fctdb->GetNchanels() ;
538 fhPedestals = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
539 for(Int_t i=0;i<npeds;i++){
540 Int_t raw =fctdb->AbsId2Raw(i) ;
542 fhPedestals->SetBinContent(raw-1,ped.Data(i)) ;
543 fhPedestals->SetBinError(raw-1,ped.DataCheck(i)) ;
547 //_____________________________________________________________________________
548 void AliPHOSCalibrator::ReadGains(const char * version,
550 { //Read data from file using AliPHOSCalibrManager
551 //version and range will be used to choose proper data
553 AliPHOSCalibrationData gains("Gains",version);
554 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
556 Warning("ReadGainss","Using database file 'PHOSBTCalibration.root'") ;
557 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
559 cmngr->ReadFromRoot(gains,range) ;
560 Int_t npeds=gains.NChannels() ;
561 fNch = fctdb->GetNchanels() ;
564 fhGains = new TH1F("hGainss","Gains mean",fNch,0.,fNch) ;
565 for(Int_t i=0;i<npeds;i++){
566 Int_t raw =fctdb->AbsId2Raw(i) ;
568 fhGains->SetBinContent(raw-1,gains.Data(i)) ;
569 fhGains->SetBinError(raw-1,gains.DataCheck(i)) ;
573 //_____________________________________________________________________________
574 void AliPHOSCalibrator::WriteGains(const char * version,
575 Int_t begin,Int_t end)
576 { //Write gains through AliPHOSCalibrManager
577 //version and validity range(begin-end) are used to identify data
580 Error("WriteGains","\n Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
584 AliPHOSCalibrationData gains("Gains",version);
585 for(Int_t i=0; i<fNch;i++){
586 Int_t absid=fctdb->Raw2AbsId(i) ;
587 gains.SetData(absid,fhGains->GetBinContent(i)) ;
588 gains.SetDataCheck(absid,fhGainsWid->GetBinContent(i)) ;
591 TIter next(fRunList) ;
595 while((file=((TObjString*)next()))){
596 TString s = file->GetString() ;
597 TSubString ss = s(s.Last('_'),s.Last('.'));
599 if(sscanf(ss.Data(),"%d",&tmp)){
606 gains.SetValidityRange(ibegin,iend) ;
609 gains.SetValidityRange(begin,end) ;
610 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
612 Warning("WriteGains","Using database file 'PHOSBTCalibration.root'") ;
613 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
615 cmngr->WriteData(&gains) ;
617 //_____________________________________________________________________________
618 void AliPHOSCalibrator::Print(const Option_t * option)const {
619 printf("--------------AliPHOSCalibrator-----------------\n") ;
620 printf("Files to handle:\n") ;
621 TIter next(fRunList) ;
623 while((r=(TObjString *)(next())))
624 printf(" %s\n",r->GetName()) ;
626 printf("Name of ConTableDB:.....................%s\n",fConTableDB.Data()) ;
627 printf("File of ConTableDB:.....................%s\n",fConTableDBFile.Data() ) ;
628 printf("Maximal deviation from mean Gain (factor):.%f\n",fGainAcceptCorr) ;
629 printf("Maximal deviation of Pedestal from mean:...%f\n",fAcceptCorr) ;
630 printf("Range used in Gain histos:..............%f\n",fGainMax) ;
631 printf("Number of bins in Gain histos:..........%d\n",fNGainBins) ;
632 printf("Number of channels to calibrate:........%d\n",fNch) ;
633 printf("Number of bins in pedestal histos:......%d\n",fNChan) ;
634 printf("trigger pattern for PEDESTAL events:....%d\n",fPedPat) ;
635 printf("trigger pattern for PULSER events:......%d\n",fPulPat) ;
636 printf("trigger pattern for LED events:.........%d\n",fLEDPat) ;
637 printf("trigger pattern for WIDE BEAM events:...%d\n",fWBPat) ;
638 printf("trigger pattern for NARROW BEAM events:.%d\n",fNBPat) ;
639 printf("--------------------------------------------------\n") ;