Public Member Functions | Private Attributes | List of all members
nearline::NearlineAna Class Reference
Inheritance diagram for nearline::NearlineAna:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 NearlineAna (fhicl::ParameterSet const &p)
 
 NearlineAna (NearlineAna const &)=delete
 
 NearlineAna (NearlineAna &&)=delete
 
NearlineAnaoperator= (NearlineAna const &)=delete
 
NearlineAnaoperator= (NearlineAna &&)=delete
 
void reconfigure (fhicl::ParameterSet const &p)
 
void printConfig ()
 
void analyze (art::Event const &e) override
 
void beginJob () override
 
void endJob () override
 
size_t getRawDigits (art::Event const &e, art::Handle< std::vector< raw::RawDigit >> &digitHandle)
 
size_t getHits (art::Event const &e, art::Handle< std::vector< recob::Hit >> &hitsHandle)
 
void makeHitsPerEventPlots ()
 
void fillHitsPerEventPlots (art::Event const &e)
 
void makePedestalPerEventPlots ()
 
void fillPedestalPerEventPlots (art::Event const &e)
 
void writePedestalPerEventSummaryFile (std::string fileName)
 
void makePedestalPerTickPlots ()
 
void fillPedestalPerTickPlots (art::Event const &e)
 
void writePedestalPerTickSummaryFile (std::string fileName)
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob ()
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
std::string const & processName () const
 
bool wantAllEvents () const
 
bool wantEvent (Event const &e)
 
fhicl::ParameterSetID selectorConfig () const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Attributes

bool fVerboseOutput
 
bool fUseOnlineChannels
 
std::string fChannelMapFile
 
std::map< int, int > fChannelMap
 
bool fMakeHitsPerEventPlots
 
std::vector< unsigned int > fHitsPerEventChannels
 
std::vector< TH1I * > fVecHitsPerEventPlots
 
art::InputTag fHitsTag
 
bool fMakePedestalPerEventPlots
 
bool fWritePedestalPerEventFile
 
std::string fPedestalPerEventFileName
 
std::vector< unsigned int > fPedestalPerEventChannels
 
std::vector< TH1I * > fVecPedestalPerEventPlots
 
art::InputTag fRawDigitsTag
 
bool fMakePedestalPerTickPlots
 
bool fWritePedestalPerTickFile
 
std::string fPedestalPerTickFileName
 
std::vector< unsigned int > fPedestalPerTickChannels
 
std::vector< TH1I * > fVecPedestalPerTickPlots
 
TTree * fHeader
 
unsigned int fRun
 
unsigned int fSubrun
 
int fFirstEvent
 
int fLastEvent
 
int fNevents
 
unsigned int fStartYear
 
unsigned int fEndYear
 
unsigned int fStartMonth
 
unsigned int fEndMonth
 
unsigned int fStartDay
 
unsigned int fEndDay
 
double fStartHour
 
double fEndHour
 
unsigned long long int fStartTime
 
unsigned long long int fEndTime
 
TH1I * fHistNearlineVersion
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &paths, fhicl::ParameterSet const &config)
 
detail::ProcessAndEventSelectorsprocessAndEventSelectors ()
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 45 of file NearlineAna_module.cc.

Constructor & Destructor Documentation

nearline::NearlineAna::NearlineAna ( fhicl::ParameterSet const &  p)
explicit

Definition at line 127 of file NearlineAna_module.cc.

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 }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:27
void reconfigure(fhicl::ParameterSet const &p)
unsigned long long int fStartTime
p
Definition: test.py:223
unsigned long long int fEndTime
nearline::NearlineAna::NearlineAna ( NearlineAna const &  )
delete
nearline::NearlineAna::NearlineAna ( NearlineAna &&  )
delete

Member Function Documentation

void nearline::NearlineAna::analyze ( art::Event const &  e)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 345 of file NearlineAna_module.cc.

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 }
void fillPedestalPerTickPlots(art::Event const &e)
unsigned int subrun
void fillPedestalPerEventPlots(art::Event const &e)
void fillHitsPerEventPlots(art::Event const &e)
const double e
unsigned long long int fStartTime
static constexpr Timestamp invalidTimestamp()
Definition: Timestamp.h:82
unsigned long long int fEndTime
Event finding and building.
unsigned int run
void nearline::NearlineAna::beginJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 247 of file NearlineAna_module.cc.

247  {
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 }
const int NearlineMajorVersion
const int NearlineMinorVersion
void nearline::NearlineAna::endJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 285 of file NearlineAna_module.cc.

285  {
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 }
std::string fPedestalPerTickFileName
unsigned long long int fStartTime
std::string fPedestalPerEventFileName
unsigned long long int fEndTime
void writePedestalPerTickSummaryFile(std::string fileName)
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:80
void writePedestalPerEventSummaryFile(std::string fileName)
void nearline::NearlineAna::fillHitsPerEventPlots ( art::Event const &  e)

Definition at line 523 of file NearlineAna_module.cc.

523  {
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 }
size_t getHits(art::Event const &e, art::Handle< std::vector< recob::Hit >> &hitsHandle)
const double e
std::vector< TH1I * > fVecHitsPerEventPlots
std::map< int, int > fChannelMap
std::vector< unsigned int > fHitsPerEventChannels
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:27
void nearline::NearlineAna::fillPedestalPerEventPlots ( art::Event const &  e)

Definition at line 572 of file NearlineAna_module.cc.

572  {
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 }
std::vector< TH1I * > fVecPedestalPerEventPlots
std::vector< short > ADCvector_t
Type representing a (compressed) vector of ADC counts.
Definition: RawDigit.h:73
std::vector< unsigned int > fPedestalPerEventChannels
const double e
size_t getRawDigits(art::Event const &e, art::Handle< std::vector< raw::RawDigit >> &digitHandle)
std::map< int, int > fChannelMap
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:755
void nearline::NearlineAna::fillPedestalPerTickPlots ( art::Event const &  e)

Definition at line 631 of file NearlineAna_module.cc.

631  {
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 }
std::vector< short > ADCvector_t
Type representing a (compressed) vector of ADC counts.
Definition: RawDigit.h:73
std::vector< unsigned int > fPedestalPerTickChannels
const double e
size_t getRawDigits(art::Event const &e, art::Handle< std::vector< raw::RawDigit >> &digitHandle)
std::map< int, int > fChannelMap
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
size_t nearline::NearlineAna::getHits ( art::Event const &  e,
art::Handle< std::vector< recob::Hit >> &  hitsHandle 
)

Definition at line 405 of file NearlineAna_module.cc.

405  {
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 }
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
bool isValid() const
Definition: Handle.h:183
const double e
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
size_t nearline::NearlineAna::getRawDigits ( art::Event const &  e,
art::Handle< std::vector< raw::RawDigit >> &  digitHandle 
)

Definition at line 377 of file NearlineAna_module.cc.

377  {
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 }
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
bool isValid() const
Definition: Handle.h:183
const double e
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
void nearline::NearlineAna::makeHitsPerEventPlots ( )

Definition at line 433 of file NearlineAna_module.cc.

433  {
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 }
std::string string
Definition: nybbler.cc:12
std::vector< TH1I * > fVecHitsPerEventPlots
std::map< int, int > fChannelMap
std::vector< unsigned int > fHitsPerEventChannels
std::string to_string(ModuleType const mt)
Definition: ModuleType.h:34
void nearline::NearlineAna::makePedestalPerEventPlots ( )

Definition at line 462 of file NearlineAna_module.cc.

462  {
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 }
std::string string
Definition: nybbler.cc:12
std::vector< TH1I * > fVecPedestalPerEventPlots
std::vector< unsigned int > fPedestalPerEventChannels
std::map< int, int > fChannelMap
std::string to_string(ModuleType const mt)
Definition: ModuleType.h:34
void nearline::NearlineAna::makePedestalPerTickPlots ( )

Definition at line 492 of file NearlineAna_module.cc.

492  {
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 }
std::string string
Definition: nybbler.cc:12
std::vector< unsigned int > fPedestalPerTickChannels
std::map< int, int > fChannelMap
std::vector< TH1I * > fVecPedestalPerTickPlots
std::string to_string(ModuleType const mt)
Definition: ModuleType.h:34
NearlineAna& nearline::NearlineAna::operator= ( NearlineAna const &  )
delete
NearlineAna& nearline::NearlineAna::operator= ( NearlineAna &&  )
delete
void nearline::NearlineAna::printConfig ( )

Definition at line 192 of file NearlineAna_module.cc.

192  {
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 }
std::vector< unsigned int > fPedestalPerTickChannels
std::vector< unsigned int > fPedestalPerEventChannels
std::string fPedestalPerTickFileName
std::map< int, int > fChannelMap
std::string fPedestalPerEventFileName
std::vector< unsigned int > fHitsPerEventChannels
void nearline::NearlineAna::reconfigure ( fhicl::ParameterSet const &  p)

Definition at line 160 of file NearlineAna_module.cc.

160  {
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 }
std::string string
Definition: nybbler.cc:12
std::vector< unsigned int > fPedestalPerTickChannels
std::vector< unsigned int > fPedestalPerEventChannels
std::string fPedestalPerTickFileName
std::map< int, int > fChannelMap
p
Definition: test.py:223
std::string fPedestalPerEventFileName
std::vector< unsigned int > fHitsPerEventChannels
void BuildTPCChannelMap(std::string channelMapFile, std::map< int, int > &channelMap)
void nearline::NearlineAna::writePedestalPerEventSummaryFile ( std::string  fileName)

Definition at line 683 of file NearlineAna_module.cc.

683  {
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 }
double rms(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:39
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::vector< TH1I * > fVecPedestalPerEventPlots
std::vector< unsigned int > fPedestalPerEventChannels
std::map< int, int > fChannelMap
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
double mean(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:15
void nearline::NearlineAna::writePedestalPerTickSummaryFile ( std::string  fileName)

Definition at line 711 of file NearlineAna_module.cc.

711  {
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 }
double rms(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:39
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::vector< unsigned int > fPedestalPerTickChannels
std::map< int, int > fChannelMap
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
std::vector< TH1I * > fVecPedestalPerTickPlots
double mean(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:15

Member Data Documentation

std::map<int,int> nearline::NearlineAna::fChannelMap
private

Definition at line 81 of file NearlineAna_module.cc.

std::string nearline::NearlineAna::fChannelMapFile
private

Definition at line 80 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fEndDay
private

Definition at line 114 of file NearlineAna_module.cc.

double nearline::NearlineAna::fEndHour
private

Definition at line 116 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fEndMonth
private

Definition at line 112 of file NearlineAna_module.cc.

unsigned long long int nearline::NearlineAna::fEndTime
private

Definition at line 118 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fEndYear
private

Definition at line 110 of file NearlineAna_module.cc.

int nearline::NearlineAna::fFirstEvent
private

Definition at line 106 of file NearlineAna_module.cc.

TTree* nearline::NearlineAna::fHeader
private

Definition at line 103 of file NearlineAna_module.cc.

TH1I* nearline::NearlineAna::fHistNearlineVersion
private

Definition at line 121 of file NearlineAna_module.cc.

std::vector<unsigned int> nearline::NearlineAna::fHitsPerEventChannels
private

Definition at line 84 of file NearlineAna_module.cc.

art::InputTag nearline::NearlineAna::fHitsTag
private

Definition at line 86 of file NearlineAna_module.cc.

int nearline::NearlineAna::fLastEvent
private

Definition at line 107 of file NearlineAna_module.cc.

bool nearline::NearlineAna::fMakeHitsPerEventPlots
private

Definition at line 83 of file NearlineAna_module.cc.

bool nearline::NearlineAna::fMakePedestalPerEventPlots
private

Definition at line 88 of file NearlineAna_module.cc.

bool nearline::NearlineAna::fMakePedestalPerTickPlots
private

Definition at line 95 of file NearlineAna_module.cc.

int nearline::NearlineAna::fNevents
private

Definition at line 108 of file NearlineAna_module.cc.

std::vector<unsigned int> nearline::NearlineAna::fPedestalPerEventChannels
private

Definition at line 91 of file NearlineAna_module.cc.

std::string nearline::NearlineAna::fPedestalPerEventFileName
private

Definition at line 90 of file NearlineAna_module.cc.

std::vector<unsigned int> nearline::NearlineAna::fPedestalPerTickChannels
private

Definition at line 98 of file NearlineAna_module.cc.

std::string nearline::NearlineAna::fPedestalPerTickFileName
private

Definition at line 97 of file NearlineAna_module.cc.

art::InputTag nearline::NearlineAna::fRawDigitsTag
private

Definition at line 93 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fRun
private

Definition at line 104 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fStartDay
private

Definition at line 113 of file NearlineAna_module.cc.

double nearline::NearlineAna::fStartHour
private

Definition at line 115 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fStartMonth
private

Definition at line 111 of file NearlineAna_module.cc.

unsigned long long int nearline::NearlineAna::fStartTime
private

Definition at line 117 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fStartYear
private

Definition at line 109 of file NearlineAna_module.cc.

unsigned int nearline::NearlineAna::fSubrun
private

Definition at line 105 of file NearlineAna_module.cc.

bool nearline::NearlineAna::fUseOnlineChannels
private

Definition at line 79 of file NearlineAna_module.cc.

std::vector<TH1I*> nearline::NearlineAna::fVecHitsPerEventPlots
private

Definition at line 85 of file NearlineAna_module.cc.

std::vector<TH1I*> nearline::NearlineAna::fVecPedestalPerEventPlots
private

Definition at line 92 of file NearlineAna_module.cc.

std::vector<TH1I*> nearline::NearlineAna::fVecPedestalPerTickPlots
private

Definition at line 99 of file NearlineAna_module.cc.

bool nearline::NearlineAna::fVerboseOutput
private

Definition at line 77 of file NearlineAna_module.cc.

bool nearline::NearlineAna::fWritePedestalPerEventFile
private

Definition at line 89 of file NearlineAna_module.cc.

bool nearline::NearlineAna::fWritePedestalPerTickFile
private

Definition at line 96 of file NearlineAna_module.cc.


The documentation for this class was generated from the following file: