tdaq-develop-2025-02-12
CrvDQM_module.cc
1 // ROOT-based DQM and viewer for the CRV
2 // Authors: Sam Grant and Simon Corrodi
3 // Data: Feb 2025
4 
5 // C++ includes
6 #include <thread>
7 
8 // art includes
9 #include "art/Framework/Core/EDAnalyzer.h"
10 #include "art/Framework/Core/ModuleMacros.h"
11 #include "art/Framework/Principal/Event.h"
12 #include "art/Framework/Principal/Handle.h"
13 #include "art/Framework/Principal/Run.h"
14 
15 // artdaq includes
16 #include "artdaq-core/Data/ContainerFragment.hh"
17 #include "artdaq-core/Data/Fragment.hh"
18 #include "artdaq-core-mu2e/Data/CRVDataDecoder.hh"
19 #include "artdaq-core-mu2e/Overlays/DTCEventFragment.hh"
20 #include "artdaq-core-mu2e/Overlays/FragmentType.hh"
21 
22 // ROOT includes
23 #include "TCanvas.h"
24 #include "TH1D.h"
25 #include "THttpServer.h"
26 #include "TSystem.h"
27 
28 // Custom includes
29 #include "CrvDQMStyle.hh"
30 
31 namespace demo { // what is the appropriate namespace?
32 
33 class CrvDQM : public art::EDAnalyzer {
34 public:
35  // Constructor
36  explicit CrvDQM(fhicl::ParameterSet const& ps);
37  // Destructor
38  ~CrvDQM() override;
39 
40 private:
41 
42  // Functions
43  void beginJob() override;
44  void analyze(art::Event const& e) override;
45  void endJob() override;
46 
47  // fcl parameters
48  int port_;
49  int diagLevel_;
50  float onlineRefreshPeriod_;
51  bool keepAlive_;
52  int keepAliveDuration_; // minutes
53  std::string histColor_; // "red"/"blue"/"green"
54 
55  // Member variables
56  TCanvas* canvas_;
57  std::unordered_map<std::string, TH1D*> hists_;
58  THttpServer* server_;
59  std::chrono::time_point<std::chrono::steady_clock> lastUpdate_;
60  std::size_t eventCounter_;
61  std::thread keepAliveThread_;
62 
63 };
64 
65 // Constructor implementation
66 CrvDQM::CrvDQM(fhicl::ParameterSet const& ps)
67  : art::EDAnalyzer(ps)
68  , port_(ps.get<int>("port", 8877))
69  , diagLevel_(ps.get<int>("diagLevel", 1))
70  , onlineRefreshPeriod_(ps.get<float> ("onlineRefreshPeriod", 500)) // ms
71  , keepAlive_(ps.get<int> ("keepAlive", true))
72  , keepAliveDuration_(ps.get<int> ("keepAliveDuration", 5)) // minutes
73  , histColor_(ps.get<std::string> ("histColor", "blue")) // minutes
74 {
75  // Initialise non-fcl member variables
76  eventCounter_ = 0;
77 }
78 
79 // Destructor implementation
80 CrvDQM::~CrvDQM() {
81  // Make sure the server thread is stopped
82  if(keepAliveThread_.joinable()) {
83  keepAliveThread_.join();
84  }
85  // Then clean up ROOT objects
86  if (server_) {
87  delete server_;
88  server_ = nullptr;
89  }
90  if (canvas_) {
91  delete canvas_;
92  canvas_ = nullptr;
93  }
94  for (auto& hist : hists_) {
95  if (hist.second) delete hist.second;
96  hist.second = nullptr;
97  }
98  hists_.clear();
99 }
100 
101 void CrvDQM::beginJob() {
102  // Create HTTP server
103  server_ = new THttpServer(Form("http:%d", port_));
104 
105  // Set global plot style
106  CrvDQMStyle::SetStyle();
107 
108  // Create canvas
109  std::string canvasName = "WebDisplay";
110  canvas_ = new TCanvas(canvasName.c_str(), "CRV web display");
111  canvas_->Divide(2,2); // divide into 4x4 grid
112 
113  // Create histograms
114  hists_["channels"] = new TH1D("Channels", ";Channel;Counts", 64, -0.5, 63.5);
115  hists_["timestamps"] = new TH1D("Timestamps", ";Timestamp;Counts", 256, 0, 255); // 0-0xff
116  hists_["nhits"] = new TH1D("Hits", ";Hits / block;Counts", 61, 0, 60);
117  hists_["adc"] = new TH1D("ADC",";ADC;Counts", 201, -100, 100);
118 
119  // Format and draw
120  int canvasIdx = 1;
121  for (auto& hist : hists_) {
122  // Get pad
123  canvas_->cd(canvasIdx);
124  ++canvasIdx;
125  // Consistent formatting
126  CrvDQMStyle::FormatHist(hist.second, histColor_);
127  // Draw
128  hist.second->Draw();
129  }
130 
131  // Register with server
132  server_->Register("/", canvas_);
133 
134  // Set item defaults
135  server_->SetItemField("/","_monitoring", Form("%f", onlineRefreshPeriod_)); // Update period
136  server_->SetItemField("/","_sidebar", "0");
137  server_->SetItemField("/","_drawitem", canvasName.c_str()); // Set DQM canvas as default item
138  server_->SetItemField("/","_http_cache", "0"); // Disable HTTP caching
139 
140  // Last update variable
141  lastUpdate_ = std::chrono::steady_clock::now();
142 
143  // Print URL
144  printf("Server running on http://localhost:%d/\n", port_);
145 
146  // Start an independent thread for server
147  keepAliveThread_ = std::thread([this]() {
148  while(keepAlive_) {
149  gSystem->ProcessEvents(); // Process events
150  std::this_thread::sleep_for(std::chrono::milliseconds(100)); // Small sleep
151  }
152  });
153 
154 }
155 
156 void CrvDQM::analyze(art::Event const& e) {
157  // Get fragments
158  std::vector<art::Handle<artdaq::Fragments> > fragmentHandles;
159  fragmentHandles = e.getMany<std::vector<artdaq::Fragment> >();
160  artdaq::FragmentPtrs containerFragments;
161  artdaq::Fragments fragments;
162 
163  // Iterate through fragment handles
164  for (const auto& handle : fragmentHandles) {
165  // Catch invalid or empty handles
166  if (!handle.isValid() || handle->empty()) {
167  continue;
168  }
169  // Check if the first object is Container Fragment
170  // ContainerFragment
171  // ├── Fragment 1 (DTCEVT)
172  // └── Fragment n (DTCEVT)
173  if (handle->front().type() == artdaq::Fragment::ContainerFragmentType) {
174  // Iterate through containers
175  for (const auto& cont : *handle) {
176  artdaq::ContainerFragment contf(cont);
177  // Break if this is single fragment rather than a container
178  if (contf.fragment_type() != mu2e::FragmentType::DTCEVT) {
179  break;
180  }
181  // Iterate through fragments in container and fill fragments vector
182  for (size_t i = 0; i < contf.block_count(); ++i) {
183  containerFragments.push_back(contf[i]);
184  fragments.push_back(*containerFragments.back());
185  }
186  }
187  } else { // If the first object in the handle a single fragment
188  if (handle->front().type() == mu2e::FragmentType::DTCEVT) {
189  // Iterate through fragments and fill fragments vector
190  for (auto frag : *handle) {
191  fragments.emplace_back(frag);
192  }
193  }
194  }
195  }
196  if (diagLevel_ > 1) {
197  TLOG(TLVL_INFO) << "[CrvDQM::analyze] Found nFragments" << fragments.size();
198  }
199 
200  // Handle the fragments
201  for (const auto& frag : fragments) {
202  try {
203  mu2e::DTCEventFragment bb(frag);
204  auto data = bb.getData();
205  auto event = &data;
206  if (diagLevel_ > 1) {
207  TLOG(TLVL_INFO) << "Event tag:\t" << "0x" << std::hex << std::setw(4) << std::setfill('0') << event->GetEventWindowTag().GetEventWindowTag(true);
208  }
209  DTCLib::DTC_EventHeader* eventHeader = event->GetHeader();
210 
211  if (diagLevel_ > 1) {
212  TLOG(TLVL_INFO) << eventHeader->toJson() << std::endl
213  << "Subevents count: " << event->GetSubEventCount() << std::endl;
214  }
215 
216  bool plotsUpdated = false;
217  for (unsigned int i = 0; i < event->GetSubEventCount(); ++i) { // In future, use GetSubsystemData to only get CRV subevents
218  DTCLib::DTC_SubEvent& subevent = *(event->GetSubEvent(i));
219  if (diagLevel_ > 1) {
220  TLOG(TLVL_INFO) << "Subevent [" << i << "]:" << std::endl;
221  TLOG(TLVL_INFO) << subevent.GetHeader()->toJson() << std::endl;
222  TLOG(TLVL_INFO) << "Number of Data Block: " << subevent.GetDataBlockCount() << std::endl;
223  }
224 
225  for (size_t bl = 0; bl < subevent.GetDataBlockCount(); ++bl) {
226  ++eventCounter_;
227  auto block = subevent.GetDataBlock(bl);
228  auto blockheader = block->GetHeader();
229  if (diagLevel_ > 1) {
230  TLOG(TLVL_INFO) << blockheader->toJSON() << std::endl;
231  for (int ii = 0; ii < blockheader->GetPacketCount(); ++ii) {
232  TLOG(TLVL_INFO) << DTCLib::DTC_DataPacket(((uint8_t*)block->blockPointer) + ((ii + 1) * 16)).toJSON() << std::endl;
233  }
234  }
235  // Check if we want to decode this data block
236  // Make sure we only process CRV data
237  if(blockheader->GetSubsystem() == 0x2) {
238  if(blockheader->isValid()) {
239  // DQM for version 0x0 data
240  if( blockheader->GetVersion() == 0x0 ) {
241  auto crvData = mu2e::CRVDataDecoder(subevent); // reference
242  //const auto crvStatus = crvData.GetCRVROCStatusPacket(bl);
243  auto hits = crvData.GetCRVHits(bl);
244  for (auto& hit : hits) {
245  // Fill histograms
246  hists_["channels"]->Fill(hit.first.febChannel);
247  hists_["timestamps"]->Fill(hit.first.HitTime);
248  for (auto& adc : hit.second) {
249  hists_["adc"]->Fill(adc.ADC);
250  }
251  }
252  hists_["nhits"]->Fill(hits.size());
253  plotsUpdated = true;
254  }
255  } else {
256  // Iterate invalid count?
257  // ++invalidEventCounter_;
258  continue;
259  }
260  }
261  }
262  }
263 
264  if(plotsUpdated) {
265  auto currentTime = std::chrono::steady_clock::now();
266  std::chrono::duration<double, std::milli> elapsed = currentTime - lastUpdate_;
267  if(elapsed.count() >= onlineRefreshPeriod_) {
268  // Auto scale y-axis
269  for (auto& hist : hists_) {
270  double maxContent = hist.second->GetBinContent(hist.second->GetMaximumBin());
271  hist.second->GetYaxis()->SetRangeUser(0, 1.15*maxContent);
272  }
273  // Update the canvas
274  canvas_->Modified();
275  canvas_->Update();
276  gSystem->ProcessEvents(); // Update display
277  lastUpdate_ = currentTime; // Update the time
278  }
279  }
280 
281  }
282  catch (const std::exception& e) {
283  if (diagLevel_ > 0) {
284  TLOG(TLVL_WARNING) << "Error processing fragment: " << e.what();
285  }
286  continue;
287  }
288  }
289 }
290 
291 // End job printouts
292 void CrvDQM::endJob() {
293  // What am we defining as an "event"? Need to think about it.
294  printf("========================================\n");
295  printf("Processed Events : %zu\n", eventCounter_);
296  // printf("Invalid Events : %zu\n", invalidEventCounter_);
297  // printf("Valid Events : %zu\n", eventCounter_ - invalidEventCounter_);
298  // printf("Error Rate : %.2f%%\n",
299  // (eventCounter_ > 0) ? (100.0 * invalidEventCounter_ / eventCounter_) : 0.0);
300  if (keepAlive_) {
301  printf("Keeping server alive for %i minutes", keepAliveDuration_);
302  printf("\n========================================\n");
303  std::this_thread::sleep_for(std::chrono::minutes(keepAliveDuration_));
304  } else {
305  printf("Killing server");
306  printf("\n========================================\n");
307  }
308  // Server server thread is then stopped by the destructor
309 }
310 
311 DEFINE_ART_MODULE(CrvDQM)
312 }