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/Fragment.hh"
15 #include "artdaq-core/Data/ContainerFragment.hh"
26 #include "art/Framework/Services/Optional/TFileService.h"
31 class MissingDataCheck;
63 void analyze(art::Event
const& e)
override;
71 std::string raw_data_label_;
72 int expected_n_fragments_;
81 unsigned int timeHigh_;
82 unsigned int timeLow_;
84 unsigned int total_n_frags_;
85 unsigned int total_data_size_;
86 unsigned int total_n_CFs_;
87 unsigned int total_n_CFs_missing_;
88 unsigned int total_n_frags_in_CFs_;
89 unsigned int total_n_frags_broken_;
92 unsigned int frag_id_;
94 unsigned int frag_type_;
96 unsigned int frag_data_size_;
104 , raw_data_label_(pset.get<std::string>(
"raw_data_label",
"daq"))
105 , expected_n_fragments_(pset.get<int>(
"expected_n_fragments",-1))
106 , verbosity_(pset.get<int>(
"verbosity",0))
107 , make_tree_(pset.get<bool>(
"make_tree",true))
110 art::ServiceHandle<art::TFileService> tfs;
111 evtree_ = tfs->make<TTree>(
"evtree",
"MissingDataCheck Event Tree");
112 evtree_->Branch(
"run",&run_,
"run/i");
113 evtree_->Branch(
"subrun",&subrun_,
"subrun/i");
114 evtree_->Branch(
"event",&event_,
"event/i");
115 evtree_->Branch(
"timeHigh",&timeHigh_,
"timeHigh/i");
116 evtree_->Branch(
"timeLow",&timeLow_,
"timeLow/i");
118 evtree_->Branch(
"n_frag_exp",&expected_n_fragments_,
"nfrag_exp/I");
119 evtree_->Branch(
"n_frag",&total_n_frags_,
"n_frag/i");
120 evtree_->Branch(
"data_size",&total_data_size_,
"data_size/i");
121 evtree_->Branch(
"n_cont_frag",&total_n_CFs_,
"n_cont_frag/i");
122 evtree_->Branch(
"n_miss_data",&total_n_CFs_missing_,
"n_miss_data/i");
123 evtree_->Branch(
"n_frag_in_cont",&total_n_frags_in_CFs_,
"n_frag_in_cont/i");
124 evtree_->Branch(
"n_frag_broken",&total_n_frags_broken_,
"n_frag_broken/i");
126 fragtree_ = tfs->make<TTree>(
"fragtree",
"MissingDataCheck Fragment Tree");
127 fragtree_->Branch(
"run",&run_,
"run/i");
128 fragtree_->Branch(
"subrun",&subrun_,
"subrun/i");
129 fragtree_->Branch(
"event",&event_,
"event/i");
130 fragtree_->Branch(
"timeHigh",&timeHigh_,
"timeHigh/i");
131 fragtree_->Branch(
"timeLow",&timeLow_,
"timeLow/i");
133 fragtree_->Branch(
"frag_id",&frag_id_,
"frag_id/i");
134 fragtree_->Branch(
"seq_id",&seq_id_,
"seq_id/i");
135 fragtree_->Branch(
"frag_type",&frag_type_,
"frag_type/i");
136 fragtree_->Branch(
"frag_data_size",&frag_data_size_,
"frag_data_size/i");
137 fragtree_->Branch(
"contained_frags",&contained_frags_,
"contained_frags/I");
138 fragtree_->Branch(
"missing_data",&missing_data_,
"missing_data/I");
146 subrun_ = e.subRun();
148 timeHigh_ = e.time().timeHigh();
149 timeLow_ = e.time().timeLow();
153 std::cout <<
"Processing:"
154 <<
" Run " << e.run()
155 <<
", Subrun " << e.subRun()
156 <<
", Event " << e.event()
157 <<
" (Time=" << e.time().timeHigh() <<
" " << e.time().timeLow() <<
")"
161 std::vector< art::Handle< std::vector<artdaq::Fragment> > > fragmentHandles;
162 e.getManyByType(fragmentHandles);
166 std::cout <<
"\tFound " << fragmentHandles.size() <<
" fragment collections." << std::endl;
168 for(
auto const& h : fragmentHandles)
169 std::cout <<
"\t\tCollection " << h.provenance()->productInstanceName()
170 <<
":\t" << h->size() <<
" fragments." << std::endl;
177 for(
auto const& h : fragmentHandles){
178 total_n_frags_ += h->size();
179 for(
auto const& f : *h)
180 total_data_size_ += f.dataSizeBytes();
184 if(expected_n_fragments_==-1)
185 expected_n_fragments_=total_n_frags_;
188 std::cout <<
"\tTotal fragments = " << total_n_frags_
189 <<
" / " << expected_n_fragments_ << std::endl;
190 std::cout <<
"\tTotal data size in fragments = " << total_data_size_ << std::endl;
195 total_n_CFs_missing_=0;
196 total_n_frags_in_CFs_=0;
197 total_n_frags_broken_=0;
198 for(
auto const& h : fragmentHandles){
199 std::string instance_name = h.provenance()->productInstanceName();
201 if(instance_name.compare(0,9,
"Container")==0){
202 total_n_CFs_ += h->size();
204 for(
auto const& f : *h){
205 artdaq::ContainerFragment cf(f);
206 if(cf.missing_data()){
207 ++total_n_CFs_missing_;
209 total_n_frags_in_CFs_ += cf.block_count();
210 frag_id_ = f.fragmentID();
211 seq_id_ = f.sequenceID();
212 frag_type_ = f.type();
213 contained_frags_ = cf.block_count();
214 missing_data_ = cf.missing_data()?1:0;
215 frag_data_size_ = f.dataSizeBytes();
222 for(
auto const& f : *h){
224 if(instance_name.compare(0,4,
"Data")==0 ||
225 instance_name.compare(0,5,
"Error")==0 ||
226 instance_name.compare(0,6,
"Broken")==0)
227 ++total_n_frags_broken_;
229 frag_id_ = f.fragmentID();
230 seq_id_ = f.sequenceID();
231 frag_type_ = f.type();
232 contained_frags_ = -1;
234 frag_data_size_ = f.dataSizeBytes();
243 std::cout <<
"\tTotal container fragments = " << total_n_CFs_
244 <<
" of which " << total_n_CFs_missing_ <<
" have missing data."
246 std::cout <<
"\tTotal fragments in containers = " << total_n_frags_in_CFs_ << std::endl;
258 std::cout <<
"---------------------------------------------------" << std::endl;
259 std::cout <<
"----------- MISSING DATA CHECK SUMMARY ------------" << std::endl;
260 std::cout <<
"---------------------------------------------------" << std::endl;
262 std::cout <<
"Total events processed: " << evtree_->GetEntries()
263 <<
", expected fragments: " << expected_n_fragments_ << std::endl;
265 unsigned int run; evtree_->SetBranchAddress(
"run",&run);
266 unsigned int subrun; evtree_->SetBranchAddress(
"subrun",&subrun);
267 unsigned int event; evtree_->SetBranchAddress(
"event",&event);
268 int n_frag_exp; evtree_->SetBranchAddress(
"n_frag_exp",&n_frag_exp);
269 unsigned int n_frag; evtree_->SetBranchAddress(
"n_frag",&n_frag);
270 unsigned int n_frag_broken; evtree_->SetBranchAddress(
"n_frag_broken",&n_frag_broken);
271 unsigned int n_miss_data; evtree_->SetBranchAddress(
"n_miss_data",&n_miss_data);
273 std::cout <<
"Events missing fragments:\t\t"
274 << evtree_->GetEntries(
"n_frag<n_frag_exp")
275 <<
" / " << evtree_->GetEntries() << std::endl;
277 if(evtree_->GetEntries(
"n_frag<n_frag_exp")>0){
278 for(
int i=0; i<evtree_->GetEntries(); ++i)
280 evtree_->GetEvent(i);
281 if((
int)n_frag<n_frag_exp){
282 std::cout <<
"\tEvent (" << run <<
"," << subrun <<
"," <<
event <<
")"
283 <<
" is missing " << n_frag_exp-n_frag <<
" fragments."
289 std::cout <<
"Events with broken fragments:\t\t"
290 << evtree_->GetEntries(
"n_frag_broken>0")
291 <<
" / " << evtree_->GetEntries() << std::endl;
292 if(evtree_->GetEntries(
"n_frag_broken>0")>0){
293 for(
int i=0; i<evtree_->GetEntries(); ++i)
295 evtree_->GetEvent(i);
296 if(n_frag_broken > 0){
297 std::cout <<
"\tEvent (" << run <<
"," << subrun <<
"," <<
event <<
")"
298 <<
" has " << n_frag_broken <<
" fragments."
303 std::cout <<
"Events missing data in fragments:\t"
304 << evtree_->GetEntries(
"n_miss_data>0")
305 <<
" / " << evtree_->GetEntries() << std::endl;
307 if(evtree_->GetEntries(
"n_miss_data>0")>0){
308 for(
int i=0; i<evtree_->GetEntries(); ++i)
310 evtree_->GetEvent(i);
312 std::cout <<
"\tEvent (" << run <<
"," << subrun <<
"," <<
event <<
")"
313 <<
" has " << n_miss_data <<
" fragments missing data."
320 std::cout <<
"---------------------------------------------------" << std::endl;
321 std::cout <<
"---------------------------------------------------" << std::endl;
void endJob() override
This method is called at the end of the job, used to print a summary.
virtual ~MissingDataCheck()=default
Default virtual Destructor.
Check art::Event for potentially missing data.
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.