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 "AliRawReaderDateV3.h"
59 #include "AliRawEventHeaderBase.h"
60 #include "AliPHOSRawStream.h"
61 #include "AliPHOSDigit.h"
62 #include "AliPHOSGetterLight.h"
63 #include "AliPHOSClusterizerv1.h"
64 #include "AliPHOSTrackSegmentMakerv1.h"
65 #include "AliPHOSPIDv1.h"
66 #include "AliPHOSCalibrManager.h"
67 #include "AliPHOSCalibrationDB.h"
69 ClassImp(AliPHOSOnlineMonitor)
72 //____________________________________________________________________________
73 AliPHOSOnlineMonitor::AliPHOSOnlineMonitor(): TDialogCanvas("PHOS","PHOS",150,300)
81 //add this TFitPanel to the list of cleanups such that in case
82 //the referenced object is deleted, its pointer be reset
83 gROOT->GetListOfCleanups()->Add(this);
86 fRefPad = (TPad*)gROOT->GetSelectedPad();
92 //Prepare list of canvas with histograms
93 fCanvasList = new TList() ;
94 fHistosList = new TList() ;
97 fGeom = AliPHOSGeometry::GetInstance("IHEP","") ;
101 //____________________________________________________________________________
102 AliPHOSOnlineMonitor::AliPHOSOnlineMonitor(const char * inputfile): TDialogCanvas("PHOS","PHOS",150,300)
110 //add this TFitPanel to the list of cleanups such that in case
111 //the referenced object is deleted, its pointer be reset
112 gROOT->GetListOfCleanups()->Add(this);
115 fRefPad = (TPad*)gROOT->GetSelectedPad();
121 //Prepare list of canvas with histograms
122 fCanvasList = new TList() ;
123 fHistosList = new TList() ;
125 fInputFile = inputfile ;
126 fGeom = AliPHOSGeometry::GetInstance("IHEP","") ;
131 //____________________________________________________________________________
132 AliPHOSOnlineMonitor::~AliPHOSOnlineMonitor()
134 //Obvious, but unevoidable comment for destructor: cleans up everething.
135 TIter nextCanvas(fCanvasList);
137 while((c=(TCanvas*)nextCanvas()))
141 TIter nextHisto(fHistosList);
143 while((h=(TH1D*)nextHisto()))
151 //____________________________________________________________________________
152 void AliPHOSOnlineMonitor::MakeButtons(void){
153 //Make buttons on graphical menu
160 Float_t dy = (ymax-ymin)/nButtons ;
164 b = new TGroupButton("APPLY","Triggers","",xmin,y1,xmax,y2);
168 b = new TGroupButton("APPLY","Pedestals","",xmin,y1,xmax,y2);
172 b = new TGroupButton("APPLY","Spectrum all","",xmin,y1,xmax,y2);
176 b = new TGroupButton("APPLY","Spectrum g","",xmin,y1,xmax,y2);
180 b = new TGroupButton("APPLY","Inv Mass","",xmin,y1,xmax,y2);
182 for(Int_t i=1; i<=5; i++){
186 sprintf(name,"Edep(ADC) %d",i) ;
187 b = new TGroupButton("APPLY",name,"",xmin,y1,xmax,y2);
190 for(Int_t i=1; i<=5; i++){
194 sprintf(name,"Edep(Cal) %d",i) ;
195 b = new TGroupButton("APPLY",name,"",xmin,y1,xmax,y2);
200 b = new TGroupButton("APPLY","Go","",xmin,y1,xmax,y2);
204 //____________________________________________________________________________
205 void AliPHOSOnlineMonitor::Apply(const char *action){
206 //Function to handle button actions
208 TDialogCanvas::Apply() ;
211 TGroupButton *button;
212 TIter next(fPrimitives);
214 if (!strcmp(action,"Triggers")) {
218 if (!strcmp(action,"Pedestals")) {
221 if (!strcmp(action,"Spectrum all")) {
222 DrawSpectrum("all") ;
224 if (!strcmp(action,"Spectrum g")) {
225 DrawSpectrum("gamma") ;
227 if (!strcmp(action,"Inv Mass")) {
230 if(strstr(action,"Edep")){
233 sscanf(action,"%s %d",tmp,&n) ;
235 if(strstr(action,"Cal"))
239 if (!strcmp(action,"Go")) {
243 //Mark button as pressed
244 if(strcmp(action,"Go")){ //Do not mark "Go" button
245 while ((obj = next())) {
246 if (obj->InheritsFrom(TGroupButton::Class())) {
247 button = (TGroupButton*)obj;
248 if(!strcmp(button->GetTitle(),action)){
249 if (button->GetBorderMode() > 0){
250 button->SetBorderMode(-1) ;
251 button->Modified(kTRUE);
258 //____________________________________________________________________________
259 void AliPHOSOnlineMonitor::SetInputFile(const char * filename){
260 //close previously opened
262 fInputFile = filename ;
264 //____________________________________________________________________________
265 void AliPHOSOnlineMonitor::DrawPedestals(){
266 //Prepare canvas and histograms for drawing pedestals
268 TIter nextCanvas(fCanvasList);
270 Bool_t exists = kFALSE ;
271 while((c=(TCanvas*)nextCanvas())){
272 if(!strcmp(c->GetName(),"Pedestals")){
278 c = new TDialogCanvas("Pedestals","Pedestals",300,200) ;
279 fCanvasList->AddLast(c) ;
282 TIter nextHisto(fHistosList);
285 while((h=(TH1D*)nextHisto())){
286 if(!strcmp(h->GetName(),"hPedestals")){
292 h = new TH1D("hPedestals","Pedestals per event",fGeom->GetNModules()*fGeom->GetNCristalsInModule(),0.,
293 1.*fGeom->GetNModules()*fGeom->GetNCristalsInModule()) ;
294 fHistosList->AddLast(h) ;
300 fScanPed = kTRUE ; //We will scan pedestals
303 //____________________________________________________________________________
304 void AliPHOSOnlineMonitor::DrawTriggers(){
305 //Prepare canvas and histogram for drawing triggers
307 TIter nextCanvas(fCanvasList);
309 Bool_t exists = kFALSE ;
310 while((c=(TCanvas*)nextCanvas())){
311 if(!strcmp(c->GetName(),"Triggers")){
317 c = new TDialogCanvas("Triggers","Triggers",200,200) ;
318 fCanvasList->AddLast(c) ;
321 TIter nextHisto(fHistosList);
324 while((h=(TH1D*)nextHisto())){
325 if(!strcmp(h->GetName(),"hTriggers")){
331 h = new TH1D("hTriggers","Triggers",2,0.,2.) ;
332 fHistosList->AddLast(h) ;
335 h->SetBit(TH1::kCanRebin);
336 h->Fill("LED",0.0000001) ;
337 h->Fill("PUL",0.0000001) ;
338 h->Fill("PED",0.0000001) ;
339 h->Fill("NEL",0.0000001) ;
340 h->Fill("WEL",0.0000001) ;
341 h->Fill("SOB",0.0000001) ;
342 h->Fill("EOB",0.0000001) ;
343 h->Fill("wrong",0.0000001) ;
344 h->LabelsOption("h");
350 //____________________________________________________________________________
351 void AliPHOSOnlineMonitor::DrawSpectrum(const char * opt){
352 //Prepare canvas and histograms for drawing spectra of all reconstructed particles or photons
354 TString name("Spectrum") ;
357 TIter nextCanvas(fCanvasList);
359 Bool_t exists = kFALSE ;
360 while((c=(TCanvas*)nextCanvas())){
361 if(!strcmp(c->GetName(),name.Data())){
367 c = new TDialogCanvas(name,name,250,300) ;
368 fCanvasList->AddLast(c) ;
371 TIter nextHisto(fHistosList);
375 while((h=(TH1D*)nextHisto())){
376 if(!strcmp(h->GetName(),name.Data())){
382 h = new TH1D(name,name,100,0.,100.) ;
383 fHistosList->AddLast(h) ;
389 fReconstruct = kTRUE ;
390 fScanSig = kTRUE ; //We will scan pedestals
393 //____________________________________________________________________________
394 void AliPHOSOnlineMonitor::DrawMinv(){
395 TIter nextCanvas(fCanvasList);
397 Bool_t exists = kFALSE ;
398 while((c=(TCanvas*)nextCanvas())){
399 if(!strcmp(c->GetName(),"InvMass")){
405 c = new TDialogCanvas("InvMass","Invariant mass",300,200) ;
406 fCanvasList->AddLast(c) ;
409 TIter nextHisto(fHistosList);
412 while((h=(TH1D*)nextHisto())){
413 if(!strcmp(h->GetName(),"hInvMass")){
419 h = new TH1D("hInvMass","hInvMass",1000,0.,1.0) ;
420 fHistosList->AddLast(h) ;
426 fReconstruct = kTRUE ;
427 fScanSig = kTRUE ; //We will scan pedestals
429 //____________________________________________________________________________
430 void AliPHOSOnlineMonitor::DrawEdep(Int_t mod,const char * opt){
432 sprintf(name,"Edep%s %d",opt,mod) ;
434 TIter nextCanvas(fCanvasList);
436 Bool_t exists = kFALSE ;
437 while((c=(TCanvas*)nextCanvas())){
438 if(!strcmp(c->GetName(),name)){
444 c = new TDialogCanvas(name,name,300,200) ;
445 fCanvasList->AddLast(c) ;
448 TIter nextHisto(fHistosList);
451 sprintf(name,"hEdep%s%d",opt,mod) ;
452 while((h=(TH2D*)nextHisto())){
453 if(!strcmp(h->GetName(),name)){
459 h = new TH2D(name,name,fGeom->GetNPhi(),0.,1.*fGeom->GetNPhi(),fGeom->GetNZ(),0.,1.*fGeom->GetNZ()) ;
460 fHistosList->AddLast(h) ;
466 fScanSig = kTRUE ; //We will scan signal events
467 if(strstr(opt,"Cal"))
468 fReconstruct = kTRUE ;
471 //____________________________________________________________________________
472 void AliPHOSOnlineMonitor::ScanPedestals(TClonesArray * digits){
473 //This method is called for events with PED trigger
474 //We fill bins with ADC values
476 TH1D * h = (TH1D*)gROOT->FindObjectAny("hPedestals");
478 Error("ScanPedestals","Can not fild histogram hPedestals") ;
481 for(Int_t i=0; i<digits->GetEntriesFast(); i++){
482 AliPHOSDigit * dig = static_cast<AliPHOSDigit*>(digits->At(i)) ;
483 h->AddBinContent(dig->GetId(),dig->GetAmp()) ;
486 //____________________________________________________________________________
487 void AliPHOSOnlineMonitor::ScanEdep(TClonesArray * digits){
488 //Fill 2D distribution of ADC values in NEL and WEL events
489 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
490 AliPHOSCalibrationDB *cdb = 0 ;
492 cdb = gime->CalibrationDB() ;
497 for(Int_t i=0; i<digits->GetEntriesFast(); i++){
498 AliPHOSDigit * dig = static_cast<AliPHOSDigit*>(digits->At(i)) ;
500 fGeom->AbsToRelNumbering(dig->GetId(),relId) ;
501 if(mod != relId[0]){ //new module, look for histograms
503 sprintf(name,"hEdep%d",mod) ;
504 h = (TH2D*)gROOT->FindObjectAny(name);
505 sprintf(name,"hEdepCal%d",mod) ;
506 hCal = (TH2D*)gROOT->FindObjectAny(name);
509 h->Fill(relId[2]-0.1,relId[3]-0.1,1.*dig->GetAmp()) ;
511 hCal->Fill(relId[2]-0.1,relId[3]-0.1,cdb->Calibrate(dig->GetAmp(),dig->GetId())) ;
514 //____________________________________________________________________________
515 void AliPHOSOnlineMonitor::ScanRecon(TClonesArray * recParticles){
516 if(!recParticles || recParticles->GetEntries()==0) return ;
518 TH1D* hSpectr = (TH1D*)gROOT->FindObjectAny("hSpectrumall");
519 TH1D* hSpectrGam= (TH1D*)gROOT->FindObjectAny("hSpectrumgamma");
520 TH1D* hInvMass = (TH1D*)gROOT->FindObjectAny("hInvMass");
521 for(Int_t i=0; i<recParticles->GetEntriesFast() ; i++){
522 AliPHOSRecParticle * p = (AliPHOSRecParticle *)recParticles->At(i) ;
523 if(hSpectr)hSpectr->Fill(p->Energy()) ;
524 if(hSpectrGam && p->IsPhoton())hSpectrGam->Fill(p->Energy()) ;
526 for(Int_t j=i+1; j<recParticles->GetEntriesFast() ; j++){
527 AliPHOSRecParticle * p2 = (AliPHOSRecParticle *)recParticles->At(j) ;
528 Double_t e = p->Energy() + p2->Energy() ;
529 Double_t x = p->Px() + p2->Px() ;
530 Double_t y = p->Py() + p2->Py() ;
531 Double_t z = p->Pz() + p2->Pz() ;
532 Double_t m = e*e-x*x-y*y-z*z ;
533 hInvMass->Fill(m>0?TMath::Sqrt(m): 0. ) ;
538 //____________________________________________________________________________
539 void AliPHOSOnlineMonitor::ScanTrigger(Int_t trig){
540 //Fills trigger distribution
542 TH1D * h = (TH1D*)gROOT->FindObjectAny("hTriggers");
545 case AliPHOSRawStream::kLED : h->Fill("LED",1.) ; break ;
546 case AliPHOSRawStream::kPUL : h->Fill("PUL",1.) ; break ;
547 case AliPHOSRawStream::kPED : h->Fill("PED",1.) ; break ;
548 case AliPHOSRawStream::kNEL : h->Fill("NEL",1.) ; break ;
549 case AliPHOSRawStream::kWEL : h->Fill("WEL",1.) ; break ;
550 case AliPHOSRawStream::kSOB : h->Fill("SOB",1.) ; break ;
551 case AliPHOSRawStream::kEOB : h->Fill("EOB",1.) ; break ;
552 default : h->Fill("wrong",1.) ;
555 //____________________________________________________________________________
556 void AliPHOSOnlineMonitor::SetConTableDB(const char * filename){
557 //Read ConnectionTableDB from file
558 TFile * file = new TFile(filename) ;
559 AliPHOSConTableDB * tmp = (AliPHOSConTableDB*)file->Get("AliPHOSConTableDB") ;
560 fcdb = new AliPHOSConTableDB(*tmp) ;
563 //____________________________________________________________________________
564 void AliPHOSOnlineMonitor::Go(){
565 //Perform scan of curent event
566 gBenchmark->Start("PHOSOnlineMon");
568 //First test if we need "Connection table" then open it
572 Info("Go","Read Connection table from file \"ConTableDB.root\"") ;
574 Error("Go","Please, set connection table with SetConTableDB() method") ;
579 AliPHOSGetterLight * gime = AliPHOSGetterLight::Instance("PHOS","On Flight") ;
581 //Configure CalibrManager to read data from file
582 //Create calibration database and read it
583 AliPHOSCalibrationDB * calibDB = 0 ;
584 if(fScanSig || fReconstruct){ //We will ned calibration parameters
585 AliPHOSCalibrManager::GetInstance("CalibrDB.root","root") ;
586 //If we configured manager to read from ASCII file,
587 //give him connection table. OK, it will not harm in any case.
588 AliPHOSCalibrManager::GetInstance()->SetConTable(fcdb) ;
590 calibDB = new AliPHOSCalibrationDB("OnLine") ;
591 calibDB->GetParameters() ; //Read parameters using Manager
592 gime->SetCalibrationDB(calibDB) ;
596 AliRawReaderDateV3 *rawReader = new AliRawReaderDateV3(fInputFile) ;
597 rawReader->RequireHeader(kFALSE);
598 AliPHOSRawStream *rawStream = new AliPHOSRawStream(rawReader) ;
599 rawStream->SetConTableDB(fcdb) ;
601 TClonesArray * digits = gime->Digits() ;
602 TClonesArray * recParticles = gime->RecParticles() ;
603 AliPHOSClusterizerv1* clu = 0 ;
604 AliPHOSTrackSegmentMakerv1 * tsm = 0 ;
605 AliPHOSPIDv1 * pid = 0 ;
606 if(fReconstruct){ //We will need calibation parameters
607 clu = new AliPHOSClusterizerv1("PHOS","On Flight") ;
608 clu->SetWriting(0) ; //Do not write to file
609 clu->SetEmcMinE(0.05) ; //Minimal energy of the digit
610 clu->SetEmcLocalMaxCut(0.05) ; //Height of local maximum over environment
611 clu->SetEmcClusteringThreshold(0.2) ; //Minimal energy to start cluster
612 // clu->SetUnfolding(kFALSE) ; //Do not unfold
613 tsm = new AliPHOSTrackSegmentMakerv1("PHOS","On Flight") ;
614 tsm->SetWriting(0) ; //Do not write to file
615 pid = new AliPHOSPIDv1("PHOS","On Flight") ;
616 pid->SetWriting(0) ; //Do not write to file
620 //Scan all event in file
622 while(rawReader->NextEvent()){
623 //Is it PHYSICAL event
624 if(rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent){
626 if(fNevents%100 ==0){
627 printf("\b\b\b\b\b\b%6d",fNevents) ;
629 if(rawStream->ReadDigits(digits)){
633 ScanTrigger(rawStream->GetTrigger()) ;
634 if(rawStream->IsPEDevent() && fScanPed){
635 ScanPedestals(digits) ;
637 if((rawStream->IsNELevent() || rawStream->IsWELevent()) && fScanSig){
640 gime->Clusterizer()->Exec("") ;
641 gime->TrackSegmentMaker()->Exec("") ;
642 gime->PID()->Exec("") ;
643 ScanRecon(recParticles) ;
648 if(fNevents%fNUpdate == 0 ){ //upate all histograms
649 TIter nextCanvas(fCanvasList);
651 while((c=(TCanvas*)nextCanvas())){
656 gSystem->ProcessEvents();
658 // if(fNevents>=200)break ;
661 gBenchmark->Stop("PHOSOnlineMon");
662 Float_t time = gBenchmark->GetCpuTime("PHOSOnlineMon") ;
663 printf("took %f seconds for scanning, i.e. %f seconds per event %d \n",
664 time,time/fNevents,fNevents) ;
666 //Update canvas with histograms at the end
667 TIter nextCanvas(fCanvasList);
669 while((c=(TCanvas*)nextCanvas())){
676 printf("delete 1 \n") ;
677 if(calibDB) delete calibDB ;
681 //____________________________________________________________________________
682 void AliPHOSOnlineMonitor::Clean(){
683 //Cleans content of all histograms
685 TIter nextHisto(fHistosList);
687 while((h=(TH1D*)nextHisto())){
690 TIter nextCanvas(fCanvasList);
692 while((c=(TCanvas*)nextCanvas())){
696 //____________________________________________________________________________
697 void AliPHOSOnlineMonitor::Reset(){
698 //delets all canvas and histograms,
699 //marks buttons as unpressed
701 TIter nextHisto(fHistosList);
703 while((h=(TH1D*)nextHisto())){
704 fHistosList->Remove(h) ;
707 TIter nextCanvas(fCanvasList);
709 while((c=(TCanvas*)nextCanvas())){
710 fCanvasList->Remove(c) ;
714 TGroupButton *button;
715 TIter next(fPrimitives);
717 //Mark buttons as anpressed
718 while ((obj = next())) {
719 if (obj->InheritsFrom(TGroupButton::Class())) {
720 button = (TGroupButton*)obj;
721 if (button->GetBorderMode() < 0){
722 button->SetBorderMode(1) ;
723 button->Modified(kTRUE);
729 //____________________________________________________________________________
730 void AliPHOSOnlineMonitor::WriteHistograms(const char * filename){
731 //Write filled histograms to file
732 TFile * file = new TFile(filename,"Update") ;
734 TIter nextHisto(fHistosList);
736 while((h=(TH1*)nextHisto())){
737 h->Write(0,TObject::kOverwrite) ;