artdaq_demo  v3_06_00
ToyDump_module.cc
1 // Class: ToyDump
3 // Module Type: analyzer
4 // File: ToyDump_module.cc
5 // Description: Prints out information about each event.
7 
8 #define TRACE_NAME "ToyDump"
9 
10 #include "art/Framework/Core/EDAnalyzer.h"
11 #include "art/Framework/Core/ModuleMacros.h"
12 #include "art/Framework/Principal/Event.h"
13 #include "art/Framework/Principal/Handle.h"
14 #include "canvas/Utilities/Exception.h"
15 
16 #include "artdaq-core-demo/Overlays/FragmentType.hh"
17 #include "artdaq-core-demo/Overlays/ToyFragment.hh"
18 #include "artdaq-core/Data/ContainerFragment.hh"
19 #include "artdaq-core/Data/Fragment.hh"
20 
21 #include <algorithm>
22 #include <cassert>
23 #include <cmath>
24 #include <fstream>
25 #include <iomanip>
26 #include <vector>
27 #include <iostream>
28 
29 namespace demo
30 {
31  class ToyDump;
32 }
33 
37 class demo::ToyDump : public art::EDAnalyzer
38 {
39 public:
54  explicit ToyDump(fhicl::ParameterSet const& pset);
55 
59  virtual ~ToyDump();
60 
65  virtual void analyze(art::Event const& evt);
66 
67 private:
68  std::string raw_data_label_;
69  int num_adcs_to_write_;
70  int num_adcs_to_print_;
71  bool binary_mode_;
72  uint32_t columns_to_display_on_screen_;
73  std::string output_file_name_;
74 };
75 
76 
77 demo::ToyDump::ToyDump(fhicl::ParameterSet const& pset)
78  : EDAnalyzer(pset)
79  , raw_data_label_(pset.get<std::string>("raw_data_label", "daq"))
80  , num_adcs_to_write_(pset.get<int>("num_adcs_to_write", 0))
81  , num_adcs_to_print_(pset.get<int>("num_adcs_to_print", 10))
82  , binary_mode_(pset.get<bool>("binary_mode", true))
83  , columns_to_display_on_screen_(pset.get<uint32_t>("columns_to_display_on_screen", 10))
84  , output_file_name_(pset.get<std::string>("output_file_name", "out.bin"))
85 {}
86 
88 
89 void demo::ToyDump::analyze(art::Event const& evt)
90 {
91  art::EventNumber_t eventNumber = evt.event();
92 
93  // ***********************
94  // *** Toy Fragments ***
95  // ***********************
96 
97  artdaq::Fragments fragments;
98  artdaq::FragmentPtrs containerFragments;
99  std::vector<std::string> fragment_type_labels{ "TOY1", "TOY2", "ContainerTOY1", "ContainerTOY2" };
100 
101  for (auto label : fragment_type_labels)
102  {
103  art::Handle<artdaq::Fragments> fragments_with_label;
104 
105  evt.getByLabel(raw_data_label_, label, fragments_with_label);
106  if (!fragments_with_label.isValid()) continue;
107 
108  if (label == "Container" || label == "ContainerTOY1" || label == "ContainerTOY2")
109  {
110  for (auto cont : *fragments_with_label)
111  {
112  artdaq::ContainerFragment contf(cont);
113  for (size_t ii = 0; ii < contf.block_count(); ++ii)
114  {
115  containerFragments.push_back(contf[ii]);
116  fragments.push_back(*containerFragments.back());
117  }
118  }
119  }
120  else
121  {
122  for (auto frag : *fragments_with_label)
123  {
124  fragments.emplace_back(frag);
125  }
126  }
127  }
128 
129  // look for raw Toy data
130  TLOG(TLVL_INFO) << "Run " << evt.run() << ", subrun " << evt.subRun()
131  << ", event " << eventNumber << " has " << fragments.size()
132  << " fragment(s) of type TOY1 or TOY2";
133 
134  for (const auto& frag : fragments)
135  {
136  ToyFragment bb(frag);
137 
138  TLOG(TLVL_INFO) << fragmentTypeToString(static_cast<demo::detail::FragmentType>(frag.type())) << " fragment " << frag.fragmentID()
139  << " w/ seqID " << frag.sequenceID() << " and timestamp " << frag.timestamp()
140  << " has total ADC counts = " << bb.total_adc_values() << ", trig # = " << bb.hdr_trigger_number() << ", dist_type = " << bb.hdr_distribution_type();
141 
142  if (frag.hasMetadata())
143  {
144  ToyFragment::Metadata const* md =
145  frag.metadata<ToyFragment::Metadata>();
146  TLOG(TLVL_DEBUG) << "Fragment metadata: " << std::showbase << "Board serial number = "
147  << md->board_serial_number << ", sample bits = "
148  << md->num_adc_bits
149  << " -> max ADC value = "
150  << bb.adc_range((int)md->num_adc_bits);
151  }
152 
153  if (num_adcs_to_write_ >= 0)
154  {
155  uint32_t numAdcs = num_adcs_to_write_;
156  if (num_adcs_to_write_ == 0) numAdcs = bb.total_adc_values();
157  else if (static_cast<uint32_t>(num_adcs_to_write_) > bb.total_adc_values())
158  {
159  TLOG(TLVL_WARNING) << "Asked for more ADC values to file than are in Fragment. Only writing what's here...";
160  numAdcs = bb.total_adc_values();
161  }
162  if (binary_mode_) {
163  std::ofstream output(output_file_name_, std::ios::out | std::ios::app | std::ios::binary);
164  for (uint32_t i_adc = 0; i_adc < numAdcs; ++i_adc)
165  {
166  output.write((char*)(bb.dataBeginADCs() + i_adc), sizeof(ToyFragment::adc_t));
167  }
168  output.close();
169  }
170  else {
171  std::ofstream output(output_file_name_, std::ios::out | std::ios::app);
172  output << fragmentTypeToString(static_cast<demo::detail::FragmentType>(frag.type())) << "\t" << frag.fragmentID();
173 
174  for (uint32_t i_adc = 0; i_adc < numAdcs; ++i_adc)
175  {
176  output << "\t" << std::to_string(*(bb.dataBeginADCs() + i_adc));
177  }
178  output << std::endl;
179  output.close();
180 
181  }
182  }
183 
184  if (num_adcs_to_print_ >= 0)
185  {
186  uint32_t numAdcs = num_adcs_to_print_;
187  if (num_adcs_to_print_ == 0) numAdcs = bb.total_adc_values();
188  else if (static_cast<uint32_t>(num_adcs_to_print_) > bb.total_adc_values())
189  {
190  TLOG(TLVL_WARNING) << "Asked for more ADC values to file than are in Fragment. Only writing what's here...";
191  numAdcs = bb.total_adc_values();
192  }
193 
194  TLOG(TLVL_INFO) << "First " << numAdcs << " ADC values in the fragment:";
195  int rows = 1 + (int)((num_adcs_to_print_ - 1) / columns_to_display_on_screen_);
196  uint32_t adc_counter = 0;
197  for (int idx = 0; idx < rows; ++idx)
198  {
199  std::ostringstream o;
200  o << std::right;
201  o << std::setw(4) << std::setfill('.');
202  o << (idx * columns_to_display_on_screen_) << ": ";
203  for (uint32_t jdx = 0; jdx < columns_to_display_on_screen_; ++jdx)
204  {
205  if (adc_counter >= numAdcs) { break; }
206  o << std::setw(6) << std::setfill(' ');
207  o << bb.adc_value(adc_counter);
208  ++adc_counter;
209  }
210  TLOG(TLVL_INFO) << o.str();
211  }
212  }
213  }
214 }
215 
216 DEFINE_ART_MODULE(demo::ToyDump)
std::vector< Fragment > Fragments
Header::trigger_number_t hdr_trigger_number() const
std::string fragmentTypeToString(FragmentType val)
adc_t adc_value(uint32_t index)
ToyDump(fhicl::ParameterSet const &pset)
ToyDump Constructor.
uint8_t hdr_distribution_type() const
An art::EDAnalyzer module designed to display the data from demo::ToyFragment objects.
adc_t const * dataBeginADCs() const
static size_t adc_range(int daq_adc_bits)
Metadata::count_t block_count() const
std::list< FragmentPtr > FragmentPtrs
virtual ~ToyDump()
ToyDump Destructor.
virtual void analyze(art::Event const &evt)
Analyze an event. Called by art for each event in run (based on command line options) ...
size_t total_adc_values() const