LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
EventMaker_module.cc
Go to the documentation of this file.
1 #ifndef EventMaker_h
2 #define EventMaker_h
7 
12 
13 namespace event {
14  class EventMaker;
15 }
16 
18 public:
19  explicit EventMaker(fhicl::ParameterSet const &p);
20  virtual ~EventMaker();
21 
22  virtual void produce(art::Event &e);
23 
24  virtual void reconfigure(fhicl::ParameterSet const & p);
25 
26 private:
27 
28  std::string fVertexModuleLabel;
29  double fProximity;
30 
31 };
32 #endif /* EventMaker_h */
33 
34 //--------------------------------------------------------
36 {
37  this->reconfigure(p);
38 
39  produces< std::vector<recob::Event> >();
40  produces< art::Assns<recob::Event, recob::Vertex> >();
41  produces< art::Assns<recob::Event, recob::Hit> >();
42 }
43 
44 //--------------------------------------------------------
46 {
47  // Clean up dynamic memory and other resources here.
48 }
49 
50 //--------------------------------------------------------
52 {
53 
54  // make the unique_ptr of the vector for the recob::Events
55  std::unique_ptr< std::vector<recob::Event> > eventcol(new std::vector<recob::Event>);
56  std::unique_ptr< art::Assns<recob::Event, recob::Vertex> > evassn(new art::Assns<recob::Event, recob::Vertex>);
57  std::unique_ptr< art::Assns<recob::Event, recob::Hit> > ehassn(new art::Assns<recob::Event, recob::Hit>);
58 
59  // first get the recob::Vertex objects out of the event
61  e.getByLabel(fVertexModuleLabel, vtxHandle);
62 
64 
65  // get a collection of art::Ptrs
66  std::vector< art::Ptr<recob::Vertex> > vtxs;
67  art::fill_ptr_vector(vtxs, vtxHandle);
68 
69  // if only 1 vertex in the event, life is easy
70  if(vtxs.size() == 1){
72  vtxvec.push_back(*(vtxs.begin()));
73  recob::Event evt(0);
74  eventcol->push_back(evt);
75 
76  // associate the event with its vertex
77  util::CreateAssn(*this, e, *eventcol, vtxvec, *evassn);
78 
79  // get the hits associated with the vertex and associate those with the event
80  std::vector< art::Ptr<recob::Hit> > hits = fmhv.at(0);
81  util::CreateAssn(*this, e, *eventcol, hits, *ehassn);
82 
83  e.put(std::move(eventcol));
84  e.put(std::move(ehassn));
85  e.put(std::move(evassn));
86 
87  return;
88  }
89 
90  // now the hard part, we have multiple vertex objects
91  // things to consider:
92  // 1. all particles coming from a common vertex location are in the same event
93  // 2. particles coming from decays of another particle, ie the electron from a muon decay
94  // 3. pi-zero decay vertex will be offset from the interaction vertex
95  int evtctr = 0;
96  std::vector< art::Ptr<recob::Vertex> >::iterator itr = vtxs.begin();
97  std::vector< art::Ptr<recob::Vertex> >::iterator itr2 = vtxs.begin();
98 
99  // make a map of the Ptr to the index in the original vector
100  std::map<art::Ptr<recob::Vertex>, size_t> ptrToIdx;
101  for(size_t v = 0; v < vtxs.size(); ++v) ptrToIdx[ vtxs[v] ] = v;
102 
103  while( itr != vtxs.end() ){
104 
106 
107  // get the current vertex object and put it into the vector
108  art::Ptr<recob::Vertex> curvtx = *itr;
109  vtxVec.push_back(curvtx);
110 
111  double curxyz[3] = {0.};
112  curvtx->XYZ(curxyz);
113 
114  // make itr2 the next one in the list, also remove this one from the list
115  // as it is being used
116  itr2 = vtxs.erase(itr);
117  while( itr2 != vtxs.end() ){
118  art::Ptr<recob::Vertex> nexvtx = *itr2;
119 
120  // get the xyz location of the vertex to compare to the current vertex
121  double nextxyz[3] = {0.};
122  nexvtx->XYZ(nextxyz);
123  if( std::sqrt((curxyz[0]-nextxyz[0])*(curxyz[0]-nextxyz[0]) +
124  (curxyz[1]-nextxyz[1])*(curxyz[1]-nextxyz[1]) +
125  (curxyz[2]-nextxyz[2])*(curxyz[2]-nextxyz[2]) ) <= fProximity){
126 
127  // add this one to the vector and remove it from the list as it is being used
128  vtxVec.push_back(nexvtx);
129  itr2 = vtxs.erase(itr2);
130  }// end if vertices are close enough to be considered from the same event
131  else itr2++;
132  }// end inner loop
133 
134  // make an event from these vertex objects and add them to the collection
135  recob::Event evt(++evtctr);
136  eventcol->push_back(evt);
137 
138  // associate the event with its vertices
139  util::CreateAssn(*this, e, *eventcol, vtxVec, *evassn);
140 
141  // get the hits associated with each vertex and associate those with the event
142  for(size_t p = 0; p < vtxVec.size(); ++p){
143  std::vector< art::Ptr<recob::Hit> > hits = fmhv.at( ptrToIdx[ vtxVec[p] ] );
144  util::CreateAssn(*this, e, *eventcol, hits, *ehassn);
145  }
146 
147  // move the initial iterator forward
148  itr++;
149 
150  }
151 
152  // put the collection of events in the art::Event
153  e.put(std::move(eventcol));
154  e.put(std::move(ehassn));
155  e.put(std::move(evassn));
156 
157  return;
158 
159 }
160 
161 //--------------------------------------------------------
163 {
164  fVertexModuleLabel = p.get< std::string >("VertexModuleLabel");
165  fProximity = p.get< double >("Proximity");
166 }
167 
168 
169 
void XYZ(double *xyz) const
Legacy method to access vertex position, preserved to avoid breaking code. Please try to use Vertex::...
Definition: Vertex.cxx:37
double fProximity
how close a vertex needs to be to another to be from the same event
Declaration of signal hit object.
virtual void reconfigure(fhicl::ParameterSet const &p)
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
std::string fVertexModuleLabel
label of the module making the recob::Vertex objects
void hits()
Definition: readHits.C:15
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:441
T get(std::string const &key) const
Definition: ParameterSet.h:231
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
virtual void produce(art::Event &e)
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
Float_t e
Definition: plot.C:34
Event finding and building.
EventMaker(fhicl::ParameterSet const &p)