Belle II Software development
HistoRelay2.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#include "daq/dqm/HistoRelay2.h"
10#include <framework/pcore/EvtMessage.h>
11#include <TFile.h>
12#include <TH1.h>
13#include <TText.h>
14#include <TKey.h>
15
16using namespace Belle2;
17using namespace std;
18
19HistoRelay2::HistoRelay2(string& filename, string& dest, int port)
20{
21 m_dest = dest;
22 m_port = port;
23 m_filename = filename;
24 m_sock = new EvtSocketSend(m_dest, m_port);
25}
26
27HistoRelay2::~HistoRelay2()
28{
29 delete m_sock;
30}
31
32int HistoRelay2::collect()
33{
34 // printf ( "HistoRelay2 : collect!!\n" );
35 EvtMessage* msg = StreamFile(m_filename);
36
37 if (!msg) return 1; // error indicator
38
39 auto ret = m_sock->send(msg);
40
41 delete (msg);
42
43 if (ret < 0) {
44 // Socket error, e.g server died
45 delete m_sock;
46 // try reconnect
47 printf("HistoRelay2: socket seems dead -> reconnect\n");
48 m_sock = new EvtSocketSend(m_dest, m_port);
49 }
50
51 return 0;
52}
53
54// Copy Shared Memory to local and stream
55EvtMessage* HistoRelay2::StreamFile(string& filename)
56{
57 TFile* file = new TFile((std::string("/dev/shm/") + filename).c_str(), "read");
58 if (file == NULL || file->IsZombie()) return NULL;
59 file->cd();
60 MsgHandler hdl(0);
61 int numobjs = 0;
62 StreamHistograms(gDirectory, &hdl, numobjs);
63 file->Close();
64 delete file;
65 // printf ( "DqmMemFile::StreamMemFile : streamed %d histograms in EvtMessage\n", numobjs );
66 EvtMessage* msg = hdl.encode_msg(MSG_EVENT);
67 (msg->header())->reserved[0] = 0;
68 (msg->header())->reserved[1] = numobjs;
69 return msg;
70}
71
72int HistoRelay2::StreamHistograms(TDirectory* curdir, MsgHandler* msg, int& numobjs)
73{
74 TList* keylist = curdir->GetListOfKeys();
75
76 TIter nextkey(keylist);
77 TKey* key = 0;
78 while ((key = (TKey*)nextkey())) {
79 TObject* obj = curdir->FindObjectAny(key->GetName());
80 if (obj->IsA()->InheritsFrom("TH1")) {
81 TH1* h1 = (TH1*) obj;
82 // printf ( "Key = %s, entry = %f\n", key->GetName(), h1->GetEntries() );
83 msg->add(h1, h1->GetName());
84 numobjs++;
85 } else if (obj->IsA()->InheritsFrom(TDirectory::Class())) {
86 // printf ( "New directory found %s, Go into subdir\n", obj->GetName() );
87 TDirectory* tdir = (TDirectory*) obj;
88 // m_msg->add(tdir, tdir->GetName());
89 TText subdir(0, 0, tdir->GetName());
90 msg->add(&subdir, "SUBDIR:" + string(obj->GetName())) ;
91 numobjs++;
92 tdir->cd();
93 StreamHistograms(tdir, msg, numobjs);
94 TText command(0, 0, "COMMAND:EXIT");
95 msg->add(&command, "SUBDIR:EXIT");
96 numobjs++;
97 curdir->cd();
98 }
99 }
100 return 0;
101}
Class to manage streamed object.
Definition: EvtMessage.h:59
EvtHeader * header()
Get pointer to EvtHeader.
Definition: EvtMessage.cc:161
A class to encode/decode an EvtMessage.
Definition: MsgHandler.h:103
virtual void add(const TObject *, const std::string &name)
Add an object to be streamed.
Definition: MsgHandler.cc:46
Abstract base class for different kinds of events.
STL namespace.