LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
EventCheater_module.cc
Go to the documentation of this file.
1 //
3 // EventCheater module
4 //
5 // brebel@fnal.gov
6 //
8 #ifndef EVENT_EVENTCHEATER_H
9 #define EVENT_EVENTCHEATER_H
10 #include <string>
11 #include <vector>
12 
13 // ROOT includes
14 
15 // LArSoft includes
24 
25 // Framework includes
30 #include "fhiclcpp/ParameterSet.h"
38 
40 namespace event {
41  class EventCheater : public art::EDProducer {
42  public:
43  explicit EventCheater(fhicl::ParameterSet const& pset);
44  virtual ~EventCheater();
45 
46  void produce(art::Event& evt);
47 
48  void reconfigure(fhicl::ParameterSet const& pset);
49 
50  private:
51 
52  std::string fCheatedVertexLabel;
53  std::string fG4ModuleLabel;
54 
55  };
56 }
57 
58 namespace event{
59 
60  //--------------------------------------------------------------------
62  {
63  this->reconfigure(pset);
64 
65  produces< std::vector<recob::Event> >();
66  produces< art::Assns<recob::Event, recob::Vertex> >();
67  produces< art::Assns<recob::Event, recob::Hit> >();
68  }
69 
70  //--------------------------------------------------------------------
72  {
73  }
74 
75  //--------------------------------------------------------------------
77  {
78  fCheatedVertexLabel = pset.get< std::string >("CheatedVertexLabel", "prong" );
79  fG4ModuleLabel = pset.get< std::string >("G4ModuleLabel", "largeant");
80 
81  return;
82  }
83 
84  //--------------------------------------------------------------------
86  {
87 
89  evt.getView(fG4ModuleLabel, pcol);
90 
92 
93  // make a map of the track id for each sim::Particle to its entry in the
94  // collection of sim::Particles
95  std::map<int, int> trackIDToPColEntry;
96  for(size_t p = 0; p < pcol.vals().size(); ++p) trackIDToPColEntry[pcol.vals().at(p)->TrackId()] = p;
97 
98  // grab the vertices that have been reconstructed
100  evt.getByLabel(fCheatedVertexLabel, vertexcol);
101 
103 
104  // make a vector of them - we aren't writing anything out to a file
105  // so no need for a art::PtrVector here
106  std::vector< art::Ptr<recob::Vertex> > vertices;
107  art::fill_ptr_vector(vertices, vertexcol);
108 
109  // loop over the vertices and figure out which primaries they are associated with
110  std::vector< art::Ptr<recob::Vertex> >::iterator vertexitr = vertices.begin();
111 
112  // make a map of primary product id's to collections of vertices
113  std::map<art::Ptr<simb::MCTruth>, std::vector< art::Ptr<recob::Vertex> > > vertexMap;
114  std::map<art::Ptr<simb::MCTruth>, std::vector< art::Ptr<recob::Vertex> > >::iterator vertexMapItr = vertexMap.begin();
115 
116  // loop over all prongs
117  while( vertexitr != vertices.end() ){
118 
119  size_t pcolEntry = trackIDToPColEntry.find((*vertexitr)->ID())->second;
120  const art::Ptr<simb::MCTruth> primary = fo.at(pcolEntry);
121 
122  vertexMap[primary].push_back(*vertexitr);
123 
124  vertexitr++;
125  }// end loop over vertices
126 
127  std::unique_ptr< std::vector<recob::Event> > eventcol(new std::vector<recob::Event>);
128  std::unique_ptr< art::Assns<recob::Event, recob::Vertex> > evassn(new art::Assns<recob::Event, recob::Vertex>);
129  std::unique_ptr< art::Assns<recob::Event, recob::Hit> > ehassn(new art::Assns<recob::Event, recob::Hit>);
130 
131  // loop over the map and associate all vertex objects with an event
132  for(vertexMapItr = vertexMap.begin(); vertexMapItr != vertexMap.end(); vertexMapItr++){
133 
135 
136  std::vector< art::Ptr<recob::Vertex> > verts( (*vertexMapItr).second );
137 
138  // add an event to the collection.
139  eventcol->push_back(recob::Event(eventcol->size()-1));
140 
141  // associate the event with its vertices
142  util::CreateAssn(*this, evt, *eventcol, verts, *evassn);
143 
144  // get the hits associated with each vertex and associate those with the event
145  for(size_t p = 0; p < ptrvs.size(); ++p){
146  std::vector< art::Ptr<recob::Hit> > hits = fm.at(p);
147  util::CreateAssn(*this, evt, *eventcol, hits, *ehassn);
148  }
149 
150 
151  mf::LogInfo("EventCheater") << "adding event: \n"
152  << eventcol->back()
153  << "\nto collection";
154 
155  } // end loop over the map
156 
157  evt.put(std::move(eventcol));
158  evt.put(std::move(evassn));
159  evt.put(std::move(ehassn));
160 
161  return;
162 
163  } // end produce
164 
165 } // end namespace
166 
167 namespace event{
168 
170 
171 }
172 
173 #endif
std::string fCheatedVertexLabel
label for module creating recob::Vertex objects
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Declaration of signal hit object.
void produce(art::Event &evt)
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
std::size_t getView(std::string const &moduleLabel, std::string const &productInstanceName, std::vector< ELEMENT const * > &result) const
Definition: DataViewImpl.h:474
void hits()
Definition: readHits.C:15
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void reconfigure(fhicl::ParameterSet const &pset)
std::string fG4ModuleLabel
label for module running G4 and making particles, etc
T get(std::string const &key) const
Definition: ParameterSet.h:231
Definition: fwd.h:47
bool CreateAssn(PRODUCER const &prod, art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t indx=UINT_MAX)
Creates a single one-to-one association.
size_type size() const
Definition: PtrVector.h:308
Utility object to perform functions of association.
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:464
Particle list in DetSim contains Monte Carlo particle information.
EventCheater(fhicl::ParameterSet const &pset)
collection_type & vals()
Definition: View.h:34
Event finding and building.