CMS 3D CMS Logo

VertexClassifierByProxy.h
Go to the documentation of this file.
1 
2 #ifndef VertexClassifierByProxy_h
3 #define VertexClassifierByProxy_h
4 
6 
8 
10 template <typename Collection>
12 public:
15 
18  : VertexClassifier(config, std::move(collector)),
19  proxy_(config.getUntrackedParameter<edm::InputTag>("vertexProducer")) {
20  collector.consumes<Association>(proxy_);
21  }
22 
24  void newEvent(edm::Event const &event, edm::EventSetup const &config) override {
25  // Get the association part of the proxy to the collection
26  event.getByLabel(proxy_, proxyHandler_);
27  // Call the previous new event
29  }
30 
34  return *this;
35  }
36 
39  const reco::VertexRefVector *vertexes = nullptr;
40 
41  try {
42  // Get a reference to the vector of associated vertexes
43  vertexes = &(proxyHandler_->find(vertex)->val);
44  } catch (edm::Exception &e) {
45  // If association fails define the vertex as unknown
46  reset();
47  unknownVertex();
48  return *this;
49  }
50 
51  // Evaluate the history for a given index
53 
54  return *this;
55  }
56 
59  const reco::VertexRefVector *vertexes = nullptr;
60 
61  try {
62  // Get a reference to the vector of associated vertexes
63  vertexes = &(proxyHandler_->find(vertex)->val);
64  } catch (edm::Exception &e) {
65  // If association fails define the vertex as unknown
66  reset();
67  unknownVertex();
68  return *this;
69  }
70 
71  // Loop over all the associated vertexes
72  for (std::size_t index = 0; index < vertexes->size(); ++index) {
73  // Copy the last status for all the flags
75 
76  // Evaluate the history for a given index
78 
79  // Combine OR the flag information
80  for (std::size_t i = 0; i < flags_.size(); ++i)
81  flags_[i] = flags_[i] || flags[i];
82  }
83 
84  return *this;
85  }
86 
87 private:
89 
91 };
92 
93 #endif
Get track history and classification by proxy.
Get track history and classify it in function of their .
void reset()
Reset the categories flags.
Definition: config.py:1
std::vector< bool > Flags
Main types associated to the class.
void newEvent(edm::Event const &event, edm::EventSetup const &config) override
Pre-process event information (for accessing reconstraction information).
Flags flags_
Flag containers.
VertexClassifier const & evaluate(reco::VertexBaseRef const &)
Classify the RecoVertex in categories.
VertexClassifierByProxy< Collection > const & evaluate(edm::Ref< Collection > const &vertex, std::size_t index)
Classify any vertexes in categories.
VertexClassifierByProxy(edm::ParameterSet const &config, edm::ConsumesCollector &&collector)
Constructor by ParameterSet.
edm::Handle< Association > proxyHandler_
VertexClassifierByProxy< Collection > const & evaluate(edm::Ref< Collection > const &vertex)
Classify any vertexes in categories.
VertexClassifierByProxy< Collection > const & evaluate(TrackingVertexRef const &vertex)
Classify the TrackingVertex in categories.
edm::AssociationMap< edm::OneToMany< Collection, reco::VertexCollection > > Association
Association type.
HLT enums.
virtual void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
const Flags & flags() const
Returns flags with the category descriptions.
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1