]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/amoreTPC-QA/src/ui/UIQA.cxx
419af12e29724dee2040e0cc2c87c00ebaf03e44
[u/mrichter/AliRoot.git] / TPC / amoreTPC-QA / src / ui / UIQA.cxx
1 /***************************************************************************
2  *   Copyright (C) 2007 by Filimon Roukoutakis                             *
3  *   Filimon.Roukoutakis@cern.ch                                           *
4  *                                                                         *
5  *   This program is free software; you can redistribute it and/or modify  *
6  *   it under the terms of the GNU General Public License as published by  *
7  *   the Free Software Foundation; either version 2 of the License, or     *
8  *   (at your option) any later version.                                   *
9  *                                                                         *
10  *   This program is distributed in the hope that it will be useful,       *
11  *   but WITHOUT ANY WARRANTY; without even the implied warranty of        *
12  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the         *
13  *   GNU General Public License for more details.                          *
14  *                                                                         *
15  *   You should have received a copy of the GNU General Public License     *
16  *   along with this program; if not, write to the                         *
17  *   Free Software Foundation, Inc.,                                       *
18  *   59 Temple Place - Suite 330, Boston, MA  02111-1307, USA.             *
19  ***************************************************************************/
20 #include "UIQA.h"
21 #include <AmoreDA.h>
22 #include <TROOT.h>
23 #include <TStyle.h>
24 #include <TList.h>
25 #include <TObjArray.h>
26 #include <TObjString.h>
27 #include <TString.h>
28 #include <TCollection.h>
29 #include <TGButton.h>
30 #include "AliTPCCalPad.h"
31 #include "AliTPCdataQA.h"
32 #include "AliTPCCalibCE.h"
33 #include "AliTPCCalibPulser.h"
34 #include "AliTPCCalibPedestal.h"
35 #include "AliTPCCalibViewer.h"
36 #include "AliTPCCalibViewerGUI.h"
37 #include "AliTPCPreprocessorOnline.h"
38 #include "AliCDBEntry.h"
39 ClassImp(amore::TPC::ui::UIQA)
40
41 #include <iostream>
42 #include <sstream>
43 #include <TCanvas.h>
44
45 namespace amore {
46
47 namespace TPC {
48
49 namespace ui {
50
51 using amore::subscriber::Subscribe;
52
53
54
55
56 UIQA::UIQA()  :
57   fMapCalibObjects(0x0),
58   fListCalibObjInfo(0x0), 
59   fListGuiObjects(0x0),
60   fAmoreDA(new amore::da::AmoreDA(amore::da::AmoreDA::kReceiver))
61   {
62
63  Construct(); // Temporary but important!!! Do not forget to put this call in the constructor for the time being!
64  fCycle=0; 
65 }
66
67
68 UIQA::~UIQA()
69 {
70   if ( fAmoreDA )       delete fAmoreDA;
71   if ( fMapCalibObjects )delete fMapCalibObjects;
72   if ( fListCalibObjInfo ) delete fListCalibObjInfo;
73   if ( fListGuiObjects ) delete fListGuiObjects;
74   printf("DA destructor called\n");
75 }
76
77 void UIQA::Construct() { // The custom GUI is constructed here. gRootFrame is the container of the custom widgets.
78   
79  fTab=new TGTab(amore::ui::gRootFrame);
80  amore::ui::gRootFrame->AddFrame(fTab);
81  //
82  //
83  TGCompositeFrame* tabCont1 =fTab->AddTab("Expert");
84  fExpert = tabCont1;
85  
86  fViewerGUI = new AliTPCCalibViewerGUI(tabCont1, 1000, 600, 0);
87  tabCont1->AddFrame(fViewerGUI , new TGLayoutHints(kLHintsExpandX | kLHintsExpandY, 0, 0, 0, 0));
88  //
89   TGCompositeFrame* tempFrame=fTab->AddTab("OverThreshold");
90  fEC[1]=new TRootEmbeddedCanvas("fEC0", tempFrame, 1000, 650);
91  tempFrame->AddFrame(fEC[1]);
92  fEC[1]->GetCanvas()->Divide(3,3);
93  //
94  //
95  tempFrame=fTab->AddTab("Charge");
96  fEC[2]=new TRootEmbeddedCanvas("fEC1", tempFrame, 1000, 650);
97  tempFrame->AddFrame(fEC[2]);
98  fEC[2]->GetCanvas()->Divide(2,3);
99  //
100  //
101  TGCompositeFrame* tabCalib =fTab->AddTab("DA Calib.");
102  SetupTabDACalib(tabCalib);
103  //
104  //
105  amore::ui::gRootFrame->MapSubwindows();
106  amore::ui::gRootFrame->Resize();
107  amore::ui::gRootFrame->MapWindow();
108
109  gROOT->SetStyle("Plain");
110  gStyle->SetFillColor(10);
111  gStyle->SetPadColor(10);
112  gStyle->SetCanvasColor(10);
113  gStyle->SetStatColor(10);
114  
115  gStyle->SetPalette(1,0);
116  gStyle->SetNumberContours(30);
117  gStyle->SetOptFit(111);
118  
119  gStyle->SetCanvasBorderMode(-1);
120  gStyle->SetCanvasBorderSize(1);
121  gStyle->SetCanvasColor(10);
122  
123  gStyle->SetFrameFillColor(10);
124  gStyle->SetFrameBorderSize(1);
125  gStyle->SetFrameBorderMode(-1);
126  gStyle->SetFrameLineWidth(1);
127  SubscribeMonitorObjects();
128  //Perhaps move later to a button action
129  RetrieveFromAmoreDB();
130 }
131
132 void UIQA::SetupTabDACalib(TGCompositeFrame *frame){
133   if ( !fListGuiObjects ){
134     fListGuiObjects = new TList;
135     fListGuiObjects->SetOwner();
136   }
137   // Update Button
138   TGTextButton *btnUpdate = new TGTextButton(frame, "&Update");
139   frame->AddFrame(btnUpdate, new TGLayoutHints(kLHintsExpandX, 10, 10, 0, 0));
140    //fBtnDraw->Connect("Clicked()", "AliTPCCalibViewerGUI", this, "DoTest(=\"fBtnDraw clicked\")");
141   btnUpdate->Connect("Clicked()", "UIQA", this, "UpdateAmoreDBValues()");
142   btnUpdate->SetToolTipText("Update values from the .");
143   fListGuiObjects->Add(btnUpdate);
144   //Update info boxes 
145 /*  TGGroupFrame *ldcPedestalGroup = new TGGroupFrame(ftabLeft0, "Pedestal/Noise Info", kVerticalFrame | kFitWidth | kFitHeight);
146   frame->AddFrame(ldcPedestalGroup, new TGLayoutHints(kLHintsExpandX, 0, 0, 10, 0));
147   for (Int_t ildc=0;ildc<36;++ildc){
148     char side='A';
149     if (isec>17) side='C';
150     TString amoreDAname(Form("ldc-TPC-%c%02d-%s",side,isec,"pedestals"));
151     
152     
153   }*/
154 }
155
156 void UIQA::UpdateAmoreDBValues(){
157   RetrieveFromAmoreDB();
158 }
159
160 void UIQA::SubscribeMonitorObjects() { // Before using any MonitorObject, a subscription should be made.
161
162   //std::ostringstream stringStream;
163  //amore::core::String_t sourceName="CEDA/", subscription; // The agent name acting as a source could be concatenated with all the objects it contains
164  //subscription=sourceName+"CE";
165  //Subscribe(subscription.c_str()); // Here you put a series of subscriptions where the string corresponds to the object name as published in the Publisher Module. As these names are internal to the QA framework, the recommended way of having consistency between AMORE and QA is to factor-out of QA the function that represents the histogram naming convention as a separate AliRoot class/function and use it from inside QA and AMORE.
166  //...
167   printf("UIQA::SubscribeMonitorObjects\n");
168   amore::core::String_t sourceName="TPCQA/", subscription; 
169   subscription=sourceName+"TPCRAW";
170   Subscribe(subscription.c_str());
171   printf("%s\n",subscription.c_str());
172   //subscription=sourceName+"hist";
173   //Subscribe(subscription.c_str());
174   
175 }
176
177 void UIQA::Update() { 
178   // This is executed after getting the updated contents of the subscribed MonitorObjects. Notice that the output of moInt[i] and moString[i] varies with time for a specific i because on the dqmAgent the "quality" check fails or succeeds. This is the essence of automatic data quality checks in AMORE. Try to use the moString[i] on a text widget to alert the shifter, or -depending of the value of moInt[i], 0 or 1- make part of the screen change color...
179  std::ostringstream stringStream;
180  // Example of accessing a normal TObject. The name is the name of the object in the QA framework
181
182
183  static Int_t counter0 = 0, counter1=0;
184  printf("UIQA::Update0\t%d\t%d\n",counter0,counter1); 
185
186  amore::core::MonitorObjectTObject* ptr=0;
187  AliTPCdataQA *tpcqa=0; 
188  ptr=gSubscriber->At<amore::core::MOTObj>("TPCQA/TPCRAW");
189  tpcqa=0;
190  printf("Pointer - %p\n",ptr);
191  if(ptr) {
192    tpcqa=(AliTPCdataQA*)ptr->Object();
193    printf("Pointertpcqa - %p\n",tpcqa);
194    tpcqa->Print();
195  }
196  
197  printf("Update1\t%d\t%d\n",counter0,counter1); 
198  counter0++;
199  if (!tpcqa) {
200    MakeTree(0x0);
201    return;
202  }
203  printf("Update2\t%d",counter1);
204  counter1++;
205  //fExpert->SetTitle(name);
206  //
207  // Over threshold
208  //
209
210  TCanvas *canvas  = fEC[1]->GetCanvas();
211  if (tpcqa->GetNoThreshold()){
212    canvas->cd(1);
213    tpcqa->GetNoThreshold()->MakeHisto1D()->Draw();
214    canvas->cd(2);
215    tpcqa->GetNoThreshold()->MakeHisto2D(0)->Draw("colz");
216    canvas->cd(3);
217    tpcqa->GetNoThreshold()->MakeHisto2D(1)->Draw("colz");
218  }
219  //
220  if (tpcqa->GetNTimeBins()){
221    canvas->cd(4);
222    tpcqa->GetNTimeBins()->MakeHisto1D()->Draw();
223    canvas->cd(5);
224    tpcqa->GetNTimeBins()->MakeHisto2D(0)->Draw("colz");
225    canvas->cd(6);
226    tpcqa->GetNTimeBins()->MakeHisto2D(1)->Draw("colz");
227  }
228
229  if (tpcqa->GetNPads()){
230    canvas->cd(7);
231    tpcqa->GetNPads()->MakeHisto1D()->Draw();
232    canvas->cd(8);
233    tpcqa->GetNPads()->MakeHisto2D(0)->Draw("colz");
234    canvas->cd(9);
235    tpcqa->GetNPads()->MakeHisto2D(1)->Draw("colz");
236  }
237
238  //
239  // Mean charge
240  //
241  canvas  = fEC[2]->GetCanvas();
242  if (tpcqa->GetMeanCharge()){
243    canvas->cd(1);
244    tpcqa->GetMeanCharge()->MakeHisto1D()->Draw();
245    canvas->cd(3);
246    tpcqa->GetMeanCharge()->MakeHisto2D(0)->Draw("colz");
247    canvas->cd(5);
248    tpcqa->GetMeanCharge()->MakeHisto2D(1)->Draw("colz");
249    //
250    canvas->cd(2);
251    tpcqa->GetMaxCharge()->MakeHisto1D()->Draw();
252    canvas->cd(4);
253    tpcqa->GetMaxCharge()->MakeHisto2D(0)->Draw("colz");
254    canvas->cd(6);
255    tpcqa->GetMaxCharge()->MakeHisto2D(1)->Draw("colz");
256  }
257  MakeTree(tpcqa);
258
259  char name[1000];
260  sprintf(name,"Summary Event -  %d",tpcqa->GetEventCounter());
261  fExpert->SetName(name);
262  for (Int_t i=0;i<3;i++) printf("*********************************\n");
263  printf("\n\nSummary Event -  %d\n\n",tpcqa->GetEventCounter());
264  for (Int_t i=0;i<3;i++) printf("*********************************\n");
265  
266
267 }
268
269 void UIQA::Process() {
270
271 }
272
273 void UIQA::StartOfCycle() {
274 }
275
276 void UIQA::EndOfCycle() {
277
278 }
279
280
281 void UIQA::MakeTree(AliTPCdataQA* ped){
282   //
283   //
284   AliTPCPreprocessorOnline * preprocesor = new AliTPCPreprocessorOnline;
285
286   //==============================
287   //    QA Information
288   //==============================
289   if (ped) {
290     if (ped->GetNLocalMaxima()) 
291       preprocesor->AddComponent(new AliTPCCalPad(*(ped->GetNLocalMaxima())));
292     if (ped->GetMaxCharge()) 
293       preprocesor->AddComponent(new AliTPCCalPad(*(ped->GetMaxCharge())));  
294     if (ped->GetMeanCharge()) 
295       preprocesor->AddComponent(new AliTPCCalPad(*(ped->GetMeanCharge())));  
296     if (ped->GetNoThreshold()) 
297       preprocesor->AddComponent(new AliTPCCalPad(*(ped->GetNoThreshold())));
298     if (ped->GetNTimeBins()) 
299       preprocesor->AddComponent(new AliTPCCalPad(*(ped->GetNTimeBins())));
300     if (ped->GetNPads()) 
301       preprocesor->AddComponent(new AliTPCCalPad(*(ped->GetNPads())));
302     if (ped->GetTimePosition()) 
303       preprocesor->AddComponent(new AliTPCCalPad(*(ped->GetTimePosition())));   
304   }
305   //==============================
306   //    Calibration Information
307   //==============================
308   TIter next(fMapCalibObjects);
309   while ( TObject *o=next() ) preprocesor->AddComponent(o->Clone());
310
311   char fname[10000];
312   sprintf(fname,"QAtree%d.root",fCycle);
313   preprocesor->DumpToFile(fname);
314   fCycle++;
315     //init viewer
316   AliTPCCalibViewer *viewer = fViewerGUI->GetViewer();
317   AliTPCCalibViewer *nviewer = new  AliTPCCalibViewer(fname, "calPads");
318   fViewerGUI->Initialize(nviewer);
319   delete preprocesor;
320 }
321
322 void UIQA::RetrieveFromAmoreDB(){
323   if ( fMapCalibObjects ) {
324     fMapCalibObjects->Delete();
325     delete fMapCalibObjects;
326   }
327   if ( !fListCalibObjInfo ){
328     fListCalibObjInfo = new TList;
329     fListCalibObjInfo->SetOwner();
330   }
331   fMapCalibObjects = new TList;
332   fMapCalibObjects->SetOwner();
333   //define calibration objects, type and role
334   TString calNames("Pedestals;Noise;PulserT0;PulserQ;PulserRMS;CET0;CEQ;CERMS");
335   TString daTypes("pedestals;pedestals;pulser;pulser;pulser;CE;CE;CE");
336   TString daRoles("ldc;ldc;ldc;ldc;ldc;mon-DA-TPC-0;mon-DA-TPC-0;mon-DA-TPC-0");
337   TObjArray *calNameArr=calNames.Tokenize(";");
338   TObjArray *daTypeArr=daTypes.Tokenize(";");
339   TObjArray *daRoleArr=daRoles.Tokenize(";");
340   Int_t nobj=calNameArr->GetEntries();
341   
342   for (Int_t ient=0; ient<nobj; ++ient){
343     const TString &calName = ((TObjString*)calNameArr->At(ient))->GetString();  
344     const TString &daType  = ((TObjString*)daTypeArr->At(ient))->GetString();
345     const TString &daRole  = ((TObjString*)daRoleArr->At(ient))->GetString();
346     AliTPCCalPad *calPad  = new AliTPCCalPad(calName.Data(),calName.Data());
347     fMapCalibObjects->Add(calPad);
348     if ( daRole.Contains("ldc") ) CollectFromLDCs(calPad,calName,daType);
349     else if ( daRole.Contains("mon") ) CollectFromMon(calPad,calName,daType,daRole);
350   }
351   delete calNameArr;
352   delete daTypeArr;
353   delete daRoleArr;
354 }
355
356 void UIQA::CollectFromLDCs(AliTPCCalPad *calPad, const TString &calName, const TString &daType){
357   for (int isec=0; isec<36; isec++){
358     TObject *obj=0x0;
359     Int_t res=0;
360     char side='A';
361     if (isec>17) side='C';
362     TString amoreDAname(Form("ldc-TPC-%c%02d-%s",side,isec,daType.Data()));
363     TString objName(Form("%s/%s",amoreDAname.Data(),calName.Data()));   
364     res=fAmoreDA->Receive(objName.Data(),obj);
365     if ( !res&&obj ){
366       printf ("Found Pedestal Object in '%s'\n",objName.Data());
367       TObjArray *arr=(TObjArray*)obj;
368        //IROC
369       AliTPCCalROC *rocMerge=(AliTPCCalROC*)arr->At(isec);
370       if (rocMerge) calPad->SetCalROC(rocMerge,isec);
371         //OROC
372       rocMerge=(AliTPCCalROC*)arr->At(isec+36);
373       if (rocMerge) calPad->SetCalROC(rocMerge,isec+36);
374     }else{
375       printf("Could not receive '%s'\n",objName.Data());
376     }      
377     // Get info message
378     TString infoName(Form("%s/%s",amoreDAname.Data(),"Info"));   
379     res=fAmoreDA->Receive(infoName.Data(),obj);
380     if ( !res&&obj ){
381       TObjString *ostr=(TObjString*)obj;
382       printf("Calib. object DA info from '%s': '%s'\n",objName.Data(),ostr->GetName());
383       //save to info list
384       TObjString *infoStr = (TObjString*)fListCalibObjInfo->FindObject(amoreDAname);
385       if ( !infoStr ) {
386         infoStr = new TObjString;
387         fListCalibObjInfo->Add(infoStr);
388       }
389       infoStr->GetString()=ostr->GetString();
390     }
391     
392   }    
393 }
394
395 void UIQA::CollectFromMon(AliTPCCalPad *calPad, const TString &calName, const TString &daType, const TString &mon){
396   TObject *obj=0x0;
397   Int_t res=0;
398   TString amoreDAname(Form("%s-%s",mon.Data(),daType.Data()));
399   TString objName(Form("%s/%s",amoreDAname.Data(),calName.Data()));   
400   res=fAmoreDA->Receive(objName.Data(),obj);
401   if ( !res&&obj ){
402     printf ("Found Pedestal Object in '%s'\n",objName.Data());
403     TObjArray *arr=(TObjArray*)obj;
404     for (Int_t isec=0; isec<72; isec++){
405       AliTPCCalROC *rocMerge=(AliTPCCalROC*)arr->At(isec);
406       if (rocMerge) calPad->SetCalROC(rocMerge,isec);     
407     }
408   }else{
409     printf("Could not receive '%s'\n",objName.Data());
410   }  
411     // Get info message
412   TString infoName(Form("%s/%s",amoreDAname.Data(),"Info"));   
413   res=fAmoreDA->Receive(infoName.Data(),obj);
414   if ( !res&&obj ){
415     TObjString *ostr=(TObjString*)obj;
416     printf("Calib. object DA info from '%s': '%s'\n",objName.Data(),ostr->GetName());
417     TObjString *infoStr = (TObjString*)fListCalibObjInfo->FindObject(amoreDAname);
418     if ( !infoStr ) {
419       infoStr = new TObjString;
420       fListCalibObjInfo->Add(infoStr);
421     }
422     infoStr->GetString()=ostr->GetString();
423   }
424 }
425
426 };
427
428 };
429
430 };