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

Public Member Functions

 TPCHits (fhicl::ParameterSet const &p)
 
 TPCHits (TPCHits const &)=delete
 
 TPCHits (TPCHits &&)=delete
 
TPCHitsoperator= (TPCHits const &)=delete
 
TPCHitsoperator= (TPCHits &&)=delete
 
void analyze (art::Event const &e) override
 
void reconfigure (fhicl::ParameterSet const &p)
 
void endJob ()
 
void makeHistograms ()
 
void printRawDigitTimingInfo (art::Event const &e)
 
void printTpcFragTimingInfo (art::Event const &e)
 
- 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

art::InputTag fTPCHitTag
 
art::InputTag fRawDigitsTag
 
bool fPrintRawDigitTimingInfo
 
art::EventID fEventIDLast
 
art::Timestamp fTimestampLast
 
art::EventNumber_t fEventLast
 
bool fPrintTpcFragTimingInfo
 
std::vector< unsigned int > fTotalNumHitsPerTPC
 
std::vector< TH1I * > fHistNumHitsPerTPC
 
std::vector< TH1I * > fHistNumTicksPerEvent
 
std::vector< TH1I * > fHistNumTicksPerEventPair
 
std::vector< TH2I * > fHistNumTicksPerEventVsLast
 
std::vector< unsigned int > fNumChansTPC
 
std::vector< unsigned int > fNumTicksTPC
 
std::vector< unsigned int > fNumChansTPCLastEvent
 
std::vector< unsigned int > fNumTicksTPCLastEvent
 
std::vector< unsigned int > fNumChansTPCTwoEvents
 
std::vector< unsigned int > fNumTicksTPCTwoEvents
 

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 44 of file TPCHits_module.cc.

Constructor & Destructor Documentation

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

Definition at line 96 of file TPCHits_module.cc.

97  :
98  EDAnalyzer(p) // ,
99  // More initializers here.
100 {
101  reconfigure(p);
102 
104  fTotalNumHitsPerTPC = std::vector<unsigned int>(geom->NTPC());
105  fNumChansTPC = std::vector<unsigned int>(geom->NTPC());
106  fNumTicksTPC = std::vector<unsigned int>(geom->NTPC());
107  fNumChansTPCLastEvent = std::vector<unsigned int>(geom->NTPC());
108  fNumTicksTPCLastEvent = std::vector<unsigned int>(geom->NTPC());
109  fNumChansTPCTwoEvents = std::vector<unsigned int>(geom->NTPC());
110  fNumTicksTPCTwoEvents = std::vector<unsigned int>(geom->NTPC());
111 
112 
113  makeHistograms();
114 
115 }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:27
std::vector< unsigned int > fNumChansTPCTwoEvents
void reconfigure(fhicl::ParameterSet const &p)
std::vector< unsigned int > fNumTicksTPC
std::vector< unsigned int > fNumTicksTPCTwoEvents
std::vector< unsigned int > fNumChansTPC
std::vector< unsigned int > fNumChansTPCLastEvent
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
p
Definition: test.py:223
std::vector< unsigned int > fNumTicksTPCLastEvent
std::vector< unsigned int > fTotalNumHitsPerTPC
nearline::TPCHits::TPCHits ( TPCHits const &  )
delete
nearline::TPCHits::TPCHits ( TPCHits &&  )
delete

Member Function Documentation

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

Implements art::EDAnalyzer.

Definition at line 299 of file TPCHits_module.cc.

300 {
301 
304  return;//FIXME
305 
307  std::vector<unsigned int> thisNumHits(geom->NTPC());
310 
311  bool retVal = e.getByLabel(fTPCHitTag, hitHandle);
312  if(retVal==true)
313  ;
314  //mf::LogInfo("TPCHits") << "Getting Hits SUCCESS: " << fTPCHitTag << std::endl;
315  else{
316  mf::LogWarning("TPCHits") << "Getting Hits FAIL: " << fTPCHitTag << std::endl;
317  return;
318  }
319 
320  try { hitHandle->size(); }
321  catch(std::exception e) {
322  mf::LogError("TPCHits") << "WARNING: Issue with hitHandle for TPC hits" << std::endl;
323  return;
324  }
325 
326  if(!hitHandle.isValid()){
327  mf::LogError("TPCHits") << "Run: " << e.run()
328  << ", SubRun: " << e.subRun()
329  << ", Event: " << e.event()
330  << " is NOT VALID";
331  throw cet::exception("hits NOT VALID");
332  return;
333  }
334 
335 
336  retVal = e.getByLabel(fRawDigitsTag, digitHandle);
337  if(retVal==true)
338  ;
339  //mf::LogInfo("TPCHits") << "Getting RawDigits SUCCESS: " << fRawDigitsTag << std::endl;
340  else{
341  mf::LogWarning("TPCHits") << "Getting RawDigits FAIL: " << fRawDigitsTag << std::endl;
342  return;
343  }
344 
345  try { digitHandle->size(); }
346  catch(std::exception e) {
347  mf::LogError("TPCHits") << "WARNING: Issue with digitHandle for RawDigits" << std::endl;
348  return;
349  }
350 
351  if(!digitHandle.isValid()){
352  mf::LogError("TPCHits") << "Run: " << e.run()
353  << ", SubRun: " << e.subRun()
354  << ", Event: " << e.event()
355  << " is NOT VALID";
356  throw cet::exception("RawDigit NOT VALID");
357  return;
358  }
359 
360  std::vector<art::Ptr<recob::Hit> > hitlist;
361  art::fill_ptr_vector(hitlist, hitHandle);
362 
363  // size_t numHits = hitlist.size();
364  // mf::LogInfo("TPCHits") << "NumHits: " << numHits << " hits" << std::endl;
365  for (size_t tpchit_index = 0; tpchit_index<hitlist.size(); ++tpchit_index){
366  geo::WireID wireid = hitlist[tpchit_index]->WireID();
367  thisNumHits.at(wireid.TPC)++;
368  fTotalNumHitsPerTPC.at(wireid.TPC)++;
369  }//hit
370 
371  for(size_t tpc=0; tpc < fHistNumHitsPerTPC.size();tpc++){
372  TH1I* tempHist = fHistNumHitsPerTPC.at(tpc);
373  if(thisNumHits.at(tpc)==0)
374  tempHist->Fill(-1);
375  else
376  tempHist->Fill(TMath::Log10(thisNumHits.at(tpc)));
377  }//tpc
378 
379  size_t numDigitChans = digitHandle->size();
380  // mf::LogInfo("TPCHits") << "numDigitChans: " << numDigitChans << std::endl;
381  // size_t numChans = geom->Nchannels();
382  // mf::LogInfo("TPCHits") << "numChans: " << numChans << std::endl;
383 
384  fNumChansTPC = std::vector<unsigned int> (geom->NTPC());
385  fNumTicksTPC = std::vector<unsigned int> (geom->NTPC());
386  for(size_t rdIter=0;rdIter<numDigitChans;rdIter++){
387  art::Ptr<raw::RawDigit> digitVec(digitHandle, rdIter);
388  auto channel = digitVec->Channel();
389  auto numSamples = digitVec->Samples();
390  std::vector<geo::WireID> wids = geom->ChannelToWire(channel);
391 
392  if(wids.at(0).Plane !=2 ) continue;
393  fNumChansTPC.at(wids.at(0).TPC)++;
394  fNumTicksTPC.at(wids.at(0).TPC)+=numSamples;
395  fNumChansTPCTwoEvents.at(wids.at(0).TPC)++;
396  fNumTicksTPCTwoEvents.at(wids.at(0).TPC)+=numSamples;
397 
398 
399  }//rdIter
400 
401  // mf::LogError loginfo("TPCHits");
402  // loginfo << "Numbers of channels\n";
403  // for(size_t tpc=0;tpc<fNumChansTPC.size();tpc++){
404  // if(tpc!=5) continue;
405  // loginfo << "tpc: " << tpc
406  // << " numChannels: " << fNumChansTPC.at(tpc)
407  // << " numSamples: " << fNumTicksTPC.at(tpc)
408  // << " numSamplesTwoEvents: " << fNumTicksTPCTwoEvents.at(tpc)
409  // << "\n";
410 
411  // }
412 
413  for(size_t tpc=0;tpc<fHistNumTicksPerEvent.size();tpc++){
414 
415  unsigned int numChans = fNumChansTPC.at(tpc);
416  unsigned int numTicks = fNumTicksTPC.at(tpc);
417  unsigned int numTicksLastEvent = fNumTicksTPCLastEvent.at(tpc);
418  // unsigned int numTicksEventPair = fNumTicksTPCTwoEvents.at(tpc);
419 
420  // if(numTicks>0){
421  // mf::LogInfo loginfo("TPCHits");
422  // loginfo << "tpc: " << tpc << "\n"
423  // << " numTicks: " << numTicks << "\n"
424  // << " numTicksEventPair: " << numTicksEventPair << "\n";
425  // }
426 
427  if(numChans==0) continue;
428  TH1I* tempHist = fHistNumTicksPerEvent.at(tpc);
429  tempHist->Fill(numTicks/numChans);
430 
431  tempHist = fHistNumTicksPerEventPair.at(tpc);
432  // tempHist->Fill(numTicksEventPair/numChans);
433  tempHist->Fill((numTicksLastEvent+numTicks)/numChans);
434 
435  TH2I* tempHist2D = fHistNumTicksPerEventVsLast.at(tpc);
436  tempHist2D->Fill(numTicks/numChans, numTicksLastEvent/numChans);
437 
438  }//tpc
439 
444 
445 }
std::vector< TH1I * > fHistNumTicksPerEventPair
std::vector< TH2I * > fHistNumTicksPerEventVsLast
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
std::vector< TH1I * > fHistNumHitsPerTPC
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
std::vector< unsigned int > fNumChansTPCTwoEvents
void printTpcFragTimingInfo(art::Event const &e)
bool isValid() const
Definition: Handle.h:183
const double e
void printRawDigitTimingInfo(art::Event const &e)
std::vector< unsigned int > fNumTicksTPC
std::vector< TH1I * > fHistNumTicksPerEvent
std::vector< unsigned int > fNumTicksTPCTwoEvents
std::vector< unsigned int > fNumChansTPC
std::vector< unsigned int > fNumChansTPCLastEvent
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
art::InputTag fTPCHitTag
std::vector< unsigned int > fNumTicksTPCLastEvent
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
art::InputTag fRawDigitsTag
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:291
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:402
std::vector< unsigned int > fTotalNumHitsPerTPC
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
void nearline::TPCHits::endJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 447 of file TPCHits_module.cc.

447  {
448 
449  std::ostringstream os;
450 
451  os << "====================================" << "\n"
452  << "End Job" << "\n"
453  << "====================================" << "\n";
454 
455  for(size_t i=0; i<fTotalNumHitsPerTPC.size(); i++) os << "i: " << i << " contents: " << fTotalNumHitsPerTPC.at(i) << std::endl;
456 
457  os << "====================================" << "\n";
458  mf::LogInfo("TPCHits") << os.str();
459 
460 
461 }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::vector< unsigned int > fTotalNumHitsPerTPC
QTextStream & endl(QTextStream &s)
void nearline::TPCHits::makeHistograms ( )

Definition at line 117 of file TPCHits_module.cc.

117  {
118 
119  mf::LogInfo loginfo("TPCHits");
120  loginfo << "====================================" << "\n"
121  << "Making Histograms" << "\n"
122  << "====================================" << "\n";
123 
126  auto const* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
127 
128  const double samplingRate = detprop->SamplingRate();
129  const unsinged int millisliceSize = 5e6;
130  const unsigned int maxNumTicksPerEvent = millisliceSize / samplingRate;
131 
132  fHistNumHitsPerTPC = std::vector<TH1I*> (geom->NTPC());
133  fHistNumTicksPerEvent = std::vector<TH1I*>(geom->NTPC());
134  fHistNumTicksPerEventPair = std::vector<TH1I*>(geom->NTPC());
135  fHistNumTicksPerEventVsLast = std::vector<TH2I*>(geom->NTPC());
136 
137  //Currently 5ms of data per pair of events
138  //Sampling rate of the TPC is 2MHz (0.5 ns spacing) -> can get this from the det properties
139  //That means that 5000 / 0.5 = 10,000 ticks per 2 events
140  //Max number of hits is therefore == 10,000 * nwires per TPC
141 
142  loginfo << "Making Number of Hits per TPC Histograms\n";
143 
144  for(size_t tpc = 0; tpc < fHistNumHitsPerTPC.size(); tpc++){
145  std::ostringstream histName;
146  histName << "hist_num_hits_tpc_" << tpc;
147 
148  std::ostringstream histTitle;
149  histTitle << "Number of hits per event TPC " << tpc << ";Log_{10}(Number of hits);Number of events";
150 
151  const unsigned int numPlanes = geom->Nplanes(tpc);
152  unsigned int totalNumWires = 0;
153  for(unsigned int plane=0;plane<numPlanes;plane++){
154  unsigned int numWires = geom->Nwires(plane, tpc);
155  totalNumWires += numWires;
156  }//plane
157 
158  unsigned int maxNumHits = totalNumWires * maxNumTicksPerEvent;
159 
160  int numBins = 100;
161  int binLow = -1;
162  int binHigh = TMath::Log10(maxNumHits);
163  TH1I* tempHist = tfs->make<TH1I>(histName.str().c_str(), histTitle.str().c_str(), numBins, binLow, binHigh);
164  fHistNumHitsPerTPC.at(tpc) = tempHist;
165  loginfo << "histName: " << histName.str() << "\n";
166  }
167 
168  loginfo << "\n"
169  << "Making Number of ticks / wire / event / TPC Histograms\n";
170 
171  for(size_t tpc = 0; tpc < fHistNumTicksPerEvent.size(); tpc++){
172  std::ostringstream histName;
173  histName << "hist_num_ticks_per_event_tpc_" << tpc;
174 
175  std::ostringstream histTitle;
176  histTitle << "Number of ticks / wire / event in TPC " << tpc << ";Number of ticks;Number of events";
177 
178  int numBins = 100;
179  int binLow = 0;
180  int binHigh = maxNumTicksPerEvent;
181  TH1I* tempHist = tfs->make<TH1I>(histName.str().c_str(), histTitle.str().c_str(), numBins, binLow, binHigh);
182  fHistNumTicksPerEvent.at(tpc) = tempHist;
183  loginfo << "histName: " << histName.str() << "\n";
184  }//tpc
185 
186  loginfo << "\n"
187  << "Making Number of ticks / wire / event pair / TPC Histograms\n";
188 
189  for(size_t tpc = 0; tpc < fHistNumTicksPerEventPair.size(); tpc++){
190  std::ostringstream histName;
191  histName << "hist_num_ticks_per_event_pair_tpc_" << tpc;
192 
193  std::ostringstream histTitle;
194  histTitle << "Number of ticks / wire / event pair in TPC " << tpc << ";Number of ticks;Number of event pairs";
195 
196  int numBins = 100;
197  int binLow = 0;
198  int binHigh = maxNumTicksPerEvent;
199  TH1I* tempHist = tfs->make<TH1I>(histName.str().c_str(), histTitle.str().c_str(), numBins, binLow, binHigh);
200  fHistNumTicksPerEventPair.at(tpc) = tempHist;
201  loginfo << "histName: " << histName.str() << "\n";
202  }//tpc
203 
204  loginfo << "\n"
205  << "Making Number of ticks / wire / event vs last / TPC Histograms\n";
206 
207  for(size_t tpc = 0; tpc < fHistNumTicksPerEventVsLast.size(); tpc++){
208  std::ostringstream histName;
209  histName << "hist_num_ticks_per_event_vs_last_tpc_" << tpc;
210 
211  std::ostringstream histTitle;
212  histTitle << "Number of ticks / wire / event vs last event in TPC " << tpc << ";Number of ticks;Number of ticks (last)";
213 
214  int numBins = 100;
215  int binLow = 0;
216  int binHigh = maxNumTicksPerEvent;
217  TH2I* tempHist = tfs->make<TH2I>(histName.str().c_str(), histTitle.str().c_str(), numBins, binLow, binHigh, numBins, binLow, binHigh);
218  fHistNumTicksPerEventVsLast.at(tpc) = tempHist;
219  loginfo << "histName: " << histName.str() << "\n";
220  }//tpc
221 
222  loginfo << "====================================" << "\n"
223  << "Finished Making Histograms" << "\n"
224  << "====================================" << "\n";
225 }
std::vector< TH1I * > fHistNumTicksPerEventPair
std::vector< TH2I * > fHistNumTicksPerEventVsLast
std::vector< TH1I * > fHistNumHitsPerTPC
unsigned int Nwires(unsigned int p, unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wires in the specified plane.
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
std::vector< TH1I * > fHistNumTicksPerEvent
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
TPCHits& nearline::TPCHits::operator= ( TPCHits const &  )
delete
TPCHits& nearline::TPCHits::operator= ( TPCHits &&  )
delete
void nearline::TPCHits::printRawDigitTimingInfo ( art::Event const &  e)

Definition at line 278 of file TPCHits_module.cc.

278  {
279 
280  art::EventID eventID = e.id();
281  art::Timestamp timestamp = e.time();
282  art::EventNumber_t event = e.event();
283 
284  mf::LogInfo("TPCHits::printRawDigitTimingInfo") << "EventID: " << eventID
285  << " time: " << timestamp.value()
286  << " EventNum: " << event
287  << std::endl;
288 
289  if((timestamp.value() - fTimestampLast.value() < 0) || (event - fEventLast < 0)){
290  mf::LogInfo("TPCHits::printRawDigitTimingInfo") << "ERROR - negative delta timestamp of eventnum" << std::endl;
291  }
292 
293  fEventIDLast = eventID;
294  fTimestampLast = timestamp;
295  fEventLast = event;
296 }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
constexpr TimeValue_t value() const
Definition: Timestamp.h:23
const double e
IDNumber_t< Level::Event > EventNumber_t
Definition: IDNumber.h:118
art::Timestamp fTimestampLast
art::EventNumber_t fEventLast
QTextStream & endl(QTextStream &s)
Event finding and building.
art::EventID fEventIDLast
void nearline::TPCHits::printTpcFragTimingInfo ( art::Event const &  e)

Definition at line 244 of file TPCHits_module.cc.

244  {
245 
246  mf::LogInfo("TPCHits::printTpcFragTimingInfo") << "Hello World" << std::endl;
247 
248  art::Handle<artdaq::Fragments> rawFragments;
249  art::InputTag tpcFragTag("daq:TPC");
250  bool retVal = e.getByLabel("daq:TPC", rawFragments);
251  if(retVal==true)
252  mf::LogInfo("TPCHits::printTpcFragTimingInfo") << "Getting TPC Frag SUCCESS: " << tpcFragTag << std::endl;
253  else{
254  mf::LogWarning("TPCHits::printTpcFragTimingInfo") << "Getting TPC Frag FAIL: " << tpcFragTag << std::endl;
255  return;
256  }
257  try { rawFragments->size(); }
258  catch(std::exception e) {
259  mf::LogError("TPCHits::printTpcFragTimingInfo") << "WARNING: Issue with rawFragments for TPC hits" << std::endl;
260  return;
261  }
262 
263  for(size_t fragIndex = 0; fragIndex < rawFragments->size(); fragIndex++){
264  const artdaq::Fragment &singleFragment = rawFragments->at(fragIndex);
265  lbne::TpcMilliSliceFragment msf(singleFragment);
266  auto nMicroSlices = msf.microSliceCount();
267 
268  mf::LogInfo("TPCHits::printTpcFragTimingInfo") << "fragIndex: " << fragIndex
269  << " fragmentID: " << singleFragment.fragmentID()
270  << " ms counter: " << nMicroSlices
271  << std::endl;
272 
273  }//fragIndex
274 }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
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::TPCHits::reconfigure ( fhicl::ParameterSet const &  p)

Definition at line 227 of file TPCHits_module.cc.

227  {
228 
229  fTPCHitTag = p.get<art::InputTag>("TPCHitTag", "a:b:c");
230  fRawDigitsTag = p.get<art::InputTag>("RawDigitsTag", "a:b:c");
231  fPrintRawDigitTimingInfo = p.get<bool>("fPrintRawDigitTimingInfo", true);
232  fPrintTpcFragTimingInfo = p.get<bool>("fPrintTpcFragTimingInfo", true);
233 
234  mf::LogInfo("TPCHits") << "====================================" << "\n"
235  << "Parameter Set" << "\n"
236  << "====================================" << "\n"
237  << "fTPCHitTag: " << fTPCHitTag << "\n"
238  << "fRawDigitsTag: " << fRawDigitsTag << "\n"
239  << "fPrintRawDigitTimingInfo: " << fPrintRawDigitTimingInfo << "\n"
240  << "fPrintTpcFragTimingInfo: " << fPrintTpcFragTimingInfo << "\n"
241  << "====================================" << "\n";
242 }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
p
Definition: test.py:223
art::InputTag fTPCHitTag
art::InputTag fRawDigitsTag

Member Data Documentation

art::EventID nearline::TPCHits::fEventIDLast
private

Definition at line 71 of file TPCHits_module.cc.

art::EventNumber_t nearline::TPCHits::fEventLast
private

Definition at line 73 of file TPCHits_module.cc.

std::vector<TH1I*> nearline::TPCHits::fHistNumHitsPerTPC
private

Definition at line 79 of file TPCHits_module.cc.

std::vector<TH1I*> nearline::TPCHits::fHistNumTicksPerEvent
private

Definition at line 80 of file TPCHits_module.cc.

std::vector<TH1I*> nearline::TPCHits::fHistNumTicksPerEventPair
private

Definition at line 81 of file TPCHits_module.cc.

std::vector<TH2I*> nearline::TPCHits::fHistNumTicksPerEventVsLast
private

Definition at line 82 of file TPCHits_module.cc.

std::vector<unsigned int> nearline::TPCHits::fNumChansTPC
private

Definition at line 85 of file TPCHits_module.cc.

std::vector<unsigned int> nearline::TPCHits::fNumChansTPCLastEvent
private

Definition at line 87 of file TPCHits_module.cc.

std::vector<unsigned int> nearline::TPCHits::fNumChansTPCTwoEvents
private

Definition at line 89 of file TPCHits_module.cc.

std::vector<unsigned int> nearline::TPCHits::fNumTicksTPC
private

Definition at line 86 of file TPCHits_module.cc.

std::vector<unsigned int> nearline::TPCHits::fNumTicksTPCLastEvent
private

Definition at line 88 of file TPCHits_module.cc.

std::vector<unsigned int> nearline::TPCHits::fNumTicksTPCTwoEvents
private

Definition at line 90 of file TPCHits_module.cc.

bool nearline::TPCHits::fPrintRawDigitTimingInfo
private

Definition at line 68 of file TPCHits_module.cc.

bool nearline::TPCHits::fPrintTpcFragTimingInfo
private

Definition at line 75 of file TPCHits_module.cc.

art::InputTag nearline::TPCHits::fRawDigitsTag
private

Definition at line 67 of file TPCHits_module.cc.

art::Timestamp nearline::TPCHits::fTimestampLast
private

Definition at line 72 of file TPCHits_module.cc.

std::vector<unsigned int> nearline::TPCHits::fTotalNumHitsPerTPC
private

Definition at line 78 of file TPCHits_module.cc.

art::InputTag nearline::TPCHits::fTPCHitTag
private

Definition at line 66 of file TPCHits_module.cc.


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