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"
52 // --- Standard library ---
54 // --- AliRoot header files ---
55 #include "AliPHOSGetter.h"
56 #include "AliPHOSCalibrator.h"
57 #include "AliPHOSConTableDB.h"
58 #include "AliPHOSCalibrManager.h"
59 #include "AliPHOSCalibrationData.h"
61 ClassImp(AliPHOSCalibrator)
64 //____________________________________________________________________________
65 AliPHOSCalibrator::AliPHOSCalibrator():TTask("AliPHOSCalibrator","Default")
67 //Default constuctor for root. Normally should not be used
75 fConTableDB = "Beamtest2002" ;
76 fConTableDBFile = "ConTableDB.root" ;
78 //____________________________________________________________________________
79 AliPHOSCalibrator::AliPHOSCalibrator(const char* file, const char* title):
80 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")
87 fRunList = new TList() ;
88 fRunList->SetOwner() ;
89 fRunList->Add(new TObjString(file)) ;
91 fPedPat = 257 ; //Patterns for different kind of events
100 fAcceptCorr = 10 ; //Maximal deviation from mean, considered as normal
102 fGainAcceptCorr = 5 ; //Factor for gain deviation
108 fConTableDB = "Beamtest2002" ;
109 fConTableDBFile = "ConTableDB.root" ;
112 //____________________________________________________________________________
113 AliPHOSCalibrator::~AliPHOSCalibrator()
123 delete fhPedestalsWid ;
129 //____________________________________________________________________________
130 void AliPHOSCalibrator::AddRun(const char * filename)
132 //Adds one more run to list of runs, which will be scanned in ScanXXX methods
134 TObjString * fn = new TObjString(filename) ;
136 fRunList=new TList() ;
137 fRunList->SetOwner() ;
142 TIter next(fRunList) ;
144 while((r=(TObjString *)(next()))){
145 if(fn->String().CompareTo(r->String())==0){
146 Error("Run already in list: ",filename) ;
154 //____________________________________________________________________________
155 void AliPHOSCalibrator::Exec(Option_t * option)
157 // reads parameters and does the calibration
158 ScanPedestals(option);
159 CalculatePedestals();
165 //____________________________________________________________________________
166 void AliPHOSCalibrator::Init(void)
168 // intializes everything
170 //check if ConTableDB already read
172 TFile * v = gROOT->GetFile(fConTableDBFile) ;
174 v = TFile::Open(fConTableDBFile) ;
176 Fatal("Can not open file with Connection Table DB:",fConTableDBFile) ;
179 fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
182 fNch = fctdb->GetNchanels() ;
183 fhPedestals = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
184 fhPedestalsWid= new TH1F("hPedestalsWid","Pedestals width",fNch,0.,fNch) ;
185 fhGains = new TH1F("hGains","Gains ",fNch,0.,fNch) ;
186 fhGainsWid = new TH1F("hGainsWid","Gains width",fNch,0.,fNch) ;
188 //____________________________________________________________________________
189 void AliPHOSCalibrator::SetConTableDB(const char * file,const char * name)
191 //Reads Connection Table database with name "name" from file "file"
193 if(file==0 || name == 0){
194 Error("Please, specify file with database"," and its title") ;
197 if(fctdb && strcmp(fctdb->GetTitle(),name)==0) //already read
206 TFile * v = gROOT->GetFile(fConTableDBFile) ;
208 v = TFile::Open(fConTableDBFile) ;
210 Error("Can not open file with Connection Table DB:",fConTableDBFile) ;
213 fctdb = dynamic_cast<AliPHOSConTableDB *>(v->Get("AliPHOSConTableDB")) ;
216 //____________________________________________________________________________
217 void AliPHOSCalibrator::PlotPedestal(Int_t chanel)
219 //Plot histogram for a given channel, filled in Scan method
220 if(fPedHistos && fPedHistos->GetEntriesFast()){
221 static_cast<TH1F*>(fPedHistos->At(chanel))->Draw() ;
224 printf("Histograms not created yet! \n") ;
227 //____________________________________________________________________________
228 void AliPHOSCalibrator::PlotPedestals(void)
230 // draws pedestals distribution
231 fhPedestals->Draw() ;
233 //____________________________________________________________________________
234 void AliPHOSCalibrator::PlotGain(Int_t chanel)
236 //Plot histogram for a given channel, filled in Scan method
237 if(fGainHistos && fGainHistos->GetEntriesFast()){
238 static_cast<TH1F*>(fGainHistos->At(chanel))->Draw() ;
241 printf("Histograms not created yet! \n") ;
244 //____________________________________________________________________________
245 void AliPHOSCalibrator::PlotGains(void)
247 // draws gains distribution
250 //____________________________________________________________________________
251 void AliPHOSCalibrator::ScanPedestals(Option_t * option )
253 //scan all files in list fRunList and fill pedestal hisgrams
254 //option: "clear" - clear pedestal histograms filled up to now
255 // "deb" - plot file name currently processed
260 if(fPedHistos && strstr(option,"clear"))
261 fPedHistos->Delete() ;
263 fPedHistos = new TObjArray(fNch) ;
265 //Create histos for each channel, fills them and extracts mean values.
266 //First - prepare histos
268 for(ich=0;ich<fNch ;ich++){
269 TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
273 TString name("Pedestal for channel ") ;
275 fPedHistos->AddAt(new TH1F(n,name,fNChan,0,fNChan),ich) ;
279 TIter next(fRunList) ;
281 while((file = static_cast<TObjString *>(next()))){
282 if(strstr(option,"deb"))
283 printf("Processing file %s \n ",file->String().Data()) ;
284 AliPHOSGetter * gime = AliPHOSGetter::Instance(file->String().Data(),GetTitle()) ;
286 for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
287 gime->Event(ievent,"D") ;
288 if(gime->EventPattern() == fPedPat ){
290 TClonesArray * digits = gime->Digits() ;
291 for(idigit = 0; idigit<digits->GetEntriesFast() ; idigit++){
292 AliPHOSDigit * digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
293 ich = fctdb->AbsId2Raw(digit->GetId());
296 Float_t amp = digit->GetAmp() ;
297 TH1F * hh = dynamic_cast<TH1F*>(fPedHistos->At(ich)) ;
306 //____________________________________________________________________________
307 void AliPHOSCalibrator::CalculatePedestals()
309 //Fit histograms, filled in ScanPedestals method with Gaussian
310 //find mean and width, check deviation from mean for each channel.
312 if(!fPedHistos || !fPedHistos->At(0)){
313 Error("CalculatePedestals","You should run ScanPedestals first!") ;
317 //Now fit results with Gauss
318 TF1 * gs = new TF1("gs","gaus",0.,10000.) ;
320 for(ich=0;ich<fNch ;ich++){
321 TH1F * h = static_cast<TH1F *>(fPedHistos->At(ich)) ;
322 Int_t max = h->GetMaximumBin() ;
323 Axis_t xmin = max/2. ;
324 Axis_t xmax = max*3/2 ;
325 gs->SetRange(xmin,xmax) ;
327 par[0] = h->GetBinContent(max) ;
330 gs->SetParameters(par[0],par[1],par[2]) ;
332 gs->GetParameters(par) ;
333 fhPedestals->SetBinContent(ich,par[1]) ;
334 fhPedestals->SetBinError(ich,par[2]) ;
335 fhPedestalsWid->Fill(ich,par[2]) ;
339 //now check reasonability of results
340 TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
341 fhPedestals->Fit("p0","Q") ;
343 p0->GetParameters(&meanPed);
344 for(ich=0;ich<fNch ;ich++){
345 Float_t ped = fhPedestals->GetBinContent(ich) ;
346 if(ped < 0 || ped > meanPed+fAcceptCorr){
347 TString out("Pedestal of channel ") ;
351 out+= "it is too far from mean " ;
353 Error("PHOSCalibrator",out) ;
359 //____________________________________________________________________________
360 void AliPHOSCalibrator::ScanGains(Option_t * option)
362 //Scan all runs, listed in fRunList and fill histograms for all channels
363 //options: "clear" - clean histograms, filled up to now
364 // "deb" - print current file name
365 // "narrow" - scan only narrow beam events
369 if(fGainHistos && strstr(option,"clear"))
370 fGainHistos->Delete() ;
372 if(strstr(option,"deball"))
373 printf("creating array for %d channels \n",fNch) ;
374 fGainHistos = new TObjArray(fNch) ;
377 //Create histos for each channel, fills them and extracts mean values.
378 //First - prepare histos
380 if(!fGainHistos->GetEntriesFast()){
382 for(ich=0;ich<fNch ;ich++){
385 TString name("Gains for channel ") ;
387 fGainHistos->AddAt(new TH1F(n,name,fNGainBins,0,fGainMax),ich) ;
388 // static_cast<TH1F*>(fGainHistos->At(ich))->Sumw2() ;
392 Bool_t all =!(Bool_t)strstr(option,"narrow") ;
395 TIter next(fRunList) ;
397 while((file = static_cast<TObjString *>(next()))){
399 AliPHOSGetter * gime = AliPHOSGetter::Instance(file->String().Data(),GetTitle()) ;
402 for(ievent = 0; ievent<gime->MaxEvent() ; ievent++){
403 gime->Event(ievent,"D") ;
404 if(gime->EventPattern() == fNBPat ||
405 (all && gime->EventPattern() == fWBPat)){
408 TClonesArray * digits = gime->Digits() ;
409 AliPHOSDigit * digit ;
412 for(idigit = 0; idigit<digits->GetEntriesFast() ; idigit++){
413 digit = static_cast<AliPHOSDigit *>(digits->At(idigit) ) ;
414 if(digit->GetAmp() > max){
416 max = digit->GetAmp() ;
419 digit = static_cast<AliPHOSDigit *>(digits->At(imax) ) ;
420 Int_t ich = fctdb->AbsId2Raw(digit->GetId());
422 Float_t pedestal = fhPedestals->GetBinContent(ich) ;
423 const Float_t kshowerInCrystall = 0.9 ;
424 Float_t beamEnergy = gime->BeamEnergy() ;
425 Float_t gain = beamEnergy*kshowerInCrystall/
426 (digit->GetAmp() - pedestal) ;
427 static_cast<TH1F*>(fGainHistos->At(ich))->Fill(gain) ;
431 if(strstr(option,"deb"))
432 printf("Hadled %d events \n",handled) ;
435 //____________________________________________________________________________
436 void AliPHOSCalibrator::CalculateGains(void)
440 if(!fGainHistos || !fGainHistos->GetEntriesFast()){
441 Error("CalculateGains","You should run ScanGains first!") ;
445 //Fit results with Landau
446 TF1 * gs = new TF1("gs","landau",0.,10000.) ;
448 for(ich=0;ich<fNch ;ich++){
449 TH1F * h = static_cast<TH1F *>(fGainHistos->At(ich)) ;
450 Int_t bmax = h->GetMaximumBin() ;
451 Axis_t center = h->GetBinCenter(bmax) ;
452 Axis_t xmin = center - 0.01 ;
453 Axis_t xmax = center + 0.02 ;
454 gs->SetRange(xmin,xmax) ;
456 par[0] = h->GetBinContent(bmax) ;
459 gs->SetParameters(par[0],par[1],par[2]) ;
461 gs->GetParameters(par) ;
462 fhGains->SetBinContent(ich,par[1]) ;
463 fhGains->SetBinError(ich,par[2]) ;
464 fhGainsWid->Fill(ich,par[2]) ;
468 //now check reasonability of results
469 TF1 * p0 = new TF1("p0","pol0",0.,fNch) ;
470 fhGains->Fit("p0","Q") ;
472 p0->GetParameters(&meanGain);
473 for(ich=0;ich<fNch ;ich++){
474 Float_t gain = fhGains->GetBinContent(ich) ;
475 if(gain < meanGain/fGainAcceptCorr || gain > meanGain*fGainAcceptCorr){
476 TString out("Gain of channel ") ;
480 out+= "it is too far from mean " ;
482 Error("PHOSCalibrator",out) ;
489 //_____________________________________________________________________________
490 void AliPHOSCalibrator::WritePedestals(const char * version,
491 Int_t begin,Int_t end)
493 //Write calculated data to file using AliPHOSCalibrManager
494 //version and validitirange (begin-end) will be used to identify data
497 Error("WritePedestals","\n Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
501 AliPHOSCalibrationData ped("Pedestals",version);
502 for(Int_t i=0; i<fNch;i++){
503 Int_t absid=fctdb->Raw2AbsId(i) ;
504 ped.SetData(absid,fhPedestals->GetBinContent(i)) ;
505 ped.SetDataCheck(absid,fhPedestalsWid->GetBinContent(i)) ;
508 //evaluate validity range
510 TIter next(fRunList) ;
514 while((file=((TObjString*)next()))){
515 TString s = file->GetString() ;
516 TString ss = s(s.Last('_'),s.Last('.'));
518 if(sscanf(ss.Data(),"%d",&tmp)){
525 ped.SetValidityRange(ibegin,iend) ;
528 ped.SetValidityRange(begin,end) ;
530 //check, may be Manager instance already configured?
531 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
533 Warning("Write Pedestals","Using database file 'PHOSBTCalibration.root'") ;
534 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
536 cmngr->WriteData(&ped) ;
538 //_____________________________________________________________________________
539 void AliPHOSCalibrator::ReadPedestals(const char * version,
542 //Read data from file using AliPHOSCalibrManager
543 //version and range will be used to choose proper data
545 AliPHOSCalibrationData ped("Pedestals",version);
546 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
548 Warning("ReadPedestals","Using database file 'PHOSBTCalibration.root'") ;
549 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
551 cmngr->ReadFromRoot(ped,range) ;
552 Int_t npeds=ped.NChannels() ;
553 fNch = fctdb->GetNchanels() ;
556 fhPedestals = new TH1F("hPedestals","Pedestals mean",fNch,0.,fNch) ;
557 for(Int_t i=0;i<npeds;i++){
558 Int_t raw =fctdb->AbsId2Raw(i) ;
560 fhPedestals->SetBinContent(raw-1,ped.Data(i)) ;
561 fhPedestals->SetBinError(raw-1,ped.DataCheck(i)) ;
565 //_____________________________________________________________________________
566 void AliPHOSCalibrator::ReadGains(const char * version,
569 //Read data from file using AliPHOSCalibrManager
570 //version and range will be used to choose proper data
572 AliPHOSCalibrationData gains("Gains",version);
573 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
575 Warning("ReadGainss","Using database file 'PHOSBTCalibration.root'") ;
576 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
578 cmngr->ReadFromRoot(gains,range) ;
579 Int_t npeds=gains.NChannels() ;
580 fNch = fctdb->GetNchanels() ;
583 fhGains = new TH1F("hGainss","Gains mean",fNch,0.,fNch) ;
584 for(Int_t i=0;i<npeds;i++){
585 Int_t raw =fctdb->AbsId2Raw(i) ;
587 fhGains->SetBinContent(raw-1,gains.Data(i)) ;
588 fhGains->SetBinError(raw-1,gains.DataCheck(i)) ;
592 //_____________________________________________________________________________
593 void AliPHOSCalibrator::WriteGains(const char * version,
594 Int_t begin,Int_t end)
596 //Write gains through AliPHOSCalibrManager
597 //version and validity range(begin-end) are used to identify data
600 Error("WriteGains","\n Please, supply Connection Table DB (use SetConTableDB()) \n" ) ;
604 AliPHOSCalibrationData gains("Gains",version);
605 for(Int_t i=0; i<fNch;i++){
606 Int_t absid=fctdb->Raw2AbsId(i) ;
607 gains.SetData(absid,fhGains->GetBinContent(i)) ;
608 gains.SetDataCheck(absid,fhGainsWid->GetBinContent(i)) ;
611 TIter next(fRunList) ;
615 while((file=((TObjString*)next()))){
616 TString s = file->GetString() ;
617 TSubString ss = s(s.Last('_'),s.Last('.'));
619 if(sscanf(ss.Data(),"%d",&tmp)){
626 gains.SetValidityRange(ibegin,iend) ;
629 gains.SetValidityRange(begin,end) ;
630 AliPHOSCalibrManager * cmngr = AliPHOSCalibrManager::GetInstance() ;
632 Warning("WriteGains","Using database file 'PHOSBTCalibration.root'") ;
633 cmngr = AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
635 cmngr->WriteData(&gains) ;
637 //_____________________________________________________________________________
638 void AliPHOSCalibrator::Print()const
641 printf("--------------AliPHOSCalibrator-----------------\n") ;
642 printf("Files to handle:\n") ;
643 TIter next(fRunList) ;
645 while((r=(TObjString *)(next())))
646 printf(" %s\n",r->GetName()) ;
648 printf("Name of ConTableDB:.....................%s\n",fConTableDB.Data()) ;
649 printf("File of ConTableDB:.....................%s\n",fConTableDBFile.Data() ) ;
650 printf("Maximal deviation from mean Gain (factor):.%f\n",fGainAcceptCorr) ;
651 printf("Maximal deviation of Pedestal from mean:...%f\n",fAcceptCorr) ;
652 printf("Range used in Gain histos:..............%f\n",fGainMax) ;
653 printf("Number of bins in Gain histos:..........%d\n",fNGainBins) ;
654 printf("Number of channels to calibrate:........%d\n",fNch) ;
655 printf("Number of bins in pedestal histos:......%d\n",fNChan) ;
656 printf("trigger pattern for PEDESTAL events:....%d\n",fPedPat) ;
657 printf("trigger pattern for PULSER events:......%d\n",fPulPat) ;
658 printf("trigger pattern for LED events:.........%d\n",fLEDPat) ;
659 printf("trigger pattern for WIDE BEAM events:...%d\n",fWBPat) ;
660 printf("trigger pattern for NARROW BEAM events:.%d\n",fNBPat) ;
661 printf("--------------------------------------------------\n") ;