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