Belle II Software development
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
22using namespace std;
23using namespace Belle2;
24
25//-----------------------------------------------------------------
26// Register the Module
27//-----------------------------------------------------------------
28REG_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 addParam("CanvasSendDefault", m_canvasSendDefault, "Send untagged canvases", false);
41 B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: Constructor done.");
42}
43
44
46
48{
49 if (m_sock != nullptr) delete m_sock;
50 m_sock = new TSocket(m_hostname.c_str(), m_port);
51 B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: initialized.");
52}
53
54
56{
57 B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: beginRun called.");
58}
59
60
62{
63 B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: event called.");
64 TMessage mess(kMESS_OBJECT);
65
66 TSeqCollection* seq = gROOT->GetListOfCanvases();
67 TIter nextkey(seq);
68 TObject* obj = 0;
69
70 time_t now = time(0);
71 char mbstr[100];
72 strftime(mbstr, sizeof(mbstr), "%F %T", localtime(&now));
73
74 auto& clist = getCanvasUpdatedList();
75 int sent_canvases = 0;
76
77 B2INFO("[" << mbstr << "] before sending " << seq->GetEntries() << " objects.");
78 bool first_try = true;
79 while ((obj = (TObject*)nextkey())) {
80 if (obj->IsA()->InheritsFrom("TCanvas")) {
81 TCanvas* c = (TCanvas*) obj;
82 auto process_canvas = m_canvasSendDefault;
83
84 auto it = clist.find(c->GetName());
85 if (it != clist.end()) {
86 process_canvas = it->second;
87 }
88 if (!process_canvas) continue;
89 sent_canvases++;
90
91 mess.Reset();
92 mess.WriteObject(c); // write object in message buffer
93 if (m_sock->Send(mess) < 0) {
94 if (!first_try) {
95 break;//Only try to reconnect once per event
96 } else {
97 first_try = false;
98 }
99 //The plain TSocket can't reconnect, so delete and create a new one
100 delete m_sock;
101 m_sock = new TSocket(m_hostname.c_str(), m_port);
102 //Try to send the failed message again
103 if (m_sock->Send(mess) < 0) {
104 //If failed a second time stop for this event
105 break;
106 }
107 }
108 }
109 }
110 now = time(0);
111 strftime(mbstr, sizeof(mbstr), "%F %T", localtime(&now));
112 B2INFO("[" << mbstr << "] after sending " << sent_canvases << " of " << seq->GetEntries() << " objects.");
113}
114
116{
117 B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: endRun called");
118}
119
120
122{
123 B2DEBUG(20, "DQMHistAnalysisOutputRelayMsg: terminate called");
124 delete m_sock;
125}
126
The base class for the histogram analysis module.
static const CanvasUpdatedList & getCanvasUpdatedList()
Get the list of the canvas update status.
TSocket * m_sock
The socket to the canvas server.
void terminate() override final
This method is called at the end of the event processing.
void event() override final
This method is called for each event.
bool m_canvasSendDefault
Send untagged canvas by default.
void endRun() override final
This method is called if the current run ends.
void beginRun() override final
Called when entering a new run.
std::string m_hostname
The hostname of the canvas server.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.
STL namespace.