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

Public Member Functions

 MakeSNeTProfileHistos (fhicl::ParameterSet const &p)
 
 MakeSNeTProfileHistos (MakeSNeTProfileHistos const &)=delete
 
 MakeSNeTProfileHistos (MakeSNeTProfileHistos &&)=delete
 
MakeSNeTProfileHistosoperator= (MakeSNeTProfileHistos const &)=delete
 
MakeSNeTProfileHistosoperator= (MakeSNeTProfileHistos &&)=delete
 
void analyze (art::Event const &e) override
 
void reconfigure (fhicl::ParameterSet const &p)
 
void beginJob () override
 
void endJob () 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

std::string fTruthLabel
 label for process that made the MCTruth info More...
 
std::string fHitLabel
 label for process that made the reco hits More...
 
TH1F * fEnergySpec
 1D event energy spectrum (used only for backgrounds) More...
 
TH1F * fTrigEff_numerator
 how many times a trigger was issued per true event energy More...
 
TH1F * fTrigEff
 trigger efficiency per true event energy (fTrigEff_numerator/fEnergySpec) More...
 
TH2F * fNHitsVsEnergy
 number of reconstructed hits per true event energy More...
 
TH2F * fSumADCVsEnergy
 summed ADC for all reco hits per true event energy More...
 
TH1F * fTPSplitProb
 probability of splitting a trigger primitive into two (or more) TPs per true event energy More...
 
TH1F * fNMCTruths
 number of MCTruths per art::event More...
 
TH1F * fCCNC
 CCNC for each MCTruth. More...
 
TH1F * fOrigin
 mctruth->Origin() More...
 
TH1F * fNuPDG
 mc_neutrino->Nu()->PdgCode() More...
 

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 43 of file MakeSNeTProfileHistos_module.cc.

Constructor & Destructor Documentation

MakeSNeTProfileHistos::MakeSNeTProfileHistos ( fhicl::ParameterSet const &  p)
explicit

Definition at line 95 of file MakeSNeTProfileHistos_module.cc.

96  :
97  EDAnalyzer(p)
98 {
99  this->reconfigure(p);
100 }
void reconfigure(fhicl::ParameterSet const &p)
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
p
Definition: test.py:223
MakeSNeTProfileHistos::MakeSNeTProfileHistos ( MakeSNeTProfileHistos const &  )
delete
MakeSNeTProfileHistos::MakeSNeTProfileHistos ( MakeSNeTProfileHistos &&  )
delete

Member Function Documentation

void MakeSNeTProfileHistos::analyze ( art::Event const &  e)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 212 of file MakeSNeTProfileHistos_module.cc.

213 {
214 
215  // TODO:
216  //
217  // 1. implement dividing the TrigEff histos in the endJob function
218  // 2. implement anything that has anything to do with the trigger data products
219  // (including the TrigEff histos, TP splitting histos, etc...???)
220  // 3. figure out a good way to generalize this for multiple MCTruths...
221 
222 
223 
224 
225  // ===== Notes and Assumptions ===== //
226  //
227  // This module assumes the following:
228  //
229  // 1. There is only one neutrino interaction per art::event. This is
230  // important because it assumes that all reco hits go with that one
231  // neutrino.
232 
233 
234 
235 
236 
237  //
238  // Lift out the desired data products:
239  //
240 
241  // Get the hits out of the event
242  auto hits_list = e.getHandle< std::vector< recob::Hit > >(fHitLabel);
243 
244  // Get the MCTruths out of the event
245  auto truths_list = e.getHandle< std::vector< simb::MCTruth > >(fTruthLabel);
246 
247  fNMCTruths->Fill(truths_list->size());
248 
249  // Assert assumption that there is one and only one MCTruth in the event
250  // (note: this needs to be fixed later so that it can handle multiple MCTruths.
251  // This only works for single generated SNe neutrinos for now...)
252  if(truths_list->size() == 1) {
253 
254  // calculate Summed ADC
255  float SummedADC = 0.0; // ...I know... this number is a float?! Who would have guessed that...
256  for(unsigned int i = 0; i < hits_list->size();i++){
257  recob::Hit const& hit = hits_list->at(i);
258  SummedADC += hit.SummedADC();
259  }
260 
261 
262  simb::MCTruth const& truth = truths_list->at(0);
263 
264  fOrigin->Fill(truth.Origin());
265 
266  if (truth.NeutrinoSet()){
267  simb::MCNeutrino const& mc_nu = truth.GetNeutrino();
268  double NuE = 1000.0*mc_nu.Nu().E(); // convert to MeV
269 
270  fEnergySpec->Fill(NuE);
271  fNHitsVsEnergy->Fill(NuE,hits_list->size());
272  fSumADCVsEnergy->Fill(NuE,SummedADC);
273 
274  fCCNC->Fill(mc_nu.CCNC());
275  fNuPDG->Fill(mc_nu.Nu().PdgCode());
276 
277  }
278 
279  }
280  else {
281  std::cerr << "\n\n\nWARNING: Number of MCTruths != 1\n\n\n\n";
282  }
283 
284 
285 
286 }
double E(const int i=0) const
Definition: MCParticle.h:233
int PdgCode() const
Definition: MCParticle.h:212
int CCNC() const
Definition: MCNeutrino.h:148
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:77
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
simb::Origin_t Origin() const
Definition: MCTruth.h:74
const double e
TH1F * fNuPDG
mc_neutrino->Nu()->PdgCode()
Detector simulation of raw signals on wires.
TH1F * fNMCTruths
number of MCTruths per art::event
TH2F * fSumADCVsEnergy
summed ADC for all reco hits per true event energy
float SummedADC() const
The sum of calibrated ADC counts of the hit (0. by default)
Definition: Hit.h:223
std::string fTruthLabel
label for process that made the MCTruth info
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:48
bool NeutrinoSet() const
Definition: MCTruth.h:78
TH2F * fNHitsVsEnergy
number of reconstructed hits per true event energy
Event generator information.
Definition: MCTruth.h:32
TH1F * fEnergySpec
1D event energy spectrum (used only for backgrounds)
Event generator information.
Definition: MCNeutrino.h:18
TH1F * fCCNC
CCNC for each MCTruth.
std::string fHitLabel
label for process that made the reco hits
TH1F * fOrigin
mctruth->Origin()
void MakeSNeTProfileHistos::beginJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 114 of file MakeSNeTProfileHistos_module.cc.

115 {
116 
118 
119  // binning for the energy axes:
120  int NbinsE = 200;
121  float minE = 0.0; // in MeV
122  float maxE = 100.0; // in MeV
123 
124  // binning for the NHits axes:
125  int NbinsNH = 5000;
126  float minNH = 0.0;
127  float maxNH = 5000.0; // a reasonable number determined from running over single nue events
128 
129  // binning for the ADC axes:
130  int NbinsADC = 500;
131  float minADC = 0.0;
132  float maxADC = 500.0e3; // a reasonable number determined from running over single nue events
133 
134 
135 
136  //
137  // Book required histograms for the Time-Profile simulator:
138  //
139 
140  // 2D neutrino energy spectrum as a function of time:
141  // (Note: generated separately by Kate Schol.)
142 
143  // 1D event energy spectrum: (used only for backgrounds)
144  fEnergySpec = tfs->make<TH1F>("EnergySpec",
145  "True energy spectra;true energy [MeV];",
146  NbinsE, minE, maxE);
147 
148  // Overall efficency for triggering as a function of true event energy:
149  fTrigEff_numerator = tfs->make<TH1F>("TrigEff_numerator",
150  "Total number of triggers issued;true energy [MeV];",
151  NbinsE, minE, maxE);
152  fTrigEff = tfs->make<TH1F>("TrigEff",
153  "Trigger efficiency;true energy [MeV];",
154  NbinsE, minE, maxE);
155 
156  // Number of reco hits vs. true event energy:
157  fNHitsVsEnergy = tfs->make<TH2F>("NHitsVsEnergy",
158  "Number of reco hits per true energy;true energy [MeV];NHits",
159  NbinsE, minE, maxE,
160  NbinsNH, minNH, maxNH);
161 
162  // Summed ADC vs. true event energy:
163  fSumADCVsEnergy = tfs->make<TH2F>("SumADCVsEnergy",
164  "Summed ADC over all reco hits per true energy;true energy [MeV];Summed ADC",
165  NbinsE, minE, maxE,
166  NbinsADC, minADC, maxADC);
167 
168  // Probability of splitting a TP into 2 TPs as a function of true event energy:
169  fTPSplitProb = tfs->make<TH1F>("TPSplitProb",
170  "Probability of spliting a TP into 2 or more per true energy;true energy [MeV];",
171  NbinsE, minE, maxE);
172 
173 
174 
175  //
176  // Book a few extra check histos:
177  //
178 
179  // Number of MCTruths per art::event:
180  fNMCTruths = tfs->make<TH1F>("NMCTruths",
181  "Number of MCTruths per art::event;# of MCTruths;",
182  101,-0.5,100.5);
183 
184  // CCNC:
185  fCCNC = tfs->make<TH1F>("CCNC",
186  "CCNC;CCNC;",
187  11,-5.5,5.5);
188 
189  // Origin:
190  fOrigin = tfs->make<TH1F>("Origin",
191  "Neutrino origin;origin type;",
192  16,-5.5,10.5);
193 
194  // Neutrino pdg code:
195  fNuPDG = tfs->make<TH1F>("NuPDG",
196  "Neutrino PDG code;pdg code;",
197  41,-20.5,20.5);
198 
199 }
TH1F * fNuPDG
mc_neutrino->Nu()->PdgCode()
TH1F * fTrigEff_numerator
how many times a trigger was issued per true event energy
TH1F * fTrigEff
trigger efficiency per true event energy (fTrigEff_numerator/fEnergySpec)
TH1F * fNMCTruths
number of MCTruths per art::event
TH2F * fSumADCVsEnergy
summed ADC for all reco hits per true event energy
TH1F * fTPSplitProb
probability of splitting a trigger primitive into two (or more) TPs per true event energy ...
TH2F * fNHitsVsEnergy
number of reconstructed hits per true event energy
TH1F * fEnergySpec
1D event energy spectrum (used only for backgrounds)
TH1F * fCCNC
CCNC for each MCTruth.
TH1F * fOrigin
mctruth->Origin()
void MakeSNeTProfileHistos::endJob ( )
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 204 of file MakeSNeTProfileHistos_module.cc.

205 {
206 
207 }
MakeSNeTProfileHistos& MakeSNeTProfileHistos::operator= ( MakeSNeTProfileHistos const &  )
delete
MakeSNeTProfileHistos& MakeSNeTProfileHistos::operator= ( MakeSNeTProfileHistos &&  )
delete
void MakeSNeTProfileHistos::reconfigure ( fhicl::ParameterSet const &  p)

Definition at line 105 of file MakeSNeTProfileHistos_module.cc.

106 {
107  fTruthLabel = p.get<std::string>("TruthLabel");
108  fHitLabel = p.get<std::string>("HitLabel");
109 }
std::string string
Definition: nybbler.cc:12
p
Definition: test.py:223
std::string fTruthLabel
label for process that made the MCTruth info
std::string fHitLabel
label for process that made the reco hits

Member Data Documentation

TH1F* MakeSNeTProfileHistos::fCCNC
private

CCNC for each MCTruth.

Definition at line 86 of file MakeSNeTProfileHistos_module.cc.

TH1F* MakeSNeTProfileHistos::fEnergySpec
private

1D event energy spectrum (used only for backgrounds)

Definition at line 75 of file MakeSNeTProfileHistos_module.cc.

std::string MakeSNeTProfileHistos::fHitLabel
private

label for process that made the reco hits

Definition at line 70 of file MakeSNeTProfileHistos_module.cc.

TH2F* MakeSNeTProfileHistos::fNHitsVsEnergy
private

number of reconstructed hits per true event energy

Definition at line 78 of file MakeSNeTProfileHistos_module.cc.

TH1F* MakeSNeTProfileHistos::fNMCTruths
private

number of MCTruths per art::event

Definition at line 85 of file MakeSNeTProfileHistos_module.cc.

TH1F* MakeSNeTProfileHistos::fNuPDG
private

mc_neutrino->Nu()->PdgCode()

Definition at line 88 of file MakeSNeTProfileHistos_module.cc.

TH1F* MakeSNeTProfileHistos::fOrigin
private

mctruth->Origin()

Definition at line 87 of file MakeSNeTProfileHistos_module.cc.

TH2F* MakeSNeTProfileHistos::fSumADCVsEnergy
private

summed ADC for all reco hits per true event energy

Definition at line 79 of file MakeSNeTProfileHistos_module.cc.

TH1F* MakeSNeTProfileHistos::fTPSplitProb
private

probability of splitting a trigger primitive into two (or more) TPs per true event energy

Definition at line 80 of file MakeSNeTProfileHistos_module.cc.

TH1F* MakeSNeTProfileHistos::fTrigEff
private

trigger efficiency per true event energy (fTrigEff_numerator/fEnergySpec)

Definition at line 77 of file MakeSNeTProfileHistos_module.cc.

TH1F* MakeSNeTProfileHistos::fTrigEff_numerator
private

how many times a trigger was issued per true event energy

Definition at line 76 of file MakeSNeTProfileHistos_module.cc.

std::string MakeSNeTProfileHistos::fTruthLabel
private

label for process that made the MCTruth info

Definition at line 69 of file MakeSNeTProfileHistos_module.cc.


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