All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
sbndcode/sbndcode/Calibration/HitEfficiencyAna_module.cc
Go to the documentation of this file.
1 // HitEfficiencyAna_module.cc
2 //
3 //Shamelessly stolen from ICARUS by A. Scarff
4 //
5 // A basic "skeleton" to read in art::Event records from a file,
6 // access their information, and do something with them.
7 // See
8 // <https://cdcvs.fnal.gov/redmine/projects/larsoftsvn/wiki/Using_the_Framework>
9 // for a description of the ART classes used here.
10 // Almost everything you see in the code below may have to be changed
11 // by you to suit your task. The example task is to make histograms
12 // and n-tuples related to dE/dx of particle tracks in the detector.
13 // As you try to understand why things are done a certain way in this
14 // example ("What's all this stuff about 'auto const&'?"), it will help
15 // to read ADDITIONAL_NOTES.txt in the same directory as this file.
16 #ifndef HitEfficiencyAna_module
17 #define HitEfficiencyAna_module
18 // LArSoft includes
20 #include "larcore/CoreUtils/ServiceUtil.h" // lar::providerFrom()
26 //#include "cetlib/search_path.h"
27 #include "cetlib/cpu_timer.h"
31 // Framework includes
32 #include "art/Framework/Core/EDAnalyzer.h"
33 #include "art/Framework/Principal/Event.h"
34 #include "art/Framework/Principal/Handle.h"
35 #include "art/Framework/Principal/View.h"
36 #include "art/Framework/Services/Registry/ServiceHandle.h"
37 #include "art_root_io/TFileService.h"
38 #include "art/Framework/Core/ModuleMacros.h"
39 #include "art/Utilities/make_tool.h"
40 #include "canvas/Utilities/InputTag.h"
41 #include "canvas/Persistency/Common/FindManyP.h"
42 #include "messagefacility/MessageLogger/MessageLogger.h"
43 #include "fhiclcpp/ParameterSet.h"
44 #include "cetlib_except/exception.h"
45 #include "nusimdata/SimulationBase/MCParticle.h"
47 #include "TTree.h"
48 // C++ Includes
49 #include <map>
50 #include <vector>
51 #include <tuple>
52 #include <algorithm>
53 #include <iostream>
54 #include <string>
55 #include <cmath>
56 #include <iostream>
57 #include <fstream>
58 namespace HitEfficiencyAna
59 {
60  //-----------------------------------------------------------------------
61  //-----------------------------------------------------------------------
62  // class definition
63  class HitEfficiencyAna : public art::EDAnalyzer
64  {
65  public:
66  // Standard constructor and destructor for an ART module.
67  explicit HitEfficiencyAna(fhicl::ParameterSet const& pset);
68  virtual ~HitEfficiencyAna();
69  // This method is called once, at the start of the job. In this
70  // example, it will define the histograms and n-tuples we'll write.
71  void beginJob();
72  void endJob();
73  // This method is called once, at the start of each run. It's a
74  // good place to read databases or files that may have
75  // run-dependent information.
76  void beginRun(const art::Run& run);
77  // This method reads in any parameters from the .fcl files. This
78  // method is called 'reconfigure' because it might be called in the
79  // middle of a job; e.g., if the user changes parameter values in an
80  // interactive event display.
81  void reconfigure(fhicl::ParameterSet const& pset);
82  // The analysis routine, called once per event.
83  void analyze (const art::Event& evt);
84  private:
85  // The variables that will go into the n-tuple.
86  int fEvent;
87  int fRun;
88  int fSubRun;
89  int fNumEvents;
90 
91  // Introduce a temporary event counter to handle the current issue where all input
92  // files start with the same run/subrun/event sequence
93  int fNumCalls;
94 
95  // Keep track of the hit histogramming tools here
96  std::vector<std::unique_ptr<IHitEfficiencyHistogramTool>> fHitHistogramToolVec;
97 
98  // And our tuple
99  mutable TTree* fTree;
100  // Other variables that will be shared between different methods.
101  const geo::GeometryCore* fGeometry; // pointer to Geometry service
102  //const lariov::DetPedestalProvider& fPedestalRetrievalAlg; ///< Keep track of an instance to the pedestal retrieval alg
103  }; // class HitEfficiencyAna
104  //-----------------------------------------------------------------------
105  //-----------------------------------------------------------------------
106  // class implementation
107  //-----------------------------------------------------------------------
108  // Constructor
109  HitEfficiencyAna::HitEfficiencyAna(fhicl::ParameterSet const& parameterSet)
110  : EDAnalyzer(parameterSet),
111  fNumCalls(0)
112  // fPedestalRetrievalAlg(*lar::providerFrom<lariov::DetPedestalService>())
113  {
114  fGeometry = lar::providerFrom<geo::Geometry>();
115  // Read in the parameters from the .fcl file.
116  this->reconfigure(parameterSet);
117  }
118  //-----------------------------------------------------------------------
119  // Destructor
120  HitEfficiencyAna::~HitEfficiencyAna()
121  {}
122  //-----------------------------------------------------------------------
123  void HitEfficiencyAna::beginJob()
124  {
125  // Get the detector length, to determine the maximum bin edge of one
126  // of the histograms.
127  //double detectorLength = fGeometry->DetLength();
128  // Access ART's TFileService, which will handle creating and writing
129  // histograms and n-tuples for us.
130  art::ServiceHandle<art::TFileService> tfs;
131  std::cout << "beginjob running" << std::endl;
132  // First task is to create a TTree..
133  fTree = tfs->makeAndRegister<TTree>("HitEffic_t","Hit Efficiency Tuple");
134  fTree->Branch("Run", &fRun, "Run/I");
135  fTree->Branch("SubRun", &fSubRun, "SubRun/I");
136  fTree->Branch("Event", &fEvent, "Event/I");
137  // The arguments to 'make<whatever>' are the same as those passed
138  // to the 'whatever' constructor, provided 'whatever' is a ROOT
139  // class that TFileService recognizes.
140  for (auto& hitHistTool : fHitHistogramToolVec)
141  {
142  hitHistTool->initializeHists(tfs, "HitEffic");
143  hitHistTool->initializeTuple(fTree);
144  }
145  // zero out the event counter
146  fNumEvents = 0;
147  }
148  //-----------------------------------------------------------------------
149  void HitEfficiencyAna::beginRun(const art::Run& /*run*/)
150  {
151  // How to convert from number of electrons to GeV. The ultimate
152  // source of this conversion factor is
153  // ${LARSIM_DIR}/include/SimpleTypesAndConstants/PhysicalConstants.h.
154  // art::ServiceHandle<sim::LArG4Parameters> larParameters;
155  // fElectronsToGeV = 1./larParameters->GeVToElectrons();
156  }
157  //-----------------------------------------------------------------------
158  void HitEfficiencyAna::reconfigure(fhicl::ParameterSet const& p)
159  {
160  // Read parameters from the .fcl file. The names in the arguments
161  // to p.get<TYPE> must match names in the .fcl file.
162  // Implement the tools for handling the responses
163  const std::vector<fhicl::ParameterSet>& hitHistogramToolVec = p.get<std::vector<fhicl::ParameterSet>>("HitEfficiencyHistogramToolList");
164 
165  for(auto& hitHistogramTool : hitHistogramToolVec)
166  fHitHistogramToolVec.push_back(art::make_tool<IHitEfficiencyHistogramTool>(hitHistogramTool));
167  std::cout << "reconfigure running" << std::endl;
168  return;
169  }
170  //-----------------------------------------------------------------------
171  void HitEfficiencyAna::analyze(const art::Event& event)
172  {
173  // Start by fetching some basic event information for our n-tuple.
174  fEvent = event.id().event();
175  fRun = ++fNumCalls; //event.run();
176  fSubRun = event.subRun();
177  fNumEvents++;
178 
179  // Make a pass through all hits to make contrasting plots
180  std::cout << "-- Run: " << fRun << ", SubRun: " << fSubRun << ", Event: " << fEvent << " -------" << std::endl;
181  std::cout << "HELLO!" << std::endl;
182 
183  for(auto& hitHistTool : fHitHistogramToolVec) hitHistTool->fillHistograms(event);
184 
185  fTree->Fill();
186  return;
187  }
188  void HitEfficiencyAna::endJob()
189  {
190  // Make a call to normalize histograms if so desired
191  for(auto& hitHistTool : fHitHistogramToolVec) hitHistTool->endJob(fNumEvents);
192  std::cout << "endjob running" << std::endl;
193  return;
194  }
195  // This macro has to be defined for this module to be invoked from a
196  // .fcl file; see HitEfficiencyAna.fcl for more information.
197  DEFINE_ART_MODULE(HitEfficiencyAna)
198 } // namespace HitEfficiencyAna
199 #endif // HitEfficiencyAna_module
process_name opflashCryo1 flashfilter analyze
Utilities related to art service access.
Declaration of signal hit object.
pdgs p
Definition: selectors.fcl:22
for pfile in ack l reconfigure(.*) override"` do echo "checking $
Description of geometry of one entire detector.
Definition of data types for geometry description.
std::vector< std::unique_ptr< IHitEfficiencyHistogramTool > > fHitHistogramToolVec
object containing MC truth information necessary for making RawDigits and doing back tracking ...
art::ServiceHandle< art::TFileService > tfs
TCEvent evt
Definition: DataStructs.cxx:8
art framework interface to geometry description
BEGIN_PROLOG could also be cout