All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BasicTrackAnalysis_tool.cc
Go to the documentation of this file.
2 
3 #include "fhiclcpp/ParameterSet.h"
4 #include "art/Utilities/ToolMacros.h"
5 #include "art/Framework/Services/Registry/ServiceHandle.h"
6 #include "art_root_io/TFileService.h"
7 #include "art/Framework/Core/ModuleMacros.h"
8 #include "art_root_io/TFileDirectory.h"
9 #include "canvas/Utilities/InputTag.h"
10 #include "messagefacility/MessageLogger/MessageLogger.h"
11 #include "canvas/Persistency/Common/FindManyP.h"
12 
14 #include "larcore/CoreUtils/ServiceUtil.h" // lar::providerFrom()
15 
17 
18 #include "TH1.h"
19 #include "TH2.h"
20 #include "TProfile.h"
21 #include "TProfile2D.h"
22 
23 #include <cmath>
24 #include <algorithm>
25 
27 {
28  ////////////////////////////////////////////////////////////////////////
29  //
30  // Class: BasicTrackAnalysis
31  // Module Type: producer
32  // File: BasicTrackAnalysis.h
33  //
34  // The intent of this module is to provide methods for
35  // "analyzing" hits on waveforms
36  //
37  // Configuration parameters:
38  //
39  // TruncMeanFraction - the fraction of waveform bins to discard when
40  //
41  // Created by Tracy Usher (usher@slac.stanford.edu) on February 19, 2016
42  //
43  ////////////////////////////////////////////////////////////////////////
44 
45 // The following typedefs will, obviously, be useful
46 using HitPtrVec = std::vector<art::Ptr<recob::Hit>>;
47 using ViewHitMap = std::map<size_t,HitPtrVec>;
48 using TrackViewHitMap = std::map<int,ViewHitMap>;
49 
51 {
52 public:
53  /**
54  * @brief Constructor
55  *
56  * @param pset
57  */
58  explicit BasicTrackAnalysis(fhicl::ParameterSet const & pset);
59 
60  /**
61  * @brief Destructor
62  */
64 
65  // provide for initialization
66  void configure(fhicl::ParameterSet const & pset) override;
67 
68  /**
69  * @brief Interface for initializing the histograms to be filled
70  *
71  * @param TFileService handle to the TFile service
72  * @param string subdirectory to store the hists in
73  */
74  void initializeHists(art::ServiceHandle<art::TFileService>&, const std::string&) override;
75 
76  /**
77  * @brief Interface for method to executve at the end of run processing
78  *
79  * @param int number of events to use for normalization
80  */
81  void endJob(int numEvents) override;
82 
83  /**
84  * @brief Interface for filling histograms
85  */
86  void fillHistograms(const art::Event&) const override;
87 
88 private:
89 
90  // Fcl parameters.
91  std::vector<art::InputTag> fTrackProducerLabelVec; ///< tag for finding the tracks
92  std::string fLocalDirName; ///< Directory name for histograms
93 
94  // Pointers to the histograms we'll create.
95  TH1D* fHitsByWire[3];
96  TH2D* fPulseHVsWidth[3];
97  TProfile* fPulseHVsHitNo[3];
98 
99  // Useful services, keep copies for now (we can update during begin run periods)
100  const geo::GeometryCore* fGeometry; ///< pointer to Geometry service
101 };
102 
103 //----------------------------------------------------------------------------
104 /// Constructor.
105 ///
106 /// Arguments:
107 ///
108 /// pset - Fcl parameters.
109 ///
110 BasicTrackAnalysis::BasicTrackAnalysis(fhicl::ParameterSet const & pset)
111 {
112  fGeometry = lar::providerFrom<geo::Geometry>();
113 
114  configure(pset);
115 
116  // Report.
117  mf::LogInfo("BasicTrackAnalysis") << "BasicTrackAnalysis configured\n";
118 }
119 
120 //----------------------------------------------------------------------------
121 /// Destructor.
122 BasicTrackAnalysis::~BasicTrackAnalysis()
123 {}
124 
125 //----------------------------------------------------------------------------
126 /// Reconfigure method.
127 ///
128 /// Arguments:
129 ///
130 /// pset - Fcl parameter set.
131 ///
132 void BasicTrackAnalysis::configure(fhicl::ParameterSet const & pset)
133 {
134  fTrackProducerLabelVec = pset.get<std::vector<art::InputTag>>("TrackProducerLabel", std::vector<art::InputTag>()={""});
135  fLocalDirName = pset.get<std::string >( "LocalDirName", std::string("wow"));
136 }
137 
138 //----------------------------------------------------------------------------
139 /// Begin job method.
140 void BasicTrackAnalysis::initializeHists(art::ServiceHandle<art::TFileService>& tfs, const std::string& dirName)
141 {
142  // Make a directory for these histograms
143  art::TFileDirectory dir = tfs->mkdir(dirName.c_str());
144 
145  fHitsByWire[0] = dir.make<TH1D>("HitsByWire0", ";Wire #", fGeometry->Nwires(0), 0., fGeometry->Nwires(0));
146  fHitsByWire[1] = dir.make<TH1D>("HitsByWire1", ";Wire #", fGeometry->Nwires(1), 0., fGeometry->Nwires(1));
147  fHitsByWire[2] = dir.make<TH1D>("HitsByWire2", ";Wire #", fGeometry->Nwires(2), 0., fGeometry->Nwires(2));
148 
149  fPulseHVsWidth[0] = dir.make<TH2D>("PHVsWidth0", ";PH;Width", 100, 0., 100., 100, 0., 20.);
150  fPulseHVsWidth[1] = dir.make<TH2D>("PHVsWidth1", ";PH;Width", 100, 0., 100., 100, 0., 20.);
151  fPulseHVsWidth[2] = dir.make<TH2D>("PHVsWidth2", ";PH;Width", 100, 0., 100., 100, 0., 20.);
152 
153  fPulseHVsHitNo[0] = dir.make<TProfile>("PHVsNo0", ";Hit #;PH", 1000, 0., 1000., 0., 100.);
154  fPulseHVsHitNo[1] = dir.make<TProfile>("PHVsNo1", ";Hit #;PH", 1000, 0., 1000., 0., 100.);
155  fPulseHVsHitNo[2] = dir.make<TProfile>("PHVsNo2", ";Hit #;PH", 1000, 0., 1000., 0., 100.);
156 
157  return;
158 }
159 
160 void BasicTrackAnalysis::fillHistograms(const art::Event& event) const
161 {
162  // The game plan for this module is to look at recob::Tracks and objects associated to tracks
163  // To do this we need a valid track collection for those we are hoping to look at
164  for(const auto& trackLabel : fTrackProducerLabelVec)
165  {
166  art::Handle<std::vector<recob::Track> > trackHandle;
167  event.getByLabel(trackLabel, trackHandle);
168 
169  if (trackHandle.isValid())
170  {
171  // Recover the collection of associations between tracks and hits
172  art::FindManyP<recob::Hit> trackHitAssns(trackHandle, event, trackLabel);
173 
174  for(size_t trackIdx = 0; trackIdx < trackHandle->size(); trackIdx++)
175  {
176  art::Ptr<recob::Track> track(trackHandle,trackIdx);
177 
178  }
179  }
180  }
181 
182  return;
183 }
184 
185 // Useful for normalizing histograms
186 void BasicTrackAnalysis::endJob(int numEvents)
187 {
188  // Normalize wire profiles to be hits/event
189  double normFactor(1./numEvents);
190 
191  for(size_t idx = 0; idx < 3; idx++) fHitsByWire[idx]->Scale(normFactor);
192 
193  return;
194 }
195 
196 DEFINE_ART_CLASS_TOOL(BasicTrackAnalysis)
197 }
void configure(fhicl::ParameterSet const &pset) override
Utilities related to art service access.
Declaration of signal hit object.
process_name use argoneut_mc_hitfinder track
Scale(size_t pos, T factor) -> Scale< T >
void endJob(int numEvents) override
Interface for method to executve at the end of run processing.
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.
BasicTrackAnalysis(fhicl::ParameterSet const &pset)
Constructor.
std::string fLocalDirName
Directory name for histograms.
std::map< size_t, HitPtrVec > ViewHitMap
void initializeHists(art::ServiceHandle< art::TFileService > &, const std::string &) override
Interface for initializing the histograms to be filled.
std::vector< art::InputTag > fTrackProducerLabelVec
tag for finding the tracks
std::map< int, ViewHitMap > TrackViewHitMap
Description of geometry of one entire detector.
tuple dir
Definition: dropbox.py:28
std::vector< art::Ptr< recob::Hit >> HitPtrVec
void fillHistograms(const art::Event &) const override
Interface for filling histograms.
art::ServiceHandle< art::TFileService > tfs
art framework interface to geometry description
const geo::GeometryCore * fGeometry
pointer to Geometry service