CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetVertexChecker.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: JetVertexChecker
4 // Class: JetVertexChecker
5 //
13 //
14 // Original Author: Andrea RIZZI
15 // Created: Mon Jan 16 11:19:48 CET 2012
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
34 
37 
40 
42 
43 //
44 // class declaration
45 //
46 
48  public:
49  explicit JetVertexChecker(const edm::ParameterSet&);
51 
52  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
53 
54  private:
55  virtual bool filter(edm::Event&, const edm::EventSetup&) override;
56 
57  // ----------member data ---------------------------
60  bool m_doFilter;
61  double m_cutMinPt;
62  double m_cutMinPtRatio;
63  int32_t m_maxNjets;
64 
65 };
66 
67 //
68 // constants, enums and typedefs
69 //
70 
71 //
72 // static data member definitions
73 //
74 
75 //
76 // constructors and destructor
77 //
79 {
80  //now do what ever initialization is needed
81  m_beamSpot = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"));
82  m_associator = consumes<reco::JetTracksAssociationCollection>(iConfig.getParameter<edm::InputTag>("jetTracks"));
83  m_doFilter = iConfig.getParameter<bool>("doFilter");
84  m_cutMinPt = iConfig.getParameter<double>("minPt");
85  m_cutMinPtRatio = iConfig.getParameter<double>("minPtRatio");
86  m_maxNjets = iConfig.getParameter<int32_t>("maxNJetsToCheck");
87  produces<std::vector<reco::CaloJet> >();
88  produces<reco::VertexCollection >();
89 }
90 
91 
93 {
94 
95  // do anything here that needs to be done at desctruction time
96  // (e.g. close files, deallocate resources etc.)
97 
98 }
99 
100 
101 //
102 // member functions
103 //
104 
105 // ------------ method called on each new Event ------------
106 bool
108 {
109  using namespace edm;
111  iEvent.getByToken(m_associator, jetTracksAssociation);
112  std::auto_ptr<std::vector<reco::CaloJet> > pOut(new std::vector<reco::CaloJet> );
113 
114  bool result=true;
115  int i = 0;
116  //limit to first two jets
117  for(reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin();
118  it != jetTracksAssociation->end() && i < m_maxNjets; it++, i++) {
119  if(fabs(it->first->eta()) < 2.4)
120  {
121  reco::TrackRefVector tracks = it->second;
122  math::XYZVector jetMomentum = it->first->momentum();
123  math::XYZVector trMomentum;
124  for(reco::TrackRefVector::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack)
125  {
126  trMomentum += (*itTrack)->momentum();
127  }
128  if(trMomentum.rho()/jetMomentum.rho() < m_cutMinPtRatio || trMomentum.rho() < m_cutMinPt)
129  {
130 // std::cout << "bad jet " << it->first->pt() << std::endl;
131  pOut->push_back(* dynamic_cast<const reco::CaloJet *>(&(*it->first)));
132  result=false;
133  }
134  }
135  }
136 
137  iEvent.put(pOut);
138 
141 
143  e(0, 0) = 0.0015 * 0.0015;
144  e(1, 1) = 0.0015 * 0.0015;
145  e(2, 2) = 1.5 * 1.5;
146  reco::Vertex::Point p(beamSpot->x0(), beamSpot->y0(), beamSpot->z0());
147  reco::Vertex thePV(p, e, 0, 0, 0);
148  std::auto_ptr<reco::VertexCollection> pOut2(new reco::VertexCollection);
149  pOut2->push_back(thePV);
150  iEvent.put(pOut2);
151 // std::cout << " filter " << result << std::endl;
152  if(m_doFilter) return result;
153  else
154  return true;
155 }
156 
157 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
158 void
160  //The following says we do not know what parameters are allowed so do no validation
161  // Please change this to state exactly what you do use, even if it is no parameters
163  desc.setUnknown();
164  descriptions.addDefault(desc);
165 }
166 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
transient_vector_type::const_iterator const_iterator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
tuple result
Definition: query.py:137
edm::EDGetTokenT< reco::BeamSpot > m_beamSpot
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
virtual bool filter(edm::Event &, const edm::EventSetup &) override
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
JetVertexChecker(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::JetTracksAssociationCollection > m_associator