Belle II Software  release-05-01-25
DQMHistAnalysisOutputRelayMsg.cc
1 //+
2 // File : DQMHistAnalysisOutputRelayMsg.cc
3 // Description : DQM Output, send Canvases to jsroot server.
4 //
5 // Author : B. Spruck
6 // Date : 25 - Mar - 2017
7 // based on work from Tomoyuki Konno, Tokyo Metropolitan Univerisity
8 //-
9 
10 
11 #include <dqm/analysis/modules/DQMHistAnalysisOutputRelayMsg.h>
12 #include <TROOT.h>
13 #include <TClass.h>
14 #include <TObject.h>
15 #include <TCanvas.h>
16 #include <TMessage.h>
17 #include <ctime>
18 
19 using namespace std;
20 using namespace Belle2;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(DQMHistAnalysisOutputRelayMsg)
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
33 {
34  //Parameter definition
35  addParam("Hostname", m_hostname, "Hostname of THTTP", std::string("localhost"));
36  addParam("Port", m_port, "Port number to THTTP", 9191);
37  B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: Constructor done.");
38 }
39 
40 
41 DQMHistAnalysisOutputRelayMsgModule::~DQMHistAnalysisOutputRelayMsgModule() { }
42 
43 void DQMHistAnalysisOutputRelayMsgModule::initialize()
44 {
45  if (m_sock != nullptr) delete m_sock;
46  m_sock = new TSocket(m_hostname.c_str(), m_port);
47  B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: initialized.");
48 }
49 
50 
51 void DQMHistAnalysisOutputRelayMsgModule::beginRun()
52 {
53  B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: beginRun called.");
54 }
55 
56 
57 void DQMHistAnalysisOutputRelayMsgModule::event()
58 {
59  B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: event called.");
60  TMessage mess(kMESS_OBJECT);
61 
62  TSeqCollection* seq = gROOT->GetListOfCanvases();
63  TIter nextkey(seq);
64  TObject* obj = 0;
65 
66  time_t now = time(0);
67  char mbstr[100];
68  strftime(mbstr, sizeof(mbstr), "%c", localtime(&now));
69 
70  B2INFO("[" << mbstr << "] before sending " << seq->GetEntries() << " objects.");
71  bool first_try = true;
72  while ((obj = (TObject*)nextkey())) {
73  if (obj->IsA()->InheritsFrom("TCanvas")) {
74  TCanvas* c = (TCanvas*) obj;
75  mess.Reset();
76  mess.WriteObject(c); // write object in message buffer
77  if (m_sock->Send(mess) < 0) {
78  if (!first_try) {
79  break;//Only try to reconnect once per event
80  } else {
81  first_try = false;
82  }
83  //The plain TSocket can't reconnect, so delete and create a new one
84  delete m_sock;
85  m_sock = new TSocket(m_hostname.c_str(), m_port);
86  //Try to send the failed message again
87  if (m_sock->Send(mess) < 0) {
88  //If failed a second time stop for this event
89  break;
90  }
91  }
92  }
93  }
94  now = time(0);
95  strftime(mbstr, sizeof(mbstr), "%c", localtime(&now));
96  B2INFO("[" << mbstr << "] after sending " << seq->GetEntries() << " objects.");
97 }
98 
99 void DQMHistAnalysisOutputRelayMsgModule::endRun()
100 {
101  B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: endRun called");
102 }
103 
104 
105 void DQMHistAnalysisOutputRelayMsgModule::terminate()
106 {
107  B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: terminate called");
108  delete m_sock;
109 }
110 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisOutputRelayMsgModule
Class definition for the output module of Sequential ROOT I/O.
Definition: DQMHistAnalysisOutputRelayMsg.h:23
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27