FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Ex2Sampler.h
Go to the documentation of this file.
1 #ifndef EX2SAMPLER_H
2 #define EX2SAMPLER_H
3 
4 #include "MyDigi.h"
5 #include "RootSerializer.h"
6 #include "SerializerExample2.h"
7 
8 #include <FairMQDevice.h>
9 #include <Rtypes.h>
10 #include <TFile.h>
11 #include <TTree.h>
12 #include <chrono>
13 
14 class Ex2Sampler : public FairMQDevice
15 {
16  public:
18  : FairMQDevice()
19  , fInput(nullptr)
20  , fTree(nullptr)
21  , fFileName()
22  , fInputFile(nullptr)
23  {}
24 
25  Ex2Sampler(const Ex2Sampler&);
27 
28  virtual ~Ex2Sampler()
29  {
30  if (fInputFile) {
31  fInputFile->Close();
32  delete fInputFile;
33  }
34  }
35 
36  protected:
37  virtual void Init()
38  {
39  fFileName = fConfig->GetValue<std::string>("input-file");
40  fInputFile = TFile::Open(fFileName.c_str(), "READ");
41  if (fInputFile) {
42  fTree = static_cast<TTree*>(fInputFile->Get("cbmsim"));
43  if (fTree) {
44  fTree->SetBranchAddress("MyDigi", &fInput);
45  } else {
46  LOG(error) << "Could not find tree 'MyDigi'";
47  }
48  } else {
49  LOG(error) << "Could not open file " << fFileName << " in SimpleTreeReader::InitSource()";
50  }
51  }
52 
53  virtual void Run()
54  {
55  int64_t sentMsgs = 0;
56  const int64_t numEvents = fTree->GetEntries();
57  LOG(info) << "Number of events to process: " << numEvents;
58 
59  auto tStart = std::chrono::high_resolution_clock::now();
60 
61  for (int64_t idx = 0; idx < numEvents; idx++) {
62  fTree->GetEntry(idx);
63  Ex2Header* header = new Ex2Header();
64  header->EventNumber = idx;
65 
66  FairMQMessagePtr msgHeader(NewMessage());
67  Serialize<SerializerEx2>(*msgHeader, header);
68 
69  FairMQMessagePtr msg(NewMessage());
70  Serialize<RootSerializer>(*msg, fInput);
71 
72  FairMQParts parts;
73  parts.AddPart(std::move(msgHeader));
74  parts.AddPart(std::move(msg));
75 
76  if (Send(parts, "data1") > 0) {
77  sentMsgs++;
78  }
79 
80  if (NewStatePending()) {
81  break;
82  }
83  }
84 
85  auto tEnd = std::chrono::high_resolution_clock::now();
86  LOG(info) << "Sent " << sentMsgs
87  << " messages in: " << std::chrono::duration<double, std::milli>(tEnd - tStart).count() << " ms";
88  }
89 
90  private:
91  TClonesArray* fInput;
92  TTree* fTree;
93  std::string fFileName;
94  TFile* fInputFile;
95 };
96 
97 #endif
Ex2Sampler & operator=(const Ex2Sampler &)
virtual void Init()
Definition: Ex2Sampler.h:37
virtual void Run()
Definition: Ex2Sampler.h:53
virtual ~Ex2Sampler()
Definition: Ex2Sampler.h:28