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  double m_maxTrackPt;
64  double m_maxChi2;
65  int32_t m_maxNjets;
67 
68 
70 
71 };
72 
73 //
74 // constants, enums and typedefs
75 //
76 
77 //
78 // static data member definitions
79 //
80 
81 //
82 // constructors and destructor
83 //
85 {
86  //now do what ever initialization is needed
87  m_beamSpot = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"));
88  m_associator = consumes<reco::JetTracksAssociationCollection>(iConfig.getParameter<edm::InputTag>("jetTracks"));
89  m_doFilter = iConfig.getParameter<bool>("doFilter");
90  m_cutMinPt = iConfig.getParameter<double>("minPt");
91  m_cutMinPtRatio = iConfig.getParameter<double>("minPtRatio");
92  m_maxNjets = iConfig.getParameter<int32_t>("maxNJetsToCheck");
93  m_maxNjetsOutput = iConfig.getParameter<int32_t>("maxNjetsOutput");
94  m_newMethod = iConfig.getParameter<bool>("newMethod");
95  m_maxTrackPt = iConfig.getParameter<double>("maxTrackPt");
96  m_maxChi2 = iConfig.getParameter<double>("maxChi2");
97  produces<std::vector<reco::CaloJet> >();
98  produces<reco::VertexCollection >();
99 }
100 
101 
103 {
104 
105  // do anything here that needs to be done at desctruction time
106  // (e.g. close files, deallocate resources etc.)
107 
108 }
109 
110 
111 //
112 // member functions
113 // m_maxChi2 m_maxTrackPt
114 
115 // ------------ method called on each new Event ------------
116 bool
118 {
119  using namespace edm;
121  iEvent.getByToken(m_associator, jetTracksAssociation);
122  std::auto_ptr<std::vector<reco::CaloJet> > pOut(new std::vector<reco::CaloJet> );
123 
124  bool result=true;
125  int i = 0;
126  float calopt=0;
127  float trkpt=0;
128  //limit to first two jets
129  for(reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin();
130  it != jetTracksAssociation->end() && i < m_maxNjets; it++, i++) {
131  if(std::abs(it->first->eta()) < 2.4)
132  {
133  reco::TrackRefVector tracks = it->second;
134  math::XYZVector jetMomentum = it->first->momentum();
135  math::XYZVector trMomentum;
136  for(reco::TrackRefVector::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack)
137  {
138  const reco::Track& iTrack = **itTrack;
139  if(m_newMethod && iTrack.chi2()>m_maxChi2) continue;
140  trMomentum += iTrack.momentum();
141  if(m_newMethod) trkpt += std::min(m_maxTrackPt,( iTrack.pt()));
142  else trkpt += iTrack.pt();
143  }
144  calopt += jetMomentum.rho();
145  if(trMomentum.rho()/jetMomentum.rho() < m_cutMinPtRatio || trMomentum.rho() < m_cutMinPt)
146  {
147  pOut->push_back(* dynamic_cast<const reco::CaloJet *>(&(*it->first)));
148  }
149  }
150  }
151  iEvent.put(pOut);
152 
155 
157  e(0, 0) = 0.0015 * 0.0015;
158  e(1, 1) = 0.0015 * 0.0015;
159  e(2, 2) = 1.5 * 1.5;
160  reco::Vertex::Point p(beamSpot->x0(), beamSpot->y0(), beamSpot->z0());
161  reco::Vertex thePV(p, e, 0, 0, 0);
162  std::auto_ptr<reco::VertexCollection> pOut2(new reco::VertexCollection);
163  pOut2->push_back(thePV);
164  iEvent.put(pOut2);
165 
166  if(m_doFilter) return result;
167  else
168  return true;
169 }
170 
171 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
172 void
175  desc.add<edm::InputTag> ("beamSpot",edm::InputTag("hltOnlineBeamSpot"));
176  desc.add<edm::InputTag> ("jetTracks",edm::InputTag("hltFastPVJetTracksAssociator"));
177  desc.add<double> ("minPtRatio",0.1);
178  desc.add<double> ("minPt",0.0);
179  desc.add<bool> ("doFilter",false);
180  desc.add<int> ("maxNJetsToCheck",2);
181  desc.add<int> ("maxNjetsOutput",2);
182  desc.add<double> ("maxChi2",20.0);
183  desc.add<double> ("maxTrackPt",20.0);
184  desc.add<bool> ("newMethod",false); // <---- newMethod
185  descriptions.add("jetVertexChecker",desc);
186 }
187 //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:446
#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 Vector & momentum() const
track momentum vector
Definition: TrackBase.h:723
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
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:597
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
double pt() const
track transverse momentum
Definition: TrackBase.h:669
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::JetTracksAssociationCollection > m_associator