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"
14 #include "artdaq-core/Data/ContainerFragment.hh"
15 #include "artdaq-core/Data/Fragment.hh"
26 #include "art_root_io/TFileService.h"
31 class MissingDataCheck;
63 void analyze(art::Event
const& e)
override;
76 std::string raw_data_label_;
77 int expected_n_fragments_;
86 unsigned int timeHigh_;
87 unsigned int timeLow_;
89 unsigned int total_n_frags_;
90 unsigned int total_data_size_;
91 unsigned int total_n_CFs_;
92 unsigned int total_n_CFs_missing_;
93 unsigned int total_n_frags_in_CFs_;
94 unsigned int total_n_frags_broken_;
96 unsigned int frag_id_;
98 unsigned int frag_type_;
100 unsigned int frag_data_size_;
106 , raw_data_label_(pset.get<std::string>(
"raw_data_label",
"daq"))
107 , expected_n_fragments_(pset.get<int>(
"expected_n_fragments", -1))
108 , verbosity_(pset.get<int>(
"verbosity", 0))
109 , make_tree_(pset.get<bool>(
"make_tree", true))
113 art::ServiceHandle<art::TFileService> tfs;
114 evtree_ = tfs->make<TTree>(
"evtree",
"MissingDataCheck Event Tree");
115 evtree_->Branch(
"run", &run_,
"run/i");
116 evtree_->Branch(
"subrun", &subrun_,
"subrun/i");
117 evtree_->Branch(
"event", &event_,
"event/i");
118 evtree_->Branch(
"timeHigh", &timeHigh_,
"timeHigh/i");
119 evtree_->Branch(
"timeLow", &timeLow_,
"timeLow/i");
121 evtree_->Branch(
"n_frag_exp", &expected_n_fragments_,
"nfrag_exp/I");
122 evtree_->Branch(
"n_frag", &total_n_frags_,
"n_frag/i");
123 evtree_->Branch(
"data_size", &total_data_size_,
"data_size/i");
124 evtree_->Branch(
"n_cont_frag", &total_n_CFs_,
"n_cont_frag/i");
125 evtree_->Branch(
"n_miss_data", &total_n_CFs_missing_,
"n_miss_data/i");
126 evtree_->Branch(
"n_frag_in_cont", &total_n_frags_in_CFs_,
"n_frag_in_cont/i");
127 evtree_->Branch(
"n_frag_broken", &total_n_frags_broken_,
"n_frag_broken/i");
129 fragtree_ = tfs->make<TTree>(
"fragtree",
"MissingDataCheck Fragment Tree");
130 fragtree_->Branch(
"run", &run_,
"run/i");
131 fragtree_->Branch(
"subrun", &subrun_,
"subrun/i");
132 fragtree_->Branch(
"event", &event_,
"event/i");
133 fragtree_->Branch(
"timeHigh", &timeHigh_,
"timeHigh/i");
134 fragtree_->Branch(
"timeLow", &timeLow_,
"timeLow/i");
136 fragtree_->Branch(
"frag_id", &frag_id_,
"frag_id/i");
137 fragtree_->Branch(
"seq_id", &seq_id_,
"seq_id/i");
138 fragtree_->Branch(
"frag_type", &frag_type_,
"frag_type/i");
139 fragtree_->Branch(
"frag_data_size", &frag_data_size_,
"frag_data_size/i");
140 fragtree_->Branch(
"contained_frags", &contained_frags_,
"contained_frags/I");
141 fragtree_->Branch(
"missing_data", &missing_data_,
"missing_data/I");
148 subrun_ = e.subRun();
150 timeHigh_ = e.time().timeHigh();
151 timeLow_ = e.time().timeLow();
156 std::cout <<
"Processing:"
157 <<
" Run " << e.run()
158 <<
", Subrun " << e.subRun()
159 <<
", Event " << e.event()
160 <<
" (Time=" << e.time().timeHigh() <<
" " << e.time().timeLow() <<
")"
165 std::vector<art::Handle<std::vector<artdaq::Fragment>>> fragmentHandles;
166 fragmentHandles = e.getMany<std::vector<artdaq::Fragment>>();
171 std::cout <<
"\tFound " << fragmentHandles.size() <<
" fragment collections." << std::endl;
174 for (
auto const& h : fragmentHandles)
176 std::cout <<
"\t\tCollection " << h.provenance()->productInstanceName()
177 <<
":\t" << h->size() <<
" fragments." << std::endl;
184 total_data_size_ = 0;
185 for (
auto const& h : fragmentHandles)
187 total_n_frags_ += h->size();
188 for (
auto const& f : *h)
190 total_data_size_ += f.dataSizeBytes();
195 if (expected_n_fragments_ == -1)
197 expected_n_fragments_ = total_n_frags_;
202 std::cout <<
"\tTotal fragments = " << total_n_frags_
203 <<
" / " << expected_n_fragments_ << std::endl;
204 std::cout <<
"\tTotal data size in fragments = " << total_data_size_ << std::endl;
208 total_n_CFs_missing_ = 0;
209 total_n_frags_in_CFs_ = 0;
210 total_n_frags_broken_ = 0;
211 for (
auto const& h : fragmentHandles)
213 std::string instance_name = h.provenance()->productInstanceName();
215 if (instance_name.compare(0, 9,
"Container") == 0)
217 total_n_CFs_ += h->size();
219 for (
auto const& f : *h)
221 artdaq::ContainerFragment cf(f);
222 if (cf.missing_data())
224 ++total_n_CFs_missing_;
226 total_n_frags_in_CFs_ += cf.block_count();
227 frag_id_ = f.fragmentID();
228 seq_id_ = f.sequenceID();
229 frag_type_ = f.type();
230 contained_frags_ = cf.block_count();
231 missing_data_ = cf.missing_data() ? 1 : 0;
232 frag_data_size_ = f.dataSizeBytes();
239 for (
auto const& f : *h)
241 if (instance_name.compare(0, 4,
"Data") == 0 ||
242 instance_name.compare(0, 5,
"Error") == 0 ||
243 instance_name.compare(0, 6,
"Broken") == 0)
245 ++total_n_frags_broken_;
248 frag_id_ = f.fragmentID();
249 seq_id_ = f.sequenceID();
250 frag_type_ = f.type();
251 contained_frags_ = -1;
253 frag_data_size_ = f.dataSizeBytes();
261 std::cout <<
"\tTotal container fragments = " << total_n_CFs_
262 <<
" of which " << total_n_CFs_missing_ <<
" have missing data."
264 std::cout <<
"\tTotal fragments in containers = " << total_n_frags_in_CFs_ << std::endl;
277 std::cout <<
"---------------------------------------------------" << std::endl;
278 std::cout <<
"----------- MISSING DATA CHECK SUMMARY ------------" << std::endl;
279 std::cout <<
"---------------------------------------------------" << std::endl;
281 std::cout <<
"Total events processed: " << evtree_->GetEntries()
282 <<
", expected fragments: " << expected_n_fragments_ << std::endl;
285 evtree_->SetBranchAddress(
"run", &run);
287 evtree_->SetBranchAddress(
"subrun", &subrun);
289 evtree_->SetBranchAddress(
"event", &event);
291 evtree_->SetBranchAddress(
"n_frag_exp", &n_frag_exp);
293 evtree_->SetBranchAddress(
"n_frag", &n_frag);
294 unsigned int n_frag_broken;
295 evtree_->SetBranchAddress(
"n_frag_broken", &n_frag_broken);
296 unsigned int n_miss_data;
297 evtree_->SetBranchAddress(
"n_miss_data", &n_miss_data);
299 std::cout <<
"Events missing fragments:\t\t"
300 << evtree_->GetEntries(
"n_frag<n_frag_exp")
301 <<
" / " << evtree_->GetEntries() << std::endl;
304 if (evtree_->GetEntries(
"n_frag<n_frag_exp") > 0)
306 for (
int i = 0; i < evtree_->GetEntries(); ++i)
308 evtree_->GetEvent(i);
309 if (static_cast<int>(n_frag) < n_frag_exp)
311 std::cout <<
"\tEvent (" << run <<
"," << subrun <<
"," <<
event <<
")"
312 <<
" is missing " << n_frag_exp - n_frag <<
" fragments."
318 std::cout <<
"Events with broken fragments:\t\t"
319 << evtree_->GetEntries(
"n_frag_broken>0")
320 <<
" / " << evtree_->GetEntries() << std::endl;
321 if (evtree_->GetEntries(
"n_frag_broken>0") > 0)
323 for (
int i = 0; i < evtree_->GetEntries(); ++i)
325 evtree_->GetEvent(i);
326 if (n_frag_broken > 0)
328 std::cout <<
"\tEvent (" << run <<
"," << subrun <<
"," <<
event <<
")"
329 <<
" has " << n_frag_broken <<
" fragments."
334 std::cout <<
"Events missing data in fragments:\t"
335 << evtree_->GetEntries(
"n_miss_data>0")
336 <<
" / " << evtree_->GetEntries() << std::endl;
339 if (evtree_->GetEntries(
"n_miss_data>0") > 0)
341 for (
int i = 0; i < evtree_->GetEntries(); ++i)
343 evtree_->GetEvent(i);
346 std::cout <<
"\tEvent (" << run <<
"," << subrun <<
"," <<
event <<
")"
347 <<
" has " << n_miss_data <<
" fragments missing data."
354 std::cout <<
"---------------------------------------------------" << std::endl;
355 std::cout <<
"---------------------------------------------------" << std::endl;
void endJob() override
This method is called at the end of the job, used to print a summary.
Check art::Event for potentially missing data.
~MissingDataCheck() override=default
Default virtual Destructor.
MissingDataCheck(fhicl::ParameterSet const &pset)
MissingDataCheck Constructor.
void analyze(art::Event const &e) override
This method is called for each art::Event in a file or run.