Belle II Software  release-05-01-25
DQMHistAnalysisInputSrv.cc
1 //+
2 // File : DQMHistAnalysisInputSrv.cc
3 // Description :
4 //
5 // Author : B. Spruck
6 // Date : 25 - Mar - 2017
7 // based on work from Tomoyuki Konno, Tokyo Metropolitan Univerisity
8 //-
9 
10 #include <dqm/analysis/modules/DQMHistAnalysisInputSrv.h>
11 
12 #include <TKey.h>
13 #include <TSystem.h>
14 
15 using namespace std;
16 using namespace Belle2;
17 
18 //-----------------------------------------------------------------
19 // Register the Module
20 //-----------------------------------------------------------------
21 REG_MODULE(DQMHistAnalysisInputSrv)
22 
23 //-----------------------------------------------------------------
24 // Implementation
25 //-----------------------------------------------------------------
26 
29 {
30  //Parameter definition
31  addParam("HistMemoryPath", m_mempath, "Path to Input Hist memory", string(""));
32  addParam("HistMemorySize", m_memsize, "Size of Input Hist memory", 10000000);
33  addParam("RefreshInterval", m_interval, "Refresh interval of histograms in ms", 2000);
34  B2DEBUG(20, "DQMHistAnalysisInputSrv: Constructor done.");
35 }
36 
37 
38 DQMHistAnalysisInputSrvModule::~DQMHistAnalysisInputSrvModule() { }
39 
40 void DQMHistAnalysisInputSrvModule::initialize()
41 {
42  if (m_memory != nullptr) delete m_memory;
43  m_memory = new DqmMemFile(m_mempath.c_str());
44  m_eventMetaDataPtr.registerInDataStore();
45  //m_serv = new THttpServer("http:8081");
46  //m_serv->SetReadOnly(kFALSE);
47  B2DEBUG(20, "DQMHistAnalysisInputSrv: initialized.");
48 }
49 
50 
51 void DQMHistAnalysisInputSrvModule::beginRun()
52 {
53  B2DEBUG(20, "DQMHistAnalysisInputSrv: beginRun called.");
54 }
55 
56 void DQMHistAnalysisInputSrvModule::event()
57 {
58  std::vector<TH1*> h;
59  TMemFile* file = m_memory->LoadMemFile();
60  file->cd();
61  TIter next(file->GetListOfKeys());
62  TKey* key = NULL;
63  while ((key = (TKey*)next())) {
64  h.push_back((TH1*)key->ReadObj());
65  }
66  resetHist();
67  for (size_t i = 0; i < h.size(); i++) {
68  addHist("", h[i]->GetName(), h[i]);
69  B2DEBUG(2, "Found : " << h[i]->GetName() << " : " << h[i]->GetEntries());
70  std::string vname = h[i]->GetName();
71  setFloatValue(vname + ".entries", h[i]->GetEntries());
72  if (h[i]->GetDimension() == 1) {
73  setFloatValue(vname + ".rms", h[i]->GetRMS());
74  setFloatValue(vname + ".rmserr", h[i]->GetRMSError());
75  setFloatValue(vname + ".mean", h[i]->GetMean());
76  setFloatValue(vname + ".meanerr", h[i]->GetMeanError());
77  } else if (h[i]->GetDimension() == 2) {
78  setFloatValue(vname + ".xrms", h[i]->GetRMS(1));
79  setFloatValue(vname + ".xrmserr", h[i]->GetRMSError(1));
80  setFloatValue(vname + ".xmean", h[i]->GetMean(1));
81  setFloatValue(vname + ".xmeanerr", h[i]->GetMeanError(1));
82  setFloatValue(vname + ".yrms", h[i]->GetRMS(2));
83  setFloatValue(vname + ".yrmserr", h[i]->GetRMSError(2));
84  setFloatValue(vname + ".ymean", h[i]->GetMean(2));
85  setFloatValue(vname + ".ymeanerr", h[i]->GetMeanError(2));
86  }
87  }
88  m_count++;
89  m_eventMetaDataPtr.create();
90  m_eventMetaDataPtr->setExperiment(m_expno);
91  m_eventMetaDataPtr->setRun(m_runno);
92  m_eventMetaDataPtr->setEvent(m_count);
93 
94  TTimer t(m_interval, kFALSE);// in ms
95 
96  do { // call at least once!
97  //m_serv->ProcessRequests();
98  //gSystem->Sleep(10); // 10 ms sleep
99  } while (!t.CheckTimer(gSystem->Now()));
100 
101 }
102 
103 void DQMHistAnalysisInputSrvModule::endRun()
104 {
105  B2DEBUG(20, "DQMHistAnalysisInputSrv: endRun called");
106 }
107 
108 
109 void DQMHistAnalysisInputSrvModule::terminate()
110 {
111  B2DEBUG(20, "DQMHistAnalysisInputSrv: terminate called");
112 }
113 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::DQMHistAnalysisInputSrvModule
Class definition for the output module of Sequential ROOT I/O.
Definition: DQMHistAnalysisInputSrv.h:28
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DqmMemFile
Definition: DqmMemFile.h:28
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27