CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BVertexFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoBTag/SecondaryVertex
4 // Class: BVertexFilterT
5 //
13 //
14 // Original Author: Andrea RIZZI
15 // Created: Mon Dec 7 18:02:10 CET 2009
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
35 
38 //
39 // class declaration
40 //
41 
42 template<typename VTX>
44  public:
45  explicit BVertexFilterT(const edm::ParameterSet&);
47 
48  private:
49  virtual bool filter(edm::Event&, const edm::EventSetup&) override;
55 };
56 
57 template<typename VTX>
59  svFilter(params.getParameter<edm::ParameterSet>("vertexFilter")),
60  useVertexKinematicAsJetAxis(params.getParameter<bool>("useVertexKinematicAsJetAxis")),
61  minVertices(params.getParameter<int>("minVertices"))
62 
63 {
64  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
65  token_secondaryVertex = consumes<edm::View<VTX> >(params.getParameter<edm::InputTag>("secondaryVertices"));
66  produces<std::vector<VTX> >();
67 
68 }
69 
70 template<typename VTX>
72 {
73 }
74 
75 template<typename VTX>
76 bool
78 {
79  int count = 0;
81  iEvent.getByToken(token_primaryVertex, pvHandle);
82  edm::Handle<edm::View<VTX> > svHandle;
83  iEvent.getByToken(token_secondaryVertex, svHandle);
84 
85  std::auto_ptr<std::vector<VTX> > recoVertices(new std::vector<VTX>);
86 
87  if(pvHandle->size()!=0) {
88  const reco::Vertex & primary = (*pvHandle.product())[0];
89  const edm::View<VTX> & vertices = *svHandle.product();
90 
91 
92  if(! primary.isFake())
93  {
94  for(typename edm::View<VTX>::const_iterator it=vertices.begin() ; it!=vertices.end() ; ++it)
95  {
96  GlobalVector axis(0,0,0);
97  if(useVertexKinematicAsJetAxis) axis = GlobalVector(it->p4().X(),it->p4().Y(),it->p4().Z());
98  if(svFilter(primary,reco::TemplatedSecondaryVertex<VTX>(primary,*it,axis,true),axis)) {
99  count++;
100  recoVertices->push_back(*it);
101  }
102  }
103  }
104  }
105  iEvent.put(recoVertices);
106 
107  return(count >= minVertices);
108 }
109 
110 
111 // define specific instances of the templated BVertexFilterT class
114 
115 // define plugins
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
BVertexFilterT< reco::Vertex > BVertexFilter
reco::VertexFilter svFilter
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
edm::EDGetTokenT< edm::View< VTX > > token_secondaryVertex
int iEvent
Definition: GenABIO.cc:230
const_iterator begin() const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
bool useVertexKinematicAsJetAxis
BVertexFilterT(const edm::ParameterSet &)
BVertexFilterT< reco::VertexCompositePtrCandidate > CandidateBVertexFilter
bool isFake() const
Definition: Vertex.h:64
T const * product() const
Definition: Handle.h:81
virtual bool filter(edm::Event &, const edm::EventSetup &) override
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
const_iterator end() const
Global3DVector GlobalVector
Definition: GlobalVector.h:10