artdaq  v3_09_00
EventDump_module.cc
1 // Class: EventDump
3 // Module Type: analyzer
4 // File: EventDump_module.cc
5 // Description: Prints out information about each event.
7 
8 #include "art/Framework/Core/EDAnalyzer.h"
9 #include "art/Framework/Core/ModuleMacros.h"
10 #include "art/Framework/Principal/Event.h"
11 #include "art/Framework/Principal/Handle.h"
12 #include "canvas/Utilities/Exception.h"
13 
14 #include "artdaq-core/Data/ContainerFragment.hh"
15 #include "artdaq-core/Data/Fragment.hh"
16 #include "artdaq-core/Data/RawEvent.hh"
17 
18 #include <algorithm>
19 #include <cassert>
20 #include <cmath>
21 #include <fstream>
22 #include <iomanip>
23 #include <iostream>
24 #include <vector>
25 
26 namespace artdaq {
27 class EventDump;
28 }
29 
33 class artdaq::EventDump : public art::EDAnalyzer
34 {
35 public:
46  explicit EventDump(fhicl::ParameterSet const& pset);
47 
51  ~EventDump() override = default;
52 
60  void analyze(art::Event const& e) override;
61 
62 private:
63  EventDump(EventDump const&) = delete;
64  EventDump(EventDump&&) = delete;
65  EventDump& operator=(EventDump const&) = delete;
66  EventDump& operator=(EventDump&&) = delete;
67 
68  std::string raw_data_label_;
69  int verbosity_;
70 };
71 
72 artdaq::EventDump::EventDump(fhicl::ParameterSet const& pset)
73  : EDAnalyzer(pset)
74  , raw_data_label_(pset.get<std::string>("raw_data_label", "daq"))
75  , verbosity_(pset.get<int>("verbosity", 0)) {}
76 
77 void artdaq::EventDump::analyze(art::Event const& e)
78 {
79  if (verbosity_ > 0)
80  {
81  std::cout << "***** Start of EventDump for event " << e.event() << " *****" << std::endl;
82 
83  art::Handle<detail::RawEventHeader> header_handle;
84  e.getByLabel(raw_data_label_, "RawEventHeader", header_handle);
85 
86  if (header_handle.isValid())
87  {
88  std::ostringstream ostr;
89  header_handle->print(ostr);
90  std::cout << "Event Header: " << ostr.str() << std::endl;
91  }
92  else
93  {
94  std::cout << "Unable to read RawEventHeader for event " << e.event() << std::endl;
95  }
96 
97  std::vector<art::Handle<std::vector<artdaq::Fragment> > > fragmentHandles;
98  e.getManyByType(fragmentHandles);
99 
100  for (auto const& handle : fragmentHandles)
101  {
102  if (!handle->empty())
103  {
104  std::string instance_name = handle.provenance()->productInstanceName();
105  std::cout << instance_name << " fragments: " << std::endl;
106 
107  int jdx = 1;
108  for (auto const& frag : *handle)
109  {
110  std::cout << " " << jdx << ") fragment ID " << frag.fragmentID() << " has type "
111  << static_cast<int>(frag.type()) << ", timestamp " << frag.timestamp()
112  << ", has metadata " << std::boolalpha << frag.hasMetadata()
113  << ", and sizeBytes " << frag.sizeBytes()
114  << " (hdr=" << frag.headerSizeBytes()
115  << ", data=" << frag.dataSizeBytes()
116  << ", meta (calculated)=" << (frag.sizeBytes() - frag.headerSizeBytes() - frag.dataSizeBytes())
117  << ")";
118 
119  if (instance_name.compare(0, 9, "Container") == 0)
120  {
121  artdaq::ContainerFragment cf(frag);
122  std::cout << " (contents: type = " << static_cast<int>(cf.fragment_type()) << ", count = "
123  << cf.block_count() << ", missing data = " << cf.missing_data()
124  << ")" << std::endl;
125  ;
126  if (verbosity_ > 1)
127  {
128  for (size_t idx = 0; idx < cf.block_count(); ++idx)
129  {
130  auto thisFrag = cf.at(idx);
131  std::cout << " " << (idx + 1) << ") fragment type " << static_cast<int>(thisFrag->type())
132  << ", timestamp " << thisFrag->timestamp()
133  << ", has metadata " << std::boolalpha << thisFrag->hasMetadata()
134  << ", and sizeBytes " << thisFrag->sizeBytes()
135  << " (hdr=" << thisFrag->headerSizeBytes()
136  << ", data=" << thisFrag->dataSizeBytes()
137  << ", meta (calculated)=" << (thisFrag->sizeBytes() - thisFrag->headerSizeBytes() - thisFrag->dataSizeBytes())
138  << ")" << std::endl;
139  }
140  }
141  }
142  else
143  {
144  std::cout << std::endl;
145  }
146  ++jdx;
147  }
148  }
149  }
150 
151  std::cout << "***** End of EventDump for event " << e.event() << " *****" << std::endl;
152  }
153 }
154 
155 DEFINE_ART_MODULE(artdaq::EventDump)// NOLINT(performance-unnecessary-value-param)
void analyze(art::Event const &e) override
This method is called for each art::Event in a file or run.
~EventDump() override=default
Default virtual Destructor.
EventDump(fhicl::ParameterSet const &pset)
EventDump Constructor.
Write Event information to the console.