LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
HitAnaModule_module.cc
Go to the documentation of this file.
1 // Class: HitAnaModule
3 // Module Type: analyzer
4 // File: HitAnaModule_module.cc
5 //
6 // Generated at Sun Jun 22 08:40:08 2014 by Wesley Ketchum using artmod
7 // from cetpkgsupport v1_05_04.
9 
28 #include "fhiclcpp/ParameterSet.h"
29 
30 #include "art_root_io/TFileService.h"
31 
32 #include <string>
33 
39 
40 #include "TTree.h"
41 
42 namespace hit {
43  class HitAnaModule;
44 }
45 
47 public:
48  explicit HitAnaModule(fhicl::ParameterSet const& p);
49 
50 private:
51  void analyze(art::Event const& e) override;
52 
53  void beginJob() override;
54 
56 
57  void createAssocVector(HitWireAssns_t const&, std::vector<std::vector<int>>&);
58 
59  void createAssocVector(std::vector<recob::Wire> const&,
60  std::vector<recob::Hit> const&,
61  std::vector<std::vector<int>>&);
62 
63  void createMCAssocVector(std::vector<recob::Wire> const&,
64  std::vector<sim::MCHitCollection> const&,
65  std::vector<std::vector<int>>&);
66 
67  // Declare member data here.
68  std::vector<std::string> fHitModuleLabels;
69  std::string fMCHitModuleLabel;
70  std::string fWireModuleLabel;
71 
72  TTree* wireDataTree;
73  std::vector<TTree*> hitDataTree;
74 
76 };
77 
79 // More initializers here.
80 {
81  fHitModuleLabels = p.get<std::vector<std::string>>("HitModuleLabels");
82  fWireModuleLabel = p.get<std::string>("WireModuleLabel");
83  fMCHitModuleLabel = p.get<std::string>("MCHitModuleLabel");
84 }
85 
87  std::vector<std::vector<int>>& WireHitAssocVector)
88 {
89  // WireHitAssocVector: for each wire, indices of all the hits associated to it
90 
91  // the iteration to art::Assns<Hit, Wire> points to a art::Ptr pair (assn_t)
92  // with a hit as first element ("left") and a wire as the second one ("right")
93  for (HitWireAssns_t::assn_t const& assn : HitToWire)
94  WireHitAssocVector.at(assn.second.key()).push_back(assn.first.key());
95 }
96 
97 void hit::HitAnaModule::createAssocVector(std::vector<recob::Wire> const& wireVector,
98  std::vector<recob::Hit> const& hitVector,
99  std::vector<std::vector<int>>& WireHitAssocVector)
100 {
101  WireHitAssocVector.resize(wireVector.size());
102  for (size_t iwire = 0; iwire < wireVector.size(); iwire++) {
103  for (size_t ihit = 0; ihit < hitVector.size(); ihit++) {
104  if (hitVector[ihit].Channel() == wireVector[iwire].Channel())
105  WireHitAssocVector[iwire].push_back(ihit);
106  }
107  }
108 }
109 
110 void hit::HitAnaModule::createMCAssocVector(std::vector<recob::Wire> const& wireVector,
111  std::vector<sim::MCHitCollection> const& mcHitVector,
112  std::vector<std::vector<int>>& WireMCHitAssocVector)
113 {
114 
115  WireMCHitAssocVector.clear();
116  WireMCHitAssocVector.resize(wireVector.size());
117 
118  //first, store all the MCHitCollection indices in a map keyed on channel
119  //then, loop through wires, and lookup mchitcollections based on the wire's channel
120 
121  std::map<unsigned int, std::vector<int>> mcHitIndicesByChannel;
122  for (unsigned int icol = 0; icol < mcHitVector.size(); icol++)
123  mcHitIndicesByChannel[mcHitVector[icol].Channel()].push_back(icol);
124 
125  for (unsigned int iwire = 0; iwire < wireVector.size(); iwire++)
126  WireMCHitAssocVector[iwire].insert(WireMCHitAssocVector[iwire].end(),
127  mcHitIndicesByChannel[wireVector[iwire].Channel()].begin(),
128  mcHitIndicesByChannel[wireVector[iwire].Channel()].end());
129 }
130 
132 {
133  //get event and run numbers
134  unsigned int eventNumber = e.id().event();
135  unsigned int runNumber = e.run();
136 
138 
139  // get the data
140  auto const& wireVector = *e.getValidHandle<std::vector<recob::Wire>>(fWireModuleLabel);
141  auto const& mcHitVector = *e.getValidHandle<std::vector<sim::MCHitCollection>>(fMCHitModuleLabel);
142 
143  // make the association vector. First index is wire index, second is
144  // mcHitCollection index
145  std::vector<std::vector<int>> WireMCHitAssocVector;
146  createMCAssocVector(wireVector, mcHitVector, WireMCHitAssocVector);
147 
148  //get the hit data
149  size_t nHitModules = fHitModuleLabels.size();
150  std::vector<art::Handle<std::vector<recob::Hit>>> hitHandles(nHitModules);
151  // for each hit module output (first index), for each wire (second index)
152  // the list of hits associated with that wire is stored
153  std::vector<std::vector<std::vector<int>>> WireHitAssocVectors(nHitModules);
154  for (size_t iter = 0; iter < nHitModules; iter++) {
155 
156  e.getByLabel(fHitModuleLabels[iter], hitHandles[iter]);
157 
158  //create association vectors by hand for now
159  //art::ValidHandle<HitWireAssns_t> HitToWireAssns
160  //= e.getValidHandle<HitWireAssns_t>(fHitModuleLabels[iter]);
161  //WireHitAssocVectors[iter].resize(wireVector.size());
162  //createAssocVector(*HitToWireAssns,WireHitAssocVectors[iter]);
163 
164  WireHitAssocVectors[iter].resize(wireVector.size());
165  art::Handle<HitWireAssns_t> HitToWireAssns;
166  if (e.getByLabel(fHitModuleLabels[iter], HitToWireAssns))
167  createAssocVector(*HitToWireAssns, WireHitAssocVectors[iter]);
168  else
169  createAssocVector(wireVector, *(hitHandles[iter]), WireHitAssocVectors[iter]);
170 
171  //load in this hit/assoc pair
173  *(hitHandles[iter]), WireHitAssocVectors[iter], fHitModuleLabels[iter]);
174  }
175 
176  // get time service
177  auto const clock_data = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(e);
178 
179  // run the analzyer alg
181  wireVector, mcHitVector, WireMCHitAssocVector, clock_data, eventNumber, runNumber);
182 }
183 
185 {
187  wireDataTree = tfs->make<TTree>("wireDataTree", "WireDataTree");
189 
190  // The below creates a Tree with one branch - a recob::Hit branch - for each
191  // Hit module specified in the fcl file. So, don't run this module once per Hit
192  // Finder as one normally would. Just run it once and specify all HitFinders.
193  // This was the design for the WireDataTree; we follow it here for the
194  // Hit trees.
195  for (auto const& label : fHitModuleLabels) {
196  std::string firstArg("hitData_");
197  firstArg += label;
198  std::string secArg("HitDataTree_");
199  secArg += label;
200  TTree* intermediateTree = tfs->make<TTree>(firstArg.c_str(), secArg.c_str());
201  hitDataTree.push_back(intermediateTree);
202  }
204 }
205 
void LoadHitAssocPair(std::vector< recob::Hit > const &, std::vector< std::vector< int >> const &, std::string const &)
Definition: HitAnaAlg.cxx:92
void createMCAssocVector(std::vector< recob::Wire > const &, std::vector< sim::MCHitCollection > const &, std::vector< std::vector< int >> &)
void SetHitDataTree(std::vector< TTree * > &trees)
Definition: HitAnaAlg.cxx:32
std::vector< std::string > fHitModuleLabels
void beginJob() override
std::string fWireModuleLabel
std::vector< TTree * > hitDataTree
Declaration of signal hit object.
std::string fMCHitModuleLabel
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
void SetWireDataTree(TTree *)
Definition: HitAnaAlg.cxx:26
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
T get(std::string const &key) const
Definition: ParameterSet.h:314
details::FindAllP< recob::Hit, recob::Wire > HitToWire
Query object connecting a hit to a wire.
Definition: HitUtils.h:54
void AnalyzeWires(std::vector< recob::Wire > const &, std::vector< sim::MCHitCollection > const &, std::vector< std::vector< int >> const &, detinfo::DetectorClocksData const &, unsigned int, unsigned int)
Definition: HitAnaAlg.cxx:103
Detector simulation of raw signals on wires.
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
HitAnaModule(fhicl::ParameterSet const &p)
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
void createAssocVector(HitWireAssns_t const &, std::vector< std::vector< int >> &)
EventNumber_t event() const
Definition: EventID.h:116
Declaration of basic channel signal object.
void analyze(art::Event const &e) override
Float_t e
Definition: plot.C:35
RunNumber_t run() const
Definition: Event.cc:29
void ClearHitModules()
Definition: HitAnaAlg.cxx:85
EventID id() const
Definition: Event.cc:23