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 **************************************************************************/
18 //_________________________________________________________________________
19 // Class intended to perform online monitoring of PHOS beamtests
20 // Being constructed, produces menu with list of available histograms to fill
21 // Once histograms are selected, button "Go" should be pressed to start scan of data.
22 // Prepared histograms will be periodically updated during scan of the data.
24 // 1. To plot most of histograms, a "Connection Table", relating ADC signal index and AbsId
25 // of PHOS crystal, should be created beforehand. To do this, call macro
26 // $ALICE_ROOT/PHOS/macros/BeamTest/MakeConTableDB.C
27 // with apropriate number of raws and columns of prototype.
28 // 2. To perform reconstruction and e.g. invariant mass analysis, a "Calibration Database"
29 // should be created beforehand. To do this, call macro
30 // $ALICE_ROOT/PHOS/macros/BeamTest/MakeConTableDB.C
31 // to read calibration parameters from file or use AliPHOSCalibrator to calculate
32 // pedestal and gains.
33 // 3. Once histograms are filled with "Go" method, they can be written to file
34 // with WriteHisto("Filename.root") method.
36 //*-- Author : D.Peressounko (RRC KI) after A.V. Kuryakin, (Sarov)
37 //////////////////////////////////////////////////////////////////////////////
39 // --- ROOT system ---
41 #include "TClonesArray.h"
45 #include "TGroupButton.h"
50 // --- Standard library ---
51 #include "TBenchmark.h"
52 #include "Riostream.h"
54 // --- AliRoot header files ---
55 #include "AliPHOSOnlineMonitor.h"
56 #include "AliPHOSConTableDB.h"
57 #include "AliPHOSGeometry.h"
58 #include "AliPHOSRawReaderDate.h"
59 #include "AliPHOSRawStream.h"
60 #include "AliPHOSDigit.h"
61 #include "AliPHOSGetterLight.h"
62 #include "AliPHOSClusterizerv1.h"
63 #include "AliPHOSTrackSegmentMakerv1.h"
64 #include "AliPHOSPIDv1.h"
65 #include "AliPHOSCalibrManager.h"
66 #include "AliPHOSCalibrationDB.h"
68 ClassImp(AliPHOSOnlineMonitor)
71 //____________________________________________________________________________
72 AliPHOSOnlineMonitor::AliPHOSOnlineMonitor(): TDialogCanvas("PHOS","PHOS",150,300)
80 //add this TFitPanel to the list of cleanups such that in case
81 //the referenced object is deleted, its pointer be reset
82 gROOT->GetListOfCleanups()->Add(this);
85 fRefPad = (TPad*)gROOT->GetSelectedPad();
91 //Prepare list of canvas with histograms
92 fCanvasList = new TList() ;
93 fHistosList = new TList() ;
96 fGeom = AliPHOSGeometry::GetInstance("IHEP","") ;
100 //____________________________________________________________________________
101 AliPHOSOnlineMonitor::AliPHOSOnlineMonitor(const char * inputfile): TDialogCanvas("PHOS","PHOS",150,300)
109 //add this TFitPanel to the list of cleanups such that in case
110 //the referenced object is deleted, its pointer be reset
111 gROOT->GetListOfCleanups()->Add(this);
114 fRefPad = (TPad*)gROOT->GetSelectedPad();
120 //Prepare list of canvas with histograms
121 fCanvasList = new TList() ;
122 fHistosList = new TList() ;
124 fInputFile = inputfile ;
125 fGeom = AliPHOSGeometry::GetInstance("IHEP","") ;
130 //____________________________________________________________________________
131 AliPHOSOnlineMonitor::~AliPHOSOnlineMonitor()
133 //Obvious, but unevoidable comment for destructor: cleans up everething.
134 TIter nextCanvas(fCanvasList);
136 while((c=(TCanvas*)nextCanvas()))
140 TIter nextHisto(fHistosList);
142 while((h=(TH1D*)nextHisto()))
150 //____________________________________________________________________________
151 void AliPHOSOnlineMonitor::MakeButtons(void){
152 //Make buttons on graphical menu
159 Float_t dy = (ymax-ymin)/nButtons ;
163 b = new TGroupButton("APPLY","Triggers","",xmin,y1,xmax,y2);
167 b = new TGroupButton("APPLY","Pedestals","",xmin,y1,xmax,y2);
171 b = new TGroupButton("APPLY","Spectrum all","",xmin,y1,xmax,y2);
175 b = new TGroupButton("APPLY","Spectrum g","",xmin,y1,xmax,y2);
179 b = new TGroupButton("APPLY","Inv Mass","",xmin,y1,xmax,y2);
181 for(Int_t i=1; i<=5; i++){
185 sprintf(name,"Edep(ADC) %d",i) ;
186 b = new TGroupButton("APPLY",name,"",xmin,y1,xmax,y2);
189 for(Int_t i=1; i<=5; i++){
193 sprintf(name,"Edep(Cal) %d",i) ;
194 b = new TGroupButton("APPLY",name,"",xmin,y1,xmax,y2);
199 b = new TGroupButton("APPLY","Go","",xmin,y1,xmax,y2);
203 //____________________________________________________________________________
204 void AliPHOSOnlineMonitor::Apply(const char *action){
205 //Function to handle button actions
207 TDialogCanvas::Apply() ;
210 TGroupButton *button;
211 TIter next(fPrimitives);
213 if (!strcmp(action,"Triggers")) {
217 if (!strcmp(action,"Pedestals")) {
220 if (!strcmp(action,"Spectrum all")) {
221 DrawSpectrum("all") ;
223 if (!strcmp(action,"Spectrum g")) {
224 DrawSpectrum("gamma") ;
226 if (!strcmp(action,"Inv Mass")) {
229 if(strstr(action,"Edep")){
232 sscanf(action,"%s %d",tmp,&n) ;
234 if(strstr(action,"Cal"))
238 if (!strcmp(action,"Go")) {
242 //Mark button as pressed
243 if(strcmp(action,"Go")){ //Do not mark "Go" button
244 while ((obj = next())) {
245 if (obj->InheritsFrom(TGroupButton::Class())) {
246 button = (TGroupButton*)obj;
247 if(!strcmp(button->GetTitle(),action)){
248 if (button->GetBorderMode() > 0){
249 button->SetBorderMode(-1) ;
250 button->Modified(kTRUE);
257 //____________________________________________________________________________
258 void AliPHOSOnlineMonitor::SetInputFile(const char * filename){
259 //close previously opened
261 fInputFile = filename ;
263 //____________________________________________________________________________
264 void AliPHOSOnlineMonitor::DrawPedestals(){
265 //Prepare canvas and histograms for drawing pedestals
267 TIter nextCanvas(fCanvasList);
269 Bool_t exists = kFALSE ;
270 while((c=(TCanvas*)nextCanvas())){
271 if(!strcmp(c->GetName(),"Pedestals")){
277 c = new TDialogCanvas("Pedestals","Pedestals",300,200) ;
278 fCanvasList->AddLast(c) ;
281 TIter nextHisto(fHistosList);
284 while((h=(TH1D*)nextHisto())){
285 if(!strcmp(h->GetName(),"hPedestals")){
291 h = new TH1D("hPedestals","Pedestals per event",fGeom->GetNModules()*fGeom->GetNCristalsInModule(),0.,
292 1.*fGeom->GetNModules()*fGeom->GetNCristalsInModule()) ;
293 fHistosList->AddLast(h) ;
299 fScanPed = kTRUE ; //We will scan pedestals
302 //____________________________________________________________________________
303 void AliPHOSOnlineMonitor::DrawTriggers(){
304 //Prepare canvas and histogram for drawing triggers
306 TIter nextCanvas(fCanvasList);
308 Bool_t exists = kFALSE ;
309 while((c=(TCanvas*)nextCanvas())){
310 if(!strcmp(c->GetName(),"Triggers")){
316 c = new TDialogCanvas("Triggers","Triggers",200,200) ;
317 fCanvasList->AddLast(c) ;
320 TIter nextHisto(fHistosList);
323 while((h=(TH1D*)nextHisto())){
324 if(!strcmp(h->GetName(),"hTriggers")){
330 h = new TH1D("hTriggers","Triggers",2,0.,2.) ;
331 fHistosList->AddLast(h) ;
334 h->SetBit(TH1::kCanRebin);
335 h->Fill("LED",0.0000001) ;
336 h->Fill("PUL",0.0000001) ;
337 h->Fill("PED",0.0000001) ;
338 h->Fill("NEL",0.0000001) ;
339 h->Fill("WEL",0.0000001) ;
340 h->Fill("SOB",0.0000001) ;
341 h->Fill("EOB",0.0000001) ;
342 h->Fill("wrong",0.0000001) ;
343 h->LabelsOption("h");
349 //____________________________________________________________________________
350 void AliPHOSOnlineMonitor::DrawSpectrum(const char * opt){
351 //Prepare canvas and histograms for drawing spectra of all reconstructed particles or photons
353 TString name("Spectrum") ;
356 TIter nextCanvas(fCanvasList);
358 Bool_t exists = kFALSE ;
359 while((c=(TCanvas*)nextCanvas())){
360 if(!strcmp(c->GetName(),name.Data())){
366 c = new TDialogCanvas(name,name,250,300) ;
367 fCanvasList->AddLast(c) ;
370 TIter nextHisto(fHistosList);
374 while((h=(TH1D*)nextHisto())){
375 if(!strcmp(h->GetName(),name.Data())){
381 h = new TH1D(name,name,100,0.,100.) ;
382 fHistosList->AddLast(h) ;
388 fReconstruct = kTRUE ;
389 fScanSig = kTRUE ; //We will scan pedestals
392 //____________________________________________________________________________
393 void AliPHOSOnlineMonitor::DrawMinv(){
394 TIter nextCanvas(fCanvasList);
396 Bool_t exists = kFALSE ;
397 while((c=(TCanvas*)nextCanvas())){
398 if(!strcmp(c->GetName(),"InvMass")){
404 c = new TDialogCanvas("InvMass","Invariant mass",300,200) ;
405 fCanvasList->AddLast(c) ;
408 TIter nextHisto(fHistosList);
411 while((h=(TH1D*)nextHisto())){
412 if(!strcmp(h->GetName(),"hInvMass")){
418 h = new TH1D("hInvMass","hInvMass",1000,0.,1.0) ;
419 fHistosList->AddLast(h) ;
425 fReconstruct = kTRUE ;
426 fScanSig = kTRUE ; //We will scan pedestals
428 //____________________________________________________________________________
429 void AliPHOSOnlineMonitor::DrawEdep(Int_t mod,const char * opt){
431 sprintf(name,"Edep%s %d",opt,mod) ;
433 TIter nextCanvas(fCanvasList);
435 Bool_t exists = kFALSE ;
436 while((c=(TCanvas*)nextCanvas())){
437 if(!strcmp(c->GetName(),name)){
443 c = new TDialogCanvas(name,name,300,200) ;
444 fCanvasList->AddLast(c) ;
447 TIter nextHisto(fHistosList);
450 sprintf(name,"hEdep%s%d",opt,mod) ;
451 while((h=(TH2D*)nextHisto())){
452 if(!strcmp(h->GetName(),name)){
458 h = new TH2D(name,name,fGeom->GetNPhi(),0.,1.*fGeom->GetNPhi(),fGeom->GetNZ(),0.,1.*fGeom->GetNZ()) ;
459 fHistosList->AddLast(h) ;
465 fScanSig = kTRUE ; //We will scan signal events
466 if(strstr(opt,"Cal"))
467 fReconstruct = kTRUE ;
470 //____________________________________________________________________________
471 void AliPHOSOnlineMonitor::ScanPedestals(TClonesArray * digits){
472 //This method is called for events with PED trigger
473 //We fill bins with ADC values
475 TH1D * h = (TH1D*)gROOT->FindObjectAny("hPedestals");
477 Error("ScanPedestals","Can not fild histogram hPedestals") ;
480 for(Int_t i=0; i<digits->GetEntriesFast(); i++){
481 AliPHOSDigit * dig = static_cast<AliPHOSDigit*>(digits->At(i)) ;
482 h->AddBinContent(dig->GetId(),dig->GetAmp()) ;
485 //____________________________________________________________________________
486 void AliPHOSOnlineMonitor::ScanEdep(TClonesArray * digits){
487 //Fill 2D distribution of ADC values in NEL and WEL events
488 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
489 AliPHOSCalibrationDB *cdb = 0 ;
491 cdb = gime->CalibrationDB() ;
496 for(Int_t i=0; i<digits->GetEntriesFast(); i++){
497 AliPHOSDigit * dig = static_cast<AliPHOSDigit*>(digits->At(i)) ;
499 fGeom->AbsToRelNumbering(dig->GetId(),relId) ;
500 if(mod != relId[0]){ //new module, look for histograms
502 sprintf(name,"hEdep%d",mod) ;
503 h = (TH2D*)gROOT->FindObjectAny(name);
504 sprintf(name,"hEdepCal%d",mod) ;
505 hCal = (TH2D*)gROOT->FindObjectAny(name);
508 h->Fill(relId[2]-0.1,relId[3]-0.1,1.*dig->GetAmp()) ;
510 hCal->Fill(relId[2]-0.1,relId[3]-0.1,cdb->Calibrate(dig->GetAmp(),dig->GetId())) ;
513 //____________________________________________________________________________
514 void AliPHOSOnlineMonitor::ScanRecon(TClonesArray * recParticles){
515 if(!recParticles || recParticles->GetEntries()==0) return ;
517 TH1D* hSpectr = (TH1D*)gROOT->FindObjectAny("hSpectrumall");
518 TH1D* hSpectrGam= (TH1D*)gROOT->FindObjectAny("hSpectrumgamma");
519 TH1D* hInvMass = (TH1D*)gROOT->FindObjectAny("hInvMass");
520 for(Int_t i=0; i<recParticles->GetEntriesFast() ; i++){
521 AliPHOSRecParticle * p = (AliPHOSRecParticle *)recParticles->At(i) ;
522 if(hSpectr)hSpectr->Fill(p->Energy()) ;
523 if(hSpectrGam && p->IsPhoton())hSpectrGam->Fill(p->Energy()) ;
525 for(Int_t j=i+1; j<recParticles->GetEntriesFast() ; j++){
526 AliPHOSRecParticle * p2 = (AliPHOSRecParticle *)recParticles->At(j) ;
527 Double_t e = p->Energy() + p2->Energy() ;
528 Double_t x = p->Px() + p2->Px() ;
529 Double_t y = p->Py() + p2->Py() ;
530 Double_t z = p->Pz() + p2->Pz() ;
531 Double_t m = e*e-x*x-y*y-z*z ;
532 hInvMass->Fill(m>0?TMath::Sqrt(m): 0. ) ;
537 //____________________________________________________________________________
538 void AliPHOSOnlineMonitor::ScanTrigger(Int_t trig){
539 //Fills trigger distribution
541 TH1D * h = (TH1D*)gROOT->FindObjectAny("hTriggers");
544 case AliPHOSRawStream::kLED : h->Fill("LED",1.) ; break ;
545 case AliPHOSRawStream::kPUL : h->Fill("PUL",1.) ; break ;
546 case AliPHOSRawStream::kPED : h->Fill("PED",1.) ; break ;
547 case AliPHOSRawStream::kNEL : h->Fill("NEL",1.) ; break ;
548 case AliPHOSRawStream::kWEL : h->Fill("WEL",1.) ; break ;
549 case AliPHOSRawStream::kSOB : h->Fill("SOB",1.) ; break ;
550 case AliPHOSRawStream::kEOB : h->Fill("EOB",1.) ; break ;
551 default : h->Fill("wrong",1.) ;
554 //____________________________________________________________________________
555 void AliPHOSOnlineMonitor::SetConTableDB(const char * filename){
556 //Read ConnectionTableDB from file
557 TFile * file = new TFile(filename) ;
558 AliPHOSConTableDB * tmp = (AliPHOSConTableDB*)file->Get("AliPHOSConTableDB") ;
559 fcdb = new AliPHOSConTableDB(*tmp) ;
562 //____________________________________________________________________________
563 void AliPHOSOnlineMonitor::Go(){
564 //Perform scan of curent event
565 gBenchmark->Start("PHOSOnlineMon");
567 //First test if we need "Connection table" then open it
571 Info("Go","Read Connection table from file \"ConTableDB.root\"") ;
573 Error("Go","Please, set connection table with SetConTableDB() method") ;
578 AliPHOSGetterLight * gime = AliPHOSGetterLight::Instance("PHOS","On Flight") ;
580 //Configure CalibrManager to read data from file
581 //Create calibration database and read it
582 AliPHOSCalibrationDB * calibDB = 0 ;
583 if(fScanSig || fReconstruct){ //We will ned calibration parameters
584 AliPHOSCalibrManager::GetInstance("CalibrDB.root","root") ;
585 //If we configured manager to read from ASCII file,
586 //give him connection table. OK, it will not harm in any case.
587 AliPHOSCalibrManager::GetInstance()->SetConTable(fcdb) ;
589 calibDB = new AliPHOSCalibrationDB("OnLine") ;
590 calibDB->GetParameters() ; //Read parameters using Manager
591 gime->SetCalibrationDB(calibDB) ;
595 AliPHOSRawReaderDate *rawReader = new AliPHOSRawReaderDate(fInputFile) ;
596 AliPHOSRawStream *rawStream = new AliPHOSRawStream(rawReader) ;
597 rawStream->SetConTableDB(fcdb) ;
599 TClonesArray * digits = gime->Digits() ;
600 TClonesArray * recParticles = gime->RecParticles() ;
601 AliPHOSClusterizerv1* clu = 0 ;
602 AliPHOSTrackSegmentMakerv1 * tsm = 0 ;
603 AliPHOSPIDv1 * pid = 0 ;
604 if(fReconstruct){ //We will need calibation parameters
605 clu = new AliPHOSClusterizerv1("PHOS","On Flight") ;
606 clu->SetWriting(0) ; //Do not write to file
607 clu->SetEmcMinE(0.05) ; //Minimal energy of the digit
608 clu->SetEmcLocalMaxCut(0.05) ; //Height of local maximum over environment
609 clu->SetEmcClusteringThreshold(0.2) ; //Minimal energy to start cluster
610 // clu->SetUnfolding(kFALSE) ; //Do not unfold
611 tsm = new AliPHOSTrackSegmentMakerv1("PHOS","On Flight") ;
612 tsm->SetWriting(0) ; //Do not write to file
613 pid = new AliPHOSPIDv1("PHOS","On Flight") ;
614 pid->SetWriting(0) ; //Do not write to file
618 //Scan all event in file
620 while(rawReader->NextEvent()){
621 //Is it PHYSICAL event
622 if(rawReader->GetType() == PHYSICS_EVENT){
624 if(fNevents%100 ==0){
625 printf("\b\b\b\b\b\b%6d",fNevents) ;
627 if(rawStream->ReadDigits(digits)){
631 ScanTrigger(rawStream->GetTrigger()) ;
632 if(rawStream->IsPEDevent() && fScanPed){
633 ScanPedestals(digits) ;
635 if((rawStream->IsNELevent() || rawStream->IsWELevent()) && fScanSig){
638 gime->Clusterizer()->Exec("") ;
639 gime->TrackSegmentMaker()->Exec("") ;
640 gime->PID()->Exec("") ;
641 ScanRecon(recParticles) ;
646 if(fNevents%fNUpdate == 0 ){ //upate all histograms
647 TIter nextCanvas(fCanvasList);
649 while((c=(TCanvas*)nextCanvas())){
654 gSystem->ProcessEvents();
656 // if(fNevents>=200)break ;
659 gBenchmark->Stop("PHOSOnlineMon");
660 Float_t time = gBenchmark->GetCpuTime("PHOSOnlineMon") ;
661 printf("took %f seconds for scanning, i.e. %f seconds per event %d \n",
662 time,time/fNevents,fNevents) ;
664 //Update canvas with histograms at the end
665 TIter nextCanvas(fCanvasList);
667 while((c=(TCanvas*)nextCanvas())){
674 printf("delete 1 \n") ;
675 if(calibDB) delete calibDB ;
679 //____________________________________________________________________________
680 void AliPHOSOnlineMonitor::Clean(){
681 //Cleans content of all histograms
683 TIter nextHisto(fHistosList);
685 while((h=(TH1D*)nextHisto())){
688 TIter nextCanvas(fCanvasList);
690 while((c=(TCanvas*)nextCanvas())){
694 //____________________________________________________________________________
695 void AliPHOSOnlineMonitor::Reset(){
696 //delets all canvas and histograms,
697 //marks buttons as unpressed
699 TIter nextHisto(fHistosList);
701 while((h=(TH1D*)nextHisto())){
702 fHistosList->Remove(h) ;
705 TIter nextCanvas(fCanvasList);
707 while((c=(TCanvas*)nextCanvas())){
708 fCanvasList->Remove(c) ;
712 TGroupButton *button;
713 TIter next(fPrimitives);
715 //Mark buttons as anpressed
716 while ((obj = next())) {
717 if (obj->InheritsFrom(TGroupButton::Class())) {
718 button = (TGroupButton*)obj;
719 if (button->GetBorderMode() < 0){
720 button->SetBorderMode(1) ;
721 button->Modified(kTRUE);
727 //____________________________________________________________________________
728 void AliPHOSOnlineMonitor::WriteHistograms(const char * filename){
729 //Write filled histograms to file
730 TFile * file = new TFile(filename,"Update") ;
732 TIter nextHisto(fHistosList);
734 while((h=(TH1*)nextHisto())){
735 h->Write(0,TObject::kOverwrite) ;