CMS 3D CMS Logo

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 // system include files
20 #include <memory>
21 
22 // user include files
25 
28 
33 
36 
39 
41 
42 //
43 // class declaration
44 //
45 
47 public:
48  explicit JetVertexChecker(const edm::ParameterSet&);
49  ~JetVertexChecker() override;
50 
51  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
52 
53 private:
54  bool filter(edm::Event&, const edm::EventSetup&) override;
55 
56  // ----------member data ---------------------------
59  const bool m_doFilter;
60  const double m_cutMinPt;
61  const double m_cutMinPtRatio;
62  const double m_maxTrackPt;
63  const double m_maxChi2;
64  const int32_t m_maxNjets;
65  const int32_t m_maxNjetsOutput;
66 
67  const bool m_newMethod;
68 
69  const double m_maxETA;
70  const double m_pvErr_x;
71  const double m_pvErr_y;
72  const double m_pvErr_z;
73 };
74 
75 //
76 // constants, enums and typedefs
77 //
78 
79 //
80 // static data member definitions
81 //
82 
83 //
84 // constructors and destructor
85 //
87  : m_associator(consumes<reco::JetTracksAssociationCollection>(iConfig.getParameter<edm::InputTag>("jetTracks"))),
88  m_beamSpot(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))),
89  m_doFilter(iConfig.getParameter<bool>("doFilter")),
90  m_cutMinPt(iConfig.getParameter<double>("minPt")),
91  m_cutMinPtRatio(iConfig.getParameter<double>("minPtRatio")),
92  m_maxTrackPt(iConfig.getParameter<double>("maxTrackPt")),
93  m_maxChi2(iConfig.getParameter<double>("maxChi2")),
94  m_maxNjets(iConfig.getParameter<int32_t>("maxNJetsToCheck")),
95  m_maxNjetsOutput(iConfig.getParameter<int32_t>("maxNjetsOutput")),
96  m_newMethod(iConfig.getParameter<bool>("newMethod")),
97  m_maxETA(iConfig.getParameter<double>("maxETA")),
98  m_pvErr_x(iConfig.getParameter<double>("pvErr_x")),
99  m_pvErr_y(iConfig.getParameter<double>("pvErr_y")),
100  m_pvErr_z(iConfig.getParameter<double>("pvErr_z")) {
101  //now do what ever initialization is needed
102  produces<std::vector<reco::CaloJet>>();
103  produces<reco::VertexCollection>();
104 }
105 
107  // do anything here that needs to be done at desctruction time
108  // (e.g. close files, deallocate resources etc.)
109 }
110 
111 //
112 // member functions
113 // m_maxChi2 m_maxTrackPt
114 
115 // ------------ method called on each new Event ------------
117  using namespace edm;
119  iEvent.getByToken(m_associator, jetTracksAssociation);
120  auto pOut = std::make_unique<std::vector<reco::CaloJet>>();
121 
122  bool result = true;
123  int i = 0;
124  float calopt = 0;
125  float trkpt = 0;
126  //limit to first two jets
127  for (reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin(),
128  et = jetTracksAssociation->end();
129  it != et && i < m_maxNjets;
130  it++, i++) {
131  if (std::abs(it->first->eta()) < m_maxETA) {
132  reco::TrackRefVector tracks = it->second;
133  math::XYZVector jetMomentum = it->first->momentum();
134  math::XYZVector trMomentum;
135  for (reco::TrackRefVector::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack) {
136  const reco::Track& iTrack = **itTrack;
137  if (m_newMethod && iTrack.chi2() > m_maxChi2)
138  continue;
139  trMomentum += iTrack.momentum();
140  if (m_newMethod)
141  trkpt += std::min(m_maxTrackPt, (iTrack.pt()));
142  else
143  trkpt += iTrack.pt();
144  }
145  calopt += jetMomentum.rho();
146  if (trMomentum.rho() / jetMomentum.rho() < m_cutMinPtRatio || trMomentum.rho() < m_cutMinPt) {
147  pOut->push_back(*dynamic_cast<const reco::CaloJet*>(&(*it->first)));
148  }
149  }
150  }
151  iEvent.put(std::move(pOut));
152 
154  iEvent.getByToken(m_beamSpot, beamSpot);
155 
157  e(0, 0) = m_pvErr_x * m_pvErr_x;
158  e(1, 1) = m_pvErr_y * m_pvErr_y;
159  e(2, 2) = m_pvErr_z * m_pvErr_z;
160  reco::Vertex::Point p(beamSpot->x0(), beamSpot->y0(), beamSpot->z0());
161  reco::Vertex thePV(p, e, 0, 0, 0);
162  auto pOut2 = std::make_unique<reco::VertexCollection>();
163  pOut2->push_back(thePV);
164  iEvent.put(std::move(pOut2));
165 
166  if (m_doFilter)
167  return result;
168  else
169  return true;
170 }
171 
172 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
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  desc.add<double>("maxETA", 2.4);
186  desc.add<double>("pvErr_x", 0.0015);
187  desc.add<double>("pvErr_y", 0.0015);
188  desc.add<double>("pvErr_z", 1.5);
189  descriptions.add("jetVertexChecker", desc);
190 }
191 //define this as a plug-in
JetTracksAssociation.h
CaloJetCollection.h
CaloJet.h
electrons_cff.bool
bool
Definition: electrons_cff.py:366
mps_fire.i
i
Definition: mps_fire.py:428
pwdgSkimBPark_cfi.beamSpot
beamSpot
Definition: pwdgSkimBPark_cfi.py:5
JetVertexChecker::m_pvErr_y
const double m_pvErr_y
Definition: JetVertexChecker.cc:71
align::BeamSpot
Definition: StructureType.h:95
JetVertexChecker::m_beamSpot
const edm::EDGetTokenT< reco::BeamSpot > m_beamSpot
Definition: JetVertexChecker.cc:58
min
T min(T a, T b)
Definition: MathUtil.h:58
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::Vertex::Error
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
JetVertexChecker::m_newMethod
const bool m_newMethod
Definition: JetVertexChecker.cc:67
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
L1TrackObjectNtupleMaker_cfg.pOut
pOut
Definition: L1TrackObjectNtupleMaker_cfg.py:172
edm::RefVector< TrackCollection >
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
edm::Handle
Definition: AssociativeIterator.h:50
JetVertexChecker::m_maxNjetsOutput
const int32_t m_maxNjetsOutput
Definition: JetVertexChecker.cc:65
reco::TrackBase::pt
double pt() const
track transverse momentum
Definition: TrackBase.h:637
JetVertexChecker::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: JetVertexChecker.cc:173
MakerMacros.h
JetVertexChecker
Definition: JetVertexChecker.cc:46
Track.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
JetVertexChecker::m_maxNjets
const int32_t m_maxNjets
Definition: JetVertexChecker.cc:64
BeamSpot.h
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
edm::AssociationVector::begin
const_iterator begin() const
Definition: AssociationVector.h:108
reco::Track
Definition: Track.h:27
edm::AssociationVector::end
const_iterator end() const
Definition: AssociationVector.h:109
JetVertexChecker::m_maxTrackPt
const double m_maxTrackPt
Definition: JetVertexChecker.cc:62
JetVertexChecker::m_maxETA
const double m_maxETA
Definition: JetVertexChecker.cc:69
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
reco::JetTracksAssociationCollection
JetTracksAssociation::Container JetTracksAssociationCollection
typedefs for backward compatibility
Definition: JetTracksAssociation.h:60
edm::ParameterSet
Definition: ParameterSet.h:47
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
Event.h
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:176
JetVertexChecker::m_cutMinPt
const double m_cutMinPt
Definition: JetVertexChecker.cc:60
JetVertexChecker::m_doFilter
const bool m_doFilter
Definition: JetVertexChecker.cc:59
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::AssociationVector::const_iterator
transient_vector_type::const_iterator const_iterator
Definition: AssociationVector.h:106
JetVertexChecker::m_associator
const edm::EDGetTokenT< reco::JetTracksAssociationCollection > m_associator
Definition: JetVertexChecker.cc:57
EgHLTOffHistBins_cfi.et
et
Definition: EgHLTOffHistBins_cfi.py:8
edm::EventSetup
Definition: EventSetup.h:58
reco::TrackBase::chi2
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
InputTag.h
VertexFwd.h
JetVertexChecker::JetVertexChecker
JetVertexChecker(const edm::ParameterSet &)
Definition: JetVertexChecker.cc:86
JetVertexChecker::filter
bool filter(edm::Event &, const edm::EventSetup &) override
Definition: JetVertexChecker.cc:116
JetVertexChecker::m_pvErr_x
const double m_pvErr_x
Definition: JetVertexChecker.cc:70
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
reco::Vertex::Point
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
eostools.move
def move(src, dest)
Definition: eostools.py:511
Vertex.h
Frameworkfwd.h
JetVertexChecker::m_cutMinPtRatio
const double m_cutMinPtRatio
Definition: JetVertexChecker.cc:61
edm::RefVectorIterator
Definition: EDProductfwd.h:33
EDFilter.h
mps_fire.result
result
Definition: mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
JetVertexChecker::~JetVertexChecker
~JetVertexChecker() override
Definition: JetVertexChecker.cc:106
reco::TrackBase::momentum
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
edm::stream::EDFilter
Definition: EDFilter.h:36
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
JetVertexChecker::m_pvErr_z
const double m_pvErr_z
Definition: JetVertexChecker.cc:72
reco::Vertex
Definition: Vertex.h:35
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
JetVertexChecker::m_maxChi2
const double m_maxChi2
Definition: JetVertexChecker.cc:63