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

Public Member Functions

 CRTSimValidation (fhicl::ParameterSet const &p)
 
 CRTSimValidation (CRTSimValidation const &)=delete
 
 CRTSimValidation (CRTSimValidation &&)=delete
 
CRTSimValidationoperator= (CRTSimValidation const &)=delete
 
CRTSimValidationoperator= (CRTSimValidation &&)=delete
 
void analyze (art::Event const &e) override
 
void beginJob () override
 
- 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 (SharedResources const &resources)
 
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 &)
 
fhicl::ParameterSetID selectorConfig () 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 fCRTLabel
 
TH1D * fDetsWithHits
 
TH1D * fAllADCs
 
TH1D * fDevWithinTrigger
 
TH1D * fDeltaTAssoc
 
TH1D * fTriggerDeltaT
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
std::string const & processName () const
 
bool wantAllEvents () const noexcept
 
bool wantEvent (ScheduleID id, Event const &e) const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &select_paths, std::vector< std::string > const &reject_paths, fhicl::ParameterSet const &config)
 
- 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 46 of file CRTSimValidation_module.cc.

Constructor & Destructor Documentation

CRT::CRTSimValidation::CRTSimValidation ( fhicl::ParameterSet const &  p)
explicit

Definition at line 85 of file CRTSimValidation_module.cc.

86  :
87  EDAnalyzer(p), fCRTLabel(p.get<art::InputTag>("CRTLabel"))
88 {
89  consumes<std::vector<CRT::Trigger>>(fCRTLabel);
90  consumes<std::vector<art::Assns<sim::AuxDetSimChannel, CRT::Trigger>>>(fCRTLabel);
91 }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
p
Definition: test.py:223
CRT::CRTSimValidation::CRTSimValidation ( CRTSimValidation const &  )
delete
CRT::CRTSimValidation::CRTSimValidation ( CRTSimValidation &&  )
delete

Member Function Documentation

void CRT::CRTSimValidation::analyze ( art::Event const &  e)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 93 of file CRTSimValidation_module.cc.

94 {
95  //Get the data products I plan to work with
96  const auto& triggers = e.getValidHandle<std::vector<CRT::Trigger>>(fCRTLabel);
97 
98  art::FindManyP<sim::AuxDetSimChannel> trigToSim(triggers, e, fCRTLabel);
99 
100  //Get a handle to the Geometry service to look up AuxDetGeos from module numbers
102 
103  //Mapping from channel to previous Trigger time
104  std::unordered_map<size_t, double> prevTimes;
105  for(const auto& trigger: *triggers)
106  {
107  fDetsWithHits->Fill(geom->AuxDet(trigger.Channel()).Name().c_str(), 1.0);
108  const auto& hits = trigger.Hits();
109  for(const auto& hit: hits) fAllADCs->Fill(hit.ADC());
110 
111  const auto found = prevTimes.find(trigger.Channel());
112  if(found != prevTimes.end())
113  {
114  fTriggerDeltaT->Fill(trigger.Timestamp() - found->second);
115  found->second = trigger.Timestamp();
116  }
117  else prevTimes[trigger.Channel()] = trigger.Timestamp();
118  }
119 
120  for(size_t trigIt = 0; trigIt < trigToSim.size(); ++trigIt)
121  {
122  const auto& trig = (*triggers)[trigIt];
123  const auto& simChans = trigToSim.at(trigIt);
124 
125  size_t nIDEs = 0;
126  const auto sumSq = std::accumulate(simChans.begin(), simChans.end(), 0.,
127  [this, &trig, &nIDEs](double sum, const auto& channel)
128  {
129  const auto& ides = channel->AuxDetIDEs();
130  nIDEs += ides.size();
131  return sum + std::accumulate(ides.begin(), ides.end(), 0.,
132  [this, &trig](double subSum, const auto& ide)
133  {
134  const auto diff = trig.Timestamp() -
135  (ide.exitT + ide.entryT)/2.;
136  fDeltaTAssoc->Fill(diff);
137  if(fabs(diff) > 200) //TODO: Read pset from original module to use readout
138  // time here
139  {
140  mf::LogInfo("LargeDeltaT") << "Found large deltaT: " << diff << "\n"
141  << "Timestamp is " << trig.Timestamp() << "\n"
142  << "True time is " << (ide.exitT +
143  ide.entryT)/2. << ".\n";
144  }
145 
146  return subSum + diff*diff;
147  });
148  });
149  fDevWithinTrigger->Fill(std::sqrt(sumSq/(nIDEs-1))); //Root of the mean of squares
150  }
151 }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
ChannelGroupService::Name Name
uint8_t channel
Definition: CRTFragment.hh:201
const double e
AuxDetGeo const & AuxDet(unsigned int const ad=0) const
Returns the specified auxiliary detector.
Detector simulation of raw signals on wires.
void CRT::CRTSimValidation::beginJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 153 of file CRTSimValidation_module.cc.

154 {
155  //Register ROOT objects I want to write with the TFileService
157  fDetsWithHits = tfs->make<TH1D>("DetsWithHits", "AuxDetGeo Names that Had CRT::Triggers;AuxDetGeo Name;Events", 1, 0, -1);
158  //Trigger automatic binning since am plotting strings.
159  fAllADCs = tfs->make<TH1D>("AllADCs", "ADC Values from CRT::Hits on All Channels;ADC Value;Hits", 500, 0, 2000);
160  fDevWithinTrigger = tfs->make<TH1D>("DevWithinTrigger", "Deviation of True Times from Associated Trigger Timestamp;"
161  "True Time #sigma;Energy Deposits", 100, 0, 300);
162  fDeltaTAssoc = tfs->make<TH1D>("DeltaTAssoc", "Time Difference Between Associated IDEs and Trigger Timestamp;#DeltaT_{IDE};Energy Deposits",
163  200, -100, 100);
164  fTriggerDeltaT = tfs->make<TH1D>("TriggerDeltaT", "#DeltaT Between Consecutive Triggers in a Module;#DeltaT;Trigger Pairs", 20, 0, 30);
165 }
CRTSimValidation& CRT::CRTSimValidation::operator= ( CRTSimValidation const &  )
delete
CRTSimValidation& CRT::CRTSimValidation::operator= ( CRTSimValidation &&  )
delete

Member Data Documentation

TH1D* CRT::CRTSimValidation::fAllADCs
private

Definition at line 73 of file CRTSimValidation_module.cc.

art::InputTag CRT::CRTSimValidation::fCRTLabel
private

Definition at line 67 of file CRTSimValidation_module.cc.

TH1D* CRT::CRTSimValidation::fDeltaTAssoc
private

Definition at line 79 of file CRTSimValidation_module.cc.

TH1D* CRT::CRTSimValidation::fDetsWithHits
private

Definition at line 72 of file CRTSimValidation_module.cc.

TH1D* CRT::CRTSimValidation::fDevWithinTrigger
private

Definition at line 77 of file CRTSimValidation_module.cc.

TH1D* CRT::CRTSimValidation::fTriggerDeltaT
private

Definition at line 80 of file CRTSimValidation_module.cc.


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