All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
EventCheater_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // EventCheater module
4 //
5 // brebel@fnal.gov
6 //
7 ////////////////////////////////////////////////////////////////////////
8 #include <string>
9 
10 // ROOT includes
11 
12 // LArSoft includes
17 #include "nusimdata/SimulationBase/MCParticle.h"
18 #include "nusimdata/SimulationBase/MCTruth.h"
19 
20 // Framework includes
21 #include "art/Framework/Core/ModuleMacros.h"
22 #include "canvas/Persistency/Common/FindManyP.h"
23 #include "art/Framework/Principal/Event.h"
24 #include "art/Framework/Principal/View.h"
25 #include "fhiclcpp/ParameterSet.h"
26 #include "art/Framework/Principal/Handle.h"
27 #include "messagefacility/MessageLogger/MessageLogger.h"
28 #include "canvas/Persistency/Common/FindOneP.h"
29 #include "art/Framework/Core/EDProducer.h"
30 
31 ///Event finding and building
32 namespace event {
33  class EventCheater : public art::EDProducer {
34  public:
35  explicit EventCheater(fhicl::ParameterSet const& pset);
36 
37  private:
38 
39  void produce(art::Event& evt);
40 
41  std::string fCheatedVertexLabel; ///< label for module creating recob::Vertex objects
42  std::string fG4ModuleLabel; ///< label for module running G4 and making particles, etc
43 
44  };
45 }
46 
47 namespace event{
48 
49  //--------------------------------------------------------------------
50  EventCheater::EventCheater(fhicl::ParameterSet const& pset)
51  : EDProducer{pset}
52  {
53  fCheatedVertexLabel = pset.get< std::string >("CheatedVertexLabel", "prong" );
54  fG4ModuleLabel = pset.get< std::string >("G4ModuleLabel", "largeant");
55 
56  produces< std::vector<recob::Event> >();
57  produces< art::Assns<recob::Event, recob::Vertex> >();
58  produces< art::Assns<recob::Event, recob::Hit> >();
59  }
60 
61  //--------------------------------------------------------------------
62  void EventCheater::produce(art::Event& evt)
63  {
64 
65  art::View<simb::MCParticle> pcol;
66  evt.getView(fG4ModuleLabel, pcol);
67 
68  art::FindOneP<simb::MCTruth> fo(pcol, evt, fG4ModuleLabel);
69 
70  // make a map of the track id for each sim::Particle to its entry in the
71  // collection of sim::Particles
72  std::map<int, int> trackIDToPColEntry;
73  for(size_t p = 0; p < pcol.vals().size(); ++p) trackIDToPColEntry[pcol.vals().at(p)->TrackId()] = p;
74 
75  // grab the vertices that have been reconstructed
76  art::Handle< std::vector<recob::Vertex> > vertexcol;
77  evt.getByLabel(fCheatedVertexLabel, vertexcol);
78 
79  art::FindManyP<recob::Hit> fm(vertexcol, evt, fCheatedVertexLabel);
80 
81  // make a vector of them - we aren't writing anything out to a file
82  // so no need for a art::PtrVector here
83  std::vector< art::Ptr<recob::Vertex> > vertices;
84  art::fill_ptr_vector(vertices, vertexcol);
85 
86  // loop over the vertices and figure out which primaries they are associated with
87  std::vector< art::Ptr<recob::Vertex> >::iterator vertexitr = vertices.begin();
88 
89  // make a map of primary product id's to collections of vertices
90  std::map<art::Ptr<simb::MCTruth>, std::vector< art::Ptr<recob::Vertex> > > vertexMap;
91  std::map<art::Ptr<simb::MCTruth>, std::vector< art::Ptr<recob::Vertex> > >::iterator vertexMapItr = vertexMap.begin();
92 
93  // loop over all prongs
94  while( vertexitr != vertices.end() ){
95 
96  size_t pcolEntry = trackIDToPColEntry.find((*vertexitr)->ID())->second;
97  const art::Ptr<simb::MCTruth> primary = fo.at(pcolEntry);
98 
99  vertexMap[primary].push_back(*vertexitr);
100 
101  vertexitr++;
102  }// end loop over vertices
103 
104  std::unique_ptr< std::vector<recob::Event> > eventcol(new std::vector<recob::Event>);
105  std::unique_ptr< art::Assns<recob::Event, recob::Vertex> > evassn(new art::Assns<recob::Event, recob::Vertex>);
106  std::unique_ptr< art::Assns<recob::Event, recob::Hit> > ehassn(new art::Assns<recob::Event, recob::Hit>);
107 
108  // loop over the map and associate all vertex objects with an event
109  for(vertexMapItr = vertexMap.begin(); vertexMapItr != vertexMap.end(); vertexMapItr++){
110 
111  art::PtrVector<recob::Vertex> ptrvs;
112 
113  std::vector< art::Ptr<recob::Vertex> > verts( (*vertexMapItr).second );
114 
115  // add an event to the collection.
116  eventcol->push_back(recob::Event(eventcol->size()-1));
117 
118  // associate the event with its vertices
119  util::CreateAssn(*this, evt, *eventcol, verts, *evassn);
120 
121  // get the hits associated with each vertex and associate those with the event
122  for(size_t p = 0; p < ptrvs.size(); ++p){
123  std::vector< art::Ptr<recob::Hit> > hits = fm.at(p);
124  util::CreateAssn(*this, evt, *eventcol, hits, *ehassn);
125  }
126 
127 
128  mf::LogInfo("EventCheater") << "adding event: \n"
129  << eventcol->back()
130  << "\nto collection";
131 
132  } // end loop over the map
133 
134  evt.put(std::move(eventcol));
135  evt.put(std::move(evassn));
136  evt.put(std::move(ehassn));
137 
138  return;
139 
140  } // end produce
141 
142 } // end namespace
143 
144 namespace event{
145 
146  DEFINE_ART_MODULE(EventCheater)
147 
148 }
std::string fCheatedVertexLabel
label for module creating recob::Vertex objects
Declaration of signal hit object.
pdgs p
Definition: selectors.fcl:22
void produce(art::Event &evt)
std::string fG4ModuleLabel
label for module running G4 and making particles, etc
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
Represents the Event implemented as a self balancing binary search tree.
TCEvent evt
Definition: DataStructs.cxx:8
EventCheater(fhicl::ParameterSet const &pset)