NearlineAna_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // Class: NearlineAna
3 // Module Type: analyzer
4 // File: NearlineAna_module.cc
5 //
6 // Generated at Wed Dec 16 08:25:59 2015 by Jonathan Davies using artmod
7 // from cetpkgsupport v1_10_01.
8 ////////////////////////////////////////////////////////////////////////
9 
17 #include "fhiclcpp/ParameterSet.h"
19 #include "art_root_io/TFileService.h"
20 
21 //cpp
22 #include <sstream>
23 #include <fstream>
24 
25 //ROOT
26 #include "TH1I.h"
27 #include "TTimeStamp.h"
28 #include "TTree.h"
29 
30 //larsoft
35 #include "lardataobj/RawData/raw.h"
36 
37 //dunetpc
39 
40 namespace nearline {
41  class NearlineAna;
42  void BuildTPCChannelMap(std::string channelMapFile, std::map<int,int>& channelMap);
43 }
44 
46 public:
47  explicit NearlineAna(fhicl::ParameterSet const & p);
48  NearlineAna(NearlineAna const &) = delete;
49  NearlineAna(NearlineAna &&) = delete;
50  NearlineAna & operator = (NearlineAna const &) = delete;
51  NearlineAna & operator = (NearlineAna &&) = delete;
52  void reconfigure(fhicl::ParameterSet const & p);
53  void printConfig();
54  void analyze(art::Event const & e) override;
55  void beginJob() override;
56  void endJob() override;
57 
58  size_t getRawDigits(art::Event const & e, art::Handle<std::vector<raw::RawDigit>> & digitHandle);
59  size_t getHits(art::Event const & e, art::Handle<std::vector<recob::Hit>> & hitsHandle);
60 
61  //HitsPerEvent plots
62  void makeHitsPerEventPlots();
63  void fillHitsPerEventPlots(art::Event const & e);
64 
65  //PedestalPerEvent plots
67  void fillPedestalPerEventPlots(art::Event const & e);
69 
70  //PedestalPerTick plots
72  void fillPedestalPerTickPlots(art::Event const & e);
74 
75 private:
76 
78 
81  std::map<int,int> fChannelMap;
82 
84  std::vector<unsigned int> fHitsPerEventChannels;
85  std::vector<TH1I*> fVecHitsPerEventPlots;
87 
91  std::vector<unsigned int> fPedestalPerEventChannels;
92  std::vector<TH1I*> fVecPedestalPerEventPlots;
94 
98  std::vector<unsigned int> fPedestalPerTickChannels;
99  std::vector<TH1I*> fVecPedestalPerTickPlots;
100 
101 
102  // Variables needed for the header info tree:
103  TTree* fHeader;
104  unsigned int fRun;
105  unsigned int fSubrun;
108  int fNevents;
109  unsigned int fStartYear;
110  unsigned int fEndYear;
111  unsigned int fStartMonth;
112  unsigned int fEndMonth;
113  unsigned int fStartDay;
114  unsigned int fEndDay;
115  double fStartHour;
116  double fEndHour;
117  unsigned long long int fStartTime;
118  unsigned long long int fEndTime;
119 
120  //Histogram to store the
122 
123 };
124 
125 ////////////////////////////////////////////////////////////////////////////////
126 
128  :
129  EDAnalyzer(p),
130  fHeader(0),
131  fRun(0),
132  fSubrun(0),
133  fFirstEvent(1e9),
134  fLastEvent(-1),
135  fNevents(0),
136  fStartYear(0),
137  fEndYear(0),
138  fStartMonth(0),
139  fEndMonth(0),
140  fStartDay(0),
141  fEndDay(0),
142  fStartHour(0.0),
143  fEndHour(0.0),
144  fStartTime(-1), // this is an unsigned int so it will default to a huge number
145  fEndTime(0)
146 {
147 
148  {
149  mf::LogInfo logInfo("NearlineAna::NearlineAna");
150  logInfo << "NearlineAna" << "\n";
151  }
152  reconfigure(p);
153 
154  printConfig();
155 
156 }
157 
158 ////////////////////////////////////////////////////////////////////////////////
159 
161 
162  mf::LogInfo logInfo("NearlineAna::reconfigure");
163  logInfo << "reconfigure" << "\n";
164 
165  fVerboseOutput = p.get<bool>("VerboseOutput", false);
166 
167  fUseOnlineChannels = p.get<bool>("UseOnlineChannels", true);
168  fChannelMapFile = p.get<std::string>("TPCChannelMapFile");
169 
170  if(fUseOnlineChannels) BuildTPCChannelMap(fChannelMapFile, fChannelMap);
171 
172  fMakeHitsPerEventPlots = p.get<bool>("MakeHitsPerEventPlots", true);
173  fHitsPerEventChannels = p.get<std::vector<unsigned int>>("HitsPerEventChannels", {1,2,3,4});
174 
175  fHitsTag = p.get<art::InputTag>("HitsTag", "a:b:c");
176 
177  fMakePedestalPerEventPlots = p.get<bool>("MakePedestalPerEventPlots", true);
178  fPedestalPerEventChannels = p.get<std::vector<unsigned int>>("PedestalPerEventChannels", {1,2,3,4});
179  fWritePedestalPerEventFile = p.get<bool>("WritePedestalPerEventFile", false);
180  fPedestalPerEventFileName = p.get<std::string>("PedestalPerEventFileName", "PedestalPerEvent.txt");
181 
182  fMakePedestalPerTickPlots = p.get<bool>("MakePedestalPerTickPlots", true);
183  fPedestalPerTickChannels = p.get<std::vector<unsigned int>>("PedestalPerTickChannels", {1,2,3,4});
184  fWritePedestalPerTickFile = p.get<bool>("WritePedestalPerTickFile", false);
185  fPedestalPerTickFileName = p.get<std::string>("PedestalPerTickFileName", "PedestalPerTick.txt");
186 
187  fRawDigitsTag = p.get<art::InputTag>("RawDigitsTag", "a:b:c");
188 }
189 
190 ////////////////////////////////////////////////////////////////////////////////
191 
193 
194  mf::LogInfo logInfo("NearlineAna::printConfig");
195 
196  logInfo << "fVerboseOutput: " << (fVerboseOutput ? "true" : "false") << "\n";
197  logInfo << "fRawDigitsTag: " << fRawDigitsTag << "\n";
198  logInfo << "fHitsTag: " << fHitsTag << "\n";
199  logInfo << "fUseOnlineChannels: " << (fUseOnlineChannels ? "true" : "false") << "\n";
200 
201  logInfo << "fMakeHitsPerEventPlots: " << (fMakeHitsPerEventPlots ? "true" : "false") << "\n";
202  if(fHitsPerEventChannels.size()){
203  logInfo << "fHitsPerEventChannels";
204  if(fUseOnlineChannels) logInfo << " (online/offline): ";
205  else logInfo << "(offline): ";
206  for(size_t i=0;i<fHitsPerEventChannels.size();i++){
207  logInfo << " " << fHitsPerEventChannels.at(i);
208  if(fUseOnlineChannels) logInfo << "/" << fChannelMap.at(fHitsPerEventChannels.at(i));
209  }//fHitsPerEventChannels
210  logInfo << "\n";
211  }
212 
213  logInfo << "fMakePedestalPerEventPlots: " << (fMakePedestalPerEventPlots ? "true" : "false") << "\n";
214  if(fPedestalPerEventChannels.size()){
215  logInfo << "fPedestalPerEventChannels";
216  if(fUseOnlineChannels) logInfo << " (online/offline): ";
217  else logInfo << "(offline): ";
218  for(size_t i=0;i<fPedestalPerEventChannels.size();i++){
219  logInfo << " " << fPedestalPerEventChannels.at(i);
220  if(fUseOnlineChannels) logInfo << "/" << fChannelMap.at(fPedestalPerEventChannels.at(i));
221  }//fPedestalPerEventChannels
222  logInfo << "\n";
223  }
224  logInfo << "fWritePedestalPerEventFile: " << (fWritePedestalPerEventFile ? "true" : "false") << "\n";
225  if(fWritePedestalPerEventFile) logInfo << "fPedestalPerEventFileName: " << fPedestalPerEventFileName << "\n";
226 
227 
228  logInfo << "fMakePedestalPerTickPlots: " << (fMakePedestalPerTickPlots ? "true" : "false") << "\n";
229  if(fPedestalPerTickChannels.size()){
230  logInfo << "fPedestalPerTickChannels";
231  if(fUseOnlineChannels) logInfo << " (online/offline): ";
232  else logInfo << "(offline): ";
233  for(size_t i=0;i<fPedestalPerTickChannels.size();i++){
234  logInfo << " " << fPedestalPerTickChannels.at(i);
235  if(fUseOnlineChannels) logInfo << "/" << fChannelMap.at(fPedestalPerTickChannels.at(i));
236  }//fPedestalPerTickChannels
237  logInfo << "\n";
238  }
239  logInfo << "fWritePedestalPerTickFile: " << (fWritePedestalPerTickFile ? "true" : "false") << "\n";
240  if(fWritePedestalPerTickFile) logInfo << "fPedestalPerTickFileName: " << fPedestalPerTickFileName << "\n";
241 
242 
243 }
244 
245 ////////////////////////////////////////////////////////////////////////////////
246 
248 
249  mf::LogInfo logInfo("NearlineAna::beginJob");
250  logInfo << "beginJob" << "\n";
251 
252  // book the header object
254  fHeader = tfs->make<TTree>("Header","Subrun Information");
255 
256  fHeader->Branch("Run",&fRun);
257  fHeader->Branch("Subrun",&fSubrun);
258  fHeader->Branch("FirstEvent",&fFirstEvent);
259  fHeader->Branch("LastEvent",&fLastEvent);
260  fHeader->Branch("Nevents",&fNevents);
261  fHeader->Branch("StartYear",&fStartYear);
262  fHeader->Branch("StartMonth",&fStartMonth);
263  fHeader->Branch("StartDay",&fStartDay);
264  fHeader->Branch("StartHour",&fStartHour);
265  fHeader->Branch("EndYear",&fEndYear);
266  fHeader->Branch("EndMonth",&fEndMonth);
267  fHeader->Branch("EndDay",&fEndDay);
268  fHeader->Branch("EndHour",&fEndHour);
269 
270  // Set Nearline Version Number
271  fHistNearlineVersion = tfs->make<TH1I>("hist_nearline_version", "hist_nearline_version", 2, 0, 2);
272  fHistNearlineVersion->GetXaxis()->SetBinLabel(1,"NearlineMinorVersion");
273  fHistNearlineVersion->GetXaxis()->SetBinLabel(2,"NearlineMajorVersion");
274  fHistNearlineVersion->SetBinContent(1, NearlineMinorVersion);
275  fHistNearlineVersion->SetBinContent(2, NearlineMajorVersion);
276 
280 
281 }
282 
283 ////////////////////////////////////////////////////////////////////////////////
284 
286 
287  //
288  // Compute header info.
289  //
290 
291  //
292  // DISECTING the time from evt.time().value() into "human readable" format to display the date & time
293  //
294  unsigned int hour, minute, second;
295  int nano;
296 
297  // Get the time stamp. art::Timestamp::value() returns a TimeValue_t which is a typedef to unsigned long long.
298  // The conventional use is for the upper 32 bits to have the seconds since 1970 epoch and the lower 32 bits to be
299  // the number of nanoseconds with the current second.
300  //
301  // NOTE: It seems that the above is NOT the convention for the 35t events. They only seem to use the lower 32 bits
302  // for the year/month/day/hour/second. For now, I have reversed the values of lup and llo to get the time right.
303  //
304  // THESE VARIABLES WILL NEED TO BE SWITCHED BACK IF USING THE LOWER 32 BITS FOR NANOSECONDS IS IMPLEMENTED!!!
305 
306 
307  const unsigned long int mask32 = 0xFFFFFFFFUL;
308 
309  // taking start time apart
310 
311  // unsigned long int lup = ( fStartTime >> 32 ) & mask32;
312  // unsigned long int llo = fStartTime & mask32;
313  unsigned long int llo = ( fStartTime >> 32 ) & mask32; // reversed value (see above comment)
314  unsigned long int lup = fStartTime & mask32; // reversed value (see above comment)
315  TTimeStamp ts1(lup, (int)llo);
316  ts1.GetDate(kTRUE,0,&fStartYear,&fStartMonth,&fStartDay);
317  ts1.GetTime(kTRUE,0,&hour,&minute,&second);
318  nano = ts1.GetNanoSec();
319  double sec = ((double)second + (double)nano/1.0e9);
320  fStartHour = (double)hour + (double)minute/60.0 + sec/3600.0;
321 
322  // taking end time apart
323  // lup = ( fEndTime >> 32 ) & mask32;
324  // llo = fEndTime & mask32;
325  llo = ( fEndTime >> 32 ) & mask32; // reversed value (see above comment)
326  lup = fEndTime & mask32; // reversed value (see above comment)
327  TTimeStamp ts2(lup, (int)llo);
328  ts2.GetDate(kTRUE,0,&fEndYear,&fEndMonth,&fEndDay);
329  ts2.GetTime(kTRUE,0,&hour,&minute,&second);
330  nano = ts2.GetNanoSec();
331  sec = ((double)second + (double)nano/1.0e9);
332  fEndHour = (double)hour + (double)minute/60.0 + sec/3600.0;
333 
334  fHeader->Fill();
335 
336  //Write out pedestal plots summary to files
339 
340 
341 }
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 
346 {
347 
348  //
349  // Extract event info for the header...
350  //
351  unsigned int run = e.run();
352  unsigned int subrun = e.subRun();
353  unsigned int event = e.id().event();
354  unsigned long long int time = e.time().value();
355 
356  fNevents++;
357  fRun = run;
358  fSubrun = subrun;
359 
360  // Don't assume first/last events are coorelated with start/end times...
361  if(time < fStartTime && e.time() != art::Timestamp::invalidTimestamp()) fStartTime = time;
362  if((int)event < fFirstEvent) fFirstEvent = event;
363  if(time > fEndTime) fEndTime = time;
364  if((int)event > fLastEvent) fLastEvent = event;
365 
366  //
367  // Fill the desired plots.
368  //
372 
373 }
374 
375 ////////////////////////////////////////////////////////////////////////////////
376 
377 size_t nearline::NearlineAna::getRawDigits(art::Event const & e, art::Handle<std::vector<raw::RawDigit>> & digitHandle){
378 
379  bool retVal = e.getByLabel(fRawDigitsTag, digitHandle);
380  if(retVal!=true){
381  mf::LogWarning("NearlineAna::getRawDigits") << "Getting RawDigits FAIL: " << fRawDigitsTag << std::endl;
382  return 0;
383  }
384 
385  try { digitHandle->size(); }
386  catch(std::exception const&) {
387  mf::LogError("NearlineAna::getRawDigits") << "WARNING: Issue with digitHandle for RawDigits" << std::endl;
388  return 0;
389  }
390 
391  if(!digitHandle.isValid()){
392  mf::LogError("NearlineAna::getRawDigits") << "Run: " << e.run()
393  << ", SubRun: " << e.subRun()
394  << ", Event: " << e.event()
395  << " is NOT VALID";
396  throw cet::exception("RawDigit NOT VALID");
397  return 0;
398  }
399 
400  return digitHandle->size();
401 }
402 
403 ////////////////////////////////////////////////////////////////////////////////
404 
405 size_t nearline::NearlineAna::getHits(art::Event const & e, art::Handle<std::vector<recob::Hit>> & hitsHandle){
406 
407  bool retVal = e.getByLabel(fHitsTag, hitsHandle);
408  if(retVal!=true){
409  mf::LogWarning("NearlineAna::getHits") << "Getting Hits FAIL: " << fHitsTag << std::endl;
410  return 0;
411  }
412  try { hitsHandle->size(); }
413  catch(std::exception const&) {
414  mf::LogError("NearlineAna::getHits") << "WARNING: Issue with hitsHandle for Hits" << std::endl;
415  return 0;
416  }
417 
418  if(!hitsHandle.isValid()){
419  mf::LogError("NearlineAna::getHits") << "Run: " << e.run()
420  << ", SubRun: " << e.subRun()
421  << ", Event: " << e.event()
422  << " is NOT VALID";
423  throw cet::exception("Hit NOT VALID");
424  return 0;
425  }
426 
427  return hitsHandle->size();
428 
429 }
430 
431 ////////////////////////////////////////////////////////////////////////////////
432 
434 
435  mf::LogInfo logInfo("NearlineAna::makeHitsPerEventPlots");
436  if(fVerboseOutput) logInfo << "fHitsPerEventChannels:" << "\n";
437 
439 
440  for(auto channel: fHitsPerEventChannels){
441  int numBins = 100;
442  int xmin = 0;
443  int xmax = 3200;
444  std::string hist_name = "hhits_per_event_chan_" + std::to_string(channel);
445  std::string hist_title = "Hits Per Event - Channel "
446  + (fUseOnlineChannels ?
447  "(online/offline) " + std::to_string(channel) + "/" + std::to_string(fChannelMap.at(channel)) :
448  "(offline) " + std::to_string(channel));
449 
450  TH1I* histTemp = tfs->make<TH1I>(hist_name.c_str(), hist_title.c_str(), numBins, xmin, xmax);
451  histTemp->GetXaxis()->SetTitle("Hits per Event");
452  histTemp->GetYaxis()->SetTitle("Events");
453  fVecHitsPerEventPlots.push_back(histTemp);
454  if(fVerboseOutput) logInfo << "channel: " << channel << " hist_name: " << hist_name << " hist_title: " << hist_title << "\n";
455  }
456  if(fVerboseOutput) logInfo << "\n";
457 
458 }
459 
460 ////////////////////////////////////////////////////////////////////////////////
461 
463 
464 
465  mf::LogInfo logInfo("NearlineAna::makePedestalPerEventPlots");
466  if(fVerboseOutput) logInfo << "fPedestalPerEventChannels:" << "\n";
467 
469 
471  int numBins = 128;
472  int xmin = 0;
473  int xmax = 4096;
474  std::string hist_name = "hped_per_event_chan_" + std::to_string(channel);
475  std::string hist_title = "Average ADC Per Event - Channel "
476  + (fUseOnlineChannels ?
477  "(online/offline) " + std::to_string(channel) + "/" + std::to_string(fChannelMap.at(channel)) :
478  "(offline) " + std::to_string(channel));
479 
480  TH1I* histTemp = tfs->make<TH1I>(hist_name.c_str(), hist_title.c_str(), numBins, xmin, xmax);
481  histTemp->GetXaxis()->SetTitle("ADC");
482  histTemp->GetYaxis()->SetTitle("Events");
483  fVecPedestalPerEventPlots.push_back(histTemp);
484  if(fVerboseOutput) logInfo << "channel: " << channel << " hist_name: " << hist_name << " hist_title: " << hist_title << "\n";
485  }
486  if(fVerboseOutput) logInfo << "\n";
487 
488 }
489 
490 ////////////////////////////////////////////////////////////////////////////////
491 
493 
494 
495  mf::LogInfo logInfo("NearlineAna::makePedestalPerTickPlots");
496  if(fVerboseOutput) logInfo << "fPedestalPerTickChannels:" << "\n";
497 
499 
501  int numBins = 128;
502  int xmin = 0;
503  int xmax = 4096;
504  std::string hist_name = "hped_per_tick_chan_" + std::to_string(channel);
505 
506  std::string hist_title = "ADC Per Tick - Channel "
507  + (fUseOnlineChannels ?
508  "(online/offline) " + std::to_string(channel) + "/" + std::to_string(fChannelMap.at(channel)) :
509  "(offline) " + std::to_string(channel));
510 
511  TH1I* histTemp = tfs->make<TH1I>(hist_name.c_str(), hist_title.c_str(), numBins, xmin, xmax);
512  histTemp->GetXaxis()->SetTitle("ADC");
513  histTemp->GetYaxis()->SetTitle("Events");
514  fVecPedestalPerTickPlots.push_back(histTemp);
515  if(fVerboseOutput) logInfo << "channel: " << channel << " hist_name: " << hist_name << " hist_title: " << hist_title << "\n";
516  }
517  if(fVerboseOutput) logInfo << "\n";
518 
519 }
520 
521 ////////////////////////////////////////////////////////////////////////////////
522 
524 
525 
526  mf::LogInfo logInfo("NearlineAna::fillHistPerEventPlots");
528  size_t numHits = getHits(e, hitHandle);
529 
530  //Count the number of hits on channels we are watching
531  std::vector<unsigned int> numHitsPerChannel(fHitsPerEventChannels.size(), 0);
532 
533  if(numHits==0) return;
534 
535  for(size_t hitIter=0;hitIter<numHits;hitIter++){
536  art::Ptr<recob::Hit> hitVec(hitHandle,hitIter);
537  raw::ChannelID_t channel = hitVec->Channel();
538  // geo::View_t view = hitVec->View();
539  // geo::SigType_t signalType = hitVec->SignalType();
540  // geo::WireID wireID = hitVec->WireID();
541  // logInfo << "Hit Channel: " << channel
542  // << " View: " << view
543  // << " SignalType: " << signalType
544  // << " WireID: " << wireID << "\n";
545 
546  for(size_t index=0;index<fHitsPerEventChannels.size();index++){
548  if(this_channel!=channel) continue;
549  numHitsPerChannel.at(index) += 1;
550  }//channel index
551  }//hitIter
552 
553  //Now fill histogram
554  for(size_t index=0;index<fHitsPerEventChannels.size();index++){
555  TH1I* histTemp = fVecHitsPerEventPlots.at(index);
556  histTemp->Fill(numHitsPerChannel.at(index));
557  }//channel index
558 
559  if(fVerboseOutput){
560  logInfo << "Channel (Number of Hits): \n";
561  for(size_t index=0;index<fHitsPerEventChannels.size();index++){
562  logInfo << fHitsPerEventChannels.at(index) << " (" << numHitsPerChannel.at(index) << ") ";
563  }
564  logInfo << "\n";
565  }
566 
567 
568 }
569 
570 ////////////////////////////////////////////////////////////////////////////////
571 
573 
574  mf::LogInfo logInfo("NearlineAna::fillPedestalPerEventPlots");
575 
576  // for(size_t index=0;index<fVecPedestalPerEventPlots.size();index++){
577  // auto channel = fPedestalPerEventChannels.at(index);
578  // auto offline_channel = (fUseOnlineChannels ? fChannelMap.at(channel) : channel );
579  // auto hist = fVecPedestalPerEventPlots.at(index);
580  // logInfo << "channel " << (fUseOnlineChannels ? "(online/offline): " + std::to_string(channel) + "/" + std::to_string(offline_channel) : "(offline): " + std::to_string(offline_channel) );
581  // logInfo << " hist_title: " << hist->GetTitle() << "\n";
582  // }
583 
584 
586  size_t numDigitChans = getRawDigits(e, digitHandle);
587 
588  //Loop through the vector of rawDigits and pick out channels that match our list of channels
589 
590  for(size_t rdIter=0;rdIter<numDigitChans;rdIter++){
591  art::Ptr<raw::RawDigit> digitVec(digitHandle, rdIter);
592  auto channel = digitVec->Channel();
593  for(size_t index=0;index<fPedestalPerEventChannels.size();index++){
595 
596  //Only proceed if this is a channel of interest
597  if(this_channel!=channel) continue;
598 
599  //DEBUG
600  if(fUseOnlineChannels && fVerboseOutput) logInfo << "this_channel (online/offline): " << this_channel << " (" << fPedestalPerEventChannels.at(index) << "/" << fChannelMap.at(fPedestalPerEventChannels.at(index)) << ")\n";
601 
602  auto numSamples = digitVec->Samples();
603  auto compression = digitVec->Compression();
604 
605  //Only proceed if there is a non-zero number of samples
606  if(numSamples==0) continue;
607 
608  //We should uncompress the samples in case there is some form of compression applied
609  raw::RawDigit::ADCvector_t ADCsUncompressed(numSamples);
610  raw::Uncompress(digitVec->ADCs(), ADCsUncompressed, compression);
611 
612  //Calculate the channel ADC average
613  double averageADC=0;
614  for(unsigned int sample=0;sample<numSamples;sample++){
615  averageADC+=ADCsUncompressed[sample];
616  }//sample
617 
618  averageADC/=numSamples;
619 
620  TH1I* histTemp = fVecPedestalPerEventPlots.at(index);
621  histTemp->Fill(averageADC);
622 
623  }//index
624  }//rdIter
625 
626 }
627 
628 
629 ////////////////////////////////////////////////////////////////////////////////
630 
632 
633  mf::LogInfo logInfo("NearlineAna::fillPedestalPerTickPlots");
634 
635  // for(size_t index=0;index<fVecPedestalPerTickPlots.size();index++){
636  // auto channel = fPedestalPerTickChannels.at(index);
637  // auto offline_channel = (fUseOnlineChannels ? fChannelMap.at(channel) : channel );
638  // auto hist = fVecPedestalPerTickPlots.at(index);
639  // logInfo << "channel " << (fUseOnlineChannels ? "(online/offline): " + std::to_string(channel) + "/" + std::to_string(offline_channel) : "(offline): " + std::to_string(offline_channel) );
640  // logInfo << " hist_title: " << hist->GetTitle() << "\n";
641  // }
642 
643 
645  size_t numDigitChans = getRawDigits(e, digitHandle);
646 
647  //Loop through the vector of rawDigits and pick out channels that match our list of channels
648 
649  for(size_t rdIter=0;rdIter<numDigitChans;rdIter++){
650  art::Ptr<raw::RawDigit> digitVec(digitHandle, rdIter);
651  auto channel = digitVec->Channel();
652  for(size_t index=0;index<fPedestalPerTickChannels.size();index++){
654 
655  //Only proceed if this is a channel of interest
656  if(this_channel!=channel) continue;
657 
658  //DEBUG
659  if(fUseOnlineChannels && fVerboseOutput) logInfo << "this_channel (online/offline): " << this_channel << " (" << fPedestalPerTickChannels.at(index) << "/" << fChannelMap.at(fPedestalPerTickChannels.at(index)) << ")\n";
660 
661  auto numSamples = digitVec->Samples();
662  auto compression = digitVec->Compression();
663 
664  //Only proceed if there is a non-zero number of samples
665  if(numSamples==0) continue;
666 
667  //We should uncompress the samples in case there is some form of compression applied
668  raw::RawDigit::ADCvector_t ADCsUncompressed(numSamples);
669  raw::Uncompress(digitVec->ADCs(), ADCsUncompressed, compression);
670 
671  TH1I* histTemp = fVecPedestalPerTickPlots.at(index);
672 
673  for(unsigned int sample=0;sample<numSamples;sample++){
674  histTemp->Fill(ADCsUncompressed[sample]);
675  }//sample
676  }//index
677  }//rdIter
678 
679 }
680 
681 ////////////////////////////////////////////////////////////////////////////////
682 
684  std::ostringstream my_ostream;
685  my_ostream << "online_channel " << "offline_channel " << "pedestal_mean " << "pedestal_rms " << "\n";
686 
687  for(size_t index=0;index<fPedestalPerEventChannels.size();index++){
688  auto online_channel = -1;
689  auto offline_channel = fPedestalPerEventChannels.at(index);
690  if(fUseOnlineChannels){
691  online_channel = fPedestalPerEventChannels.at(index);
692  offline_channel = fChannelMap.at(fPedestalPerEventChannels.at(index));
693  }
694  TH1I* histTemp = fVecPedestalPerEventPlots.at(index);
695  double mean = histTemp->GetMean();
696  double rms = histTemp->GetRMS();
697  my_ostream << online_channel << " " << offline_channel << " " << mean << " " << rms << "\n";
698  }//index
699 
700  std::ofstream outFile(fileName);
701  if(outFile.is_open()) outFile << my_ostream.str();
702  else mf::LogWarning("writePedestalPerEventSummaryFile") << "FAILED to open file: " << fileName;
703  outFile.close();
704 
705  //DEBUG
706  mf::LogInfo("writePedestalPerEventSummaryFile") << my_ostream.str();
707 }
708 
709 ////////////////////////////////////////////////////////////////////////////////
710 
712  std::ostringstream my_ostream;
713  my_ostream << "online_channel " << "offline_channel " << "pedestal_mean " << "pedestal_rms " << "\n";
714 
715  for(size_t index=0;index<fPedestalPerTickChannels.size();index++){
716  auto online_channel = -1;
717  auto offline_channel = fPedestalPerTickChannels.at(index);
718  if(fUseOnlineChannels){
719  online_channel = fPedestalPerTickChannels.at(index);
720  offline_channel = fChannelMap.at(fPedestalPerTickChannels.at(index));
721  }
722  TH1I* histTemp = fVecPedestalPerTickPlots.at(index);
723  double mean = histTemp->GetMean();
724  double rms = histTemp->GetRMS();
725  my_ostream << online_channel << " " << offline_channel << " " << mean << " " << rms << "\n";
726  }//index
727 
728  std::ofstream outFile(fileName);
729  if(outFile.is_open()) outFile << my_ostream.str();
730  else mf::LogWarning("writePedestalPerEventSummaryFile") << "FAILED to open file: " << fileName;
731  outFile.close();
732 
733  //DEBUG
734  mf::LogInfo("writePedestalPerTickSummaryFile") << my_ostream.str();
735 
736 }
737 
738 ////////////////////////////////////////////////////////////////////////////////
739 
740 
741 void nearline::BuildTPCChannelMap(std::string channelMapFile, std::map<int,int>& channelMap) {
742 
743  /// Builds TPC channel map from the map txt file
744 
745  channelMap.clear();
746 
747  int onlineChannel;
748  int offlineChannel;
749 
750  std::string fullname;
751  cet::search_path sp("FW_SEARCH_PATH");
752  sp.find_file(channelMapFile, fullname);
753 
754  if (fullname.empty())
755  mf::LogWarning("DAQToOffline") << "Input TPC channel map file " << channelMapFile << " not found in FW_SEARCH_PATH. Using online channel numbers!" << std::endl;
756 
757  else {
758  mf::LogVerbatim("DAQToOffline") << "Build TPC Online->Offline channel Map from " << fullname;
759  std::ifstream infile(fullname);
760  while (infile.good()) {
761  infile >> onlineChannel >> offlineChannel;
762  channelMap.insert(std::make_pair(onlineChannel,offlineChannel));
763  mf::LogVerbatim("DAQToOffline") << " " << onlineChannel << " -> " << offlineChannel;
764  }
765  std::cout << "channelMap has size " << channelMap.size() << ". If this is 2048, then it's fine even if the above lines skipped a 'few' channels..." << std::endl;
766  }
767 
768 }
769 
void analyze(art::Event const &e) override
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
const ADCvector_t & ADCs() const
Reference to the compressed ADC count vector.
Definition: RawDigit.h:210
EventNumber_t event() const
Definition: DataViewImpl.cc:96
void fillPedestalPerTickPlots(art::Event const &e)
ULong64_t Samples() const
Number of samples in the uncompressed ADC data.
Definition: RawDigit.h:213
unsigned int subrun
const int NearlineMajorVersion
double rms(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:39
size_t getHits(art::Event const &e, art::Handle< std::vector< recob::Hit >> &hitsHandle)
NearlineAna & operator=(NearlineAna const &)=delete
std::string string
Definition: nybbler.cc:12
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Declaration of signal hit object.
void fillPedestalPerEventPlots(art::Event const &e)
std::vector< TH1I * > fVecPedestalPerEventPlots
ChannelID_t Channel() const
DAQ channel this raw data was read from.
Definition: RawDigit.h:212
std::vector< short > ADCvector_t
Type representing a (compressed) vector of ADC counts.
Definition: RawDigit.h:73
Definition of basic raw digits.
std::vector< unsigned int > fPedestalPerTickChannels
void fillHitsPerEventPlots(art::Event const &e)
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:27
constexpr TimeValue_t value() const
Definition: Timestamp.h:23
std::vector< unsigned int > fPedestalPerEventChannels
void reconfigure(fhicl::ParameterSet const &p)
std::string fPedestalPerTickFileName
const int NearlineMinorVersion
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:446
const double e
Timestamp time() const
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:69
string infile
std::vector< TH1I * > fVecHitsPerEventPlots
size_t getRawDigits(art::Event const &e, art::Handle< std::vector< raw::RawDigit >> &digitHandle)
Collect all the RawData header files together.
T get(std::string const &key) const
Definition: ParameterSet.h:231
unsigned long long int fStartTime
SubRunNumber_t subRun() const
Definition: DataViewImpl.cc:89
RunNumber_t run() const
Definition: DataViewImpl.cc:82
std::map< int, int > fChannelMap
raw::Compress_t Compression() const
Compression algorithm used to store the ADC counts.
Definition: RawDigit.h:216
NearlineAna(fhicl::ParameterSet const &p)
static constexpr Timestamp invalidTimestamp()
Definition: Timestamp.h:82
p
Definition: test.py:223
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
std::string find_file(std::string const &filename) const
Definition: search_path.cc:94
std::string fPedestalPerEventFileName
std::vector< unsigned int > fHitsPerEventChannels
EventNumber_t event() const
Definition: EventID.h:116
unsigned long long int fEndTime
void writePedestalPerTickSummaryFile(std::string fileName)
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:27
std::vector< TH1I * > fVecPedestalPerTickPlots
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:755
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:80
void writePedestalPerEventSummaryFile(std::string fileName)
void BuildTPCChannelMap(std::string channelMapFile, std::map< int, int > &channelMap)
raw::ChannelID_t Channel() const
ID of the readout channel the hit was extracted from.
Definition: Hit.h:230
double mean(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:15
std::string to_string(ModuleType const mt)
Definition: ModuleType.h:34
EventID id() const
Definition: Event.cc:37
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
Event finding and building.
unsigned int run