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  const bool m_doFilter;
61  const double m_cutMinPt;
62  const double m_cutMinPtRatio;
63  const double m_maxTrackPt;
64  const double m_maxChi2;
65  const int32_t m_maxNjets;
66  const int32_t m_maxNjetsOutput;
67 
68  const bool m_newMethod;
69 
70  const double m_maxETA;
71  const double m_pvErr_x;
72  const double m_pvErr_y;
73  const double m_pvErr_z;
74 
75 };
76 
77 //
78 // constants, enums and typedefs
79 //
80 
81 //
82 // static data member definitions
83 //
84 
85 //
86 // constructors and destructor
87 //
89  : m_associator ( consumes<reco::JetTracksAssociationCollection>(iConfig.getParameter<edm::InputTag>("jetTracks")) )
90  , m_beamSpot ( consumes<reco::BeamSpot> (iConfig.getParameter<edm::InputTag>("beamSpot") ) )
91  , m_doFilter ( iConfig.getParameter<bool> ("doFilter") )
92  , m_cutMinPt ( iConfig.getParameter<double> ("minPt") )
93  , m_cutMinPtRatio ( iConfig.getParameter<double> ("minPtRatio") )
94  , m_maxTrackPt ( iConfig.getParameter<double> ("maxTrackPt") )
95  , m_maxChi2 ( iConfig.getParameter<double> ("maxChi2") )
96  , m_maxNjets ( iConfig.getParameter<int32_t>("maxNJetsToCheck") )
97  , m_maxNjetsOutput ( iConfig.getParameter<int32_t>("maxNjetsOutput") )
98  , m_newMethod ( iConfig.getParameter<bool> ("newMethod") )
99  , m_maxETA ( iConfig.getParameter<double> ("maxETA") )
100  , m_pvErr_x ( iConfig.getParameter<double> ("pvErr_x") )
101  , m_pvErr_y ( iConfig.getParameter<double> ("pvErr_y") )
102  , m_pvErr_z ( iConfig.getParameter<double> ("pvErr_z") )
103 {
104  //now do what ever initialization is needed
105  produces<std::vector<reco::CaloJet> >();
106  produces<reco::VertexCollection >();
107 }
108 
109 
111 {
112 
113  // do anything here that needs to be done at desctruction time
114  // (e.g. close files, deallocate resources etc.)
115 
116 }
117 
118 
119 //
120 // member functions
121 // m_maxChi2 m_maxTrackPt
122 
123 // ------------ method called on each new Event ------------
124 bool
126 {
127  using namespace edm;
129  iEvent.getByToken(m_associator, jetTracksAssociation);
130  std::auto_ptr<std::vector<reco::CaloJet> > pOut(new std::vector<reco::CaloJet> );
131 
132  bool result=true;
133  int i = 0;
134  float calopt=0;
135  float trkpt=0;
136  //limit to first two jets
137  for (reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin(), et = jetTracksAssociation->end();
138  it != et && i < m_maxNjets; it++, i++) {
139  if(std::abs(it->first->eta()) < m_maxETA)
140  {
141  reco::TrackRefVector tracks = it->second;
142  math::XYZVector jetMomentum = it->first->momentum();
143  math::XYZVector trMomentum;
144  for(reco::TrackRefVector::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack)
145  {
146  const reco::Track& iTrack = **itTrack;
147  if(m_newMethod && iTrack.chi2()>m_maxChi2) continue;
148  trMomentum += iTrack.momentum();
149  if(m_newMethod) trkpt += std::min(m_maxTrackPt,( iTrack.pt()));
150  else trkpt += iTrack.pt();
151  }
152  calopt += jetMomentum.rho();
153  if(trMomentum.rho()/jetMomentum.rho() < m_cutMinPtRatio || trMomentum.rho() < m_cutMinPt)
154  {
155  pOut->push_back(* dynamic_cast<const reco::CaloJet *>(&(*it->first)));
156  }
157  }
158  }
159  iEvent.put(pOut);
160 
163 
165  e(0, 0) = m_pvErr_x * m_pvErr_x;
166  e(1, 1) = m_pvErr_y * m_pvErr_y;
167  e(2, 2) = m_pvErr_z * m_pvErr_z;
168  reco::Vertex::Point p(beamSpot->x0(), beamSpot->y0(), beamSpot->z0());
169  reco::Vertex thePV(p, e, 0, 0, 0);
170  std::auto_ptr<reco::VertexCollection> pOut2(new reco::VertexCollection);
171  pOut2->push_back(thePV);
172  iEvent.put(pOut2);
173 
174  if(m_doFilter) return result;
175  else
176  return true;
177 }
178 
179 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
180 void
183  desc.add<edm::InputTag> ("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
184  desc.add<edm::InputTag> ("jetTracks",edm::InputTag("hltFastPVJetTracksAssociator"));
185  desc.add<double> ("minPtRatio", 0.1);
186  desc.add<double> ("minPt", 0.0);
187  desc.add<bool> ("doFilter", false);
188  desc.add<int> ("maxNJetsToCheck",2);
189  desc.add<int> ("maxNjetsOutput", 2);
190  desc.add<double> ("maxChi2", 20.0);
191  desc.add<double> ("maxTrackPt", 20.0);
192  desc.add<bool> ("newMethod", false); // <---- newMethod
193  desc.add<double> ("maxETA", 2.4 );
194  desc.add<double> ("pvErr_x", 0.0015);
195  desc.add<double> ("pvErr_y", 0.0015);
196  desc.add<double> ("pvErr_z", 1.5 );
197  descriptions.add("jetVertexChecker",desc);
198 }
199 //define this as a plug-in
int i
Definition: DBlmapReader.cc:9
const double m_maxChi2
transient_vector_type::const_iterator const_iterator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
JetTracksAssociation::Container JetTracksAssociationCollection
typedefs for backward compatibility
const double m_maxTrackPt
const bool m_newMethod
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 double m_pvErr_z
const double m_pvErr_y
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:670
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:253
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:248
tuple result
Definition: mps_fire.py:84
int iEvent
Definition: GenABIO.cc:230
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:544
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
const double m_cutMinPt
double pt() const
track transverse momentum
Definition: TrackBase.h:616
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const int32_t m_maxNjets
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
virtual bool filter(edm::Event &, const edm::EventSetup &) override
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const double m_maxETA
const edm::EDGetTokenT< reco::JetTracksAssociationCollection > m_associator
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
JetVertexChecker(const edm::ParameterSet &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const double m_pvErr_x
const edm::EDGetTokenT< reco::BeamSpot > m_beamSpot
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const int32_t m_maxNjetsOutput
const double m_cutMinPtRatio