File indexing completed on 2024-09-07 04:38:05
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDFilter.h"
0025
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "DataFormats/TrackReco/interface/Track.h"
0031 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
0032 #include "FWCore/Utilities/interface/InputTag.h"
0033
0034 #include "DataFormats/JetReco/interface/CaloJet.h"
0035 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0036
0037 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0038 #include "DataFormats/VertexReco/interface/Vertex.h"
0039
0040 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0041
0042
0043
0044
0045
0046 class JetVertexChecker : public edm::stream::EDFilter<> {
0047 public:
0048 explicit JetVertexChecker(const edm::ParameterSet&);
0049 ~JetVertexChecker() override;
0050
0051 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0052
0053 private:
0054 bool filter(edm::Event&, const edm::EventSetup&) override;
0055
0056
0057 const edm::EDGetTokenT<reco::JetTracksAssociationCollection> m_associator;
0058 const edm::EDGetTokenT<reco::BeamSpot> m_beamSpot;
0059 const bool m_doFilter;
0060 const double m_cutMinPt;
0061 const double m_cutMinPtRatio;
0062 const double m_maxTrackPt;
0063 const double m_maxChi2;
0064 const int32_t m_maxNjets;
0065 const int32_t m_maxNjetsOutput;
0066
0067 const bool m_newMethod;
0068
0069 const double m_maxETA;
0070 const double m_pvErr_x;
0071 const double m_pvErr_y;
0072 const double m_pvErr_z;
0073 };
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085
0086 JetVertexChecker::JetVertexChecker(const edm::ParameterSet& iConfig)
0087 : m_associator(consumes<reco::JetTracksAssociationCollection>(iConfig.getParameter<edm::InputTag>("jetTracks"))),
0088 m_beamSpot(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))),
0089 m_doFilter(iConfig.getParameter<bool>("doFilter")),
0090 m_cutMinPt(iConfig.getParameter<double>("minPt")),
0091 m_cutMinPtRatio(iConfig.getParameter<double>("minPtRatio")),
0092 m_maxTrackPt(iConfig.getParameter<double>("maxTrackPt")),
0093 m_maxChi2(iConfig.getParameter<double>("maxChi2")),
0094 m_maxNjets(iConfig.getParameter<int32_t>("maxNJetsToCheck")),
0095 m_maxNjetsOutput(iConfig.getParameter<int32_t>("maxNjetsOutput")),
0096 m_newMethod(iConfig.getParameter<bool>("newMethod")),
0097 m_maxETA(iConfig.getParameter<double>("maxETA")),
0098 m_pvErr_x(iConfig.getParameter<double>("pvErr_x")),
0099 m_pvErr_y(iConfig.getParameter<double>("pvErr_y")),
0100 m_pvErr_z(iConfig.getParameter<double>("pvErr_z")) {
0101
0102 produces<std::vector<reco::CaloJet>>();
0103 produces<reco::VertexCollection>();
0104 }
0105
0106 JetVertexChecker::~JetVertexChecker() {
0107
0108
0109 }
0110
0111
0112
0113
0114
0115
0116 bool JetVertexChecker::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0117 using namespace edm;
0118 Handle<reco::JetTracksAssociationCollection> jetTracksAssociation;
0119 iEvent.getByToken(m_associator, jetTracksAssociation);
0120 auto pOut = std::make_unique<std::vector<reco::CaloJet>>();
0121
0122 bool result = true;
0123 int i = 0;
0124
0125 for (reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin(),
0126 et = jetTracksAssociation->end();
0127 it != et && i < m_maxNjets;
0128 it++, i++) {
0129 if (std::abs(it->first->eta()) < m_maxETA) {
0130 reco::TrackRefVector tracks = it->second;
0131 math::XYZVector jetMomentum = it->first->momentum();
0132 math::XYZVector trMomentum;
0133 for (reco::TrackRefVector::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack) {
0134 const reco::Track& iTrack = **itTrack;
0135 if (m_newMethod && iTrack.chi2() > m_maxChi2)
0136 continue;
0137 trMomentum += iTrack.momentum();
0138 }
0139 if (trMomentum.rho() / jetMomentum.rho() < m_cutMinPtRatio || trMomentum.rho() < m_cutMinPt) {
0140 pOut->push_back(*dynamic_cast<const reco::CaloJet*>(&(*it->first)));
0141 }
0142 }
0143 }
0144 iEvent.put(std::move(pOut));
0145
0146 edm::Handle<reco::BeamSpot> beamSpot;
0147 iEvent.getByToken(m_beamSpot, beamSpot);
0148
0149 reco::Vertex::Error e;
0150 e(0, 0) = m_pvErr_x * m_pvErr_x;
0151 e(1, 1) = m_pvErr_y * m_pvErr_y;
0152 e(2, 2) = m_pvErr_z * m_pvErr_z;
0153 reco::Vertex::Point p(beamSpot->x0(), beamSpot->y0(), beamSpot->z0());
0154 reco::Vertex thePV(p, e, 0, 0, 0);
0155 auto pOut2 = std::make_unique<reco::VertexCollection>();
0156 pOut2->push_back(thePV);
0157 iEvent.put(std::move(pOut2));
0158
0159 if (m_doFilter)
0160 return result;
0161 else
0162 return true;
0163 }
0164
0165
0166 void JetVertexChecker::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0167 edm::ParameterSetDescription desc;
0168 desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
0169 desc.add<edm::InputTag>("jetTracks", edm::InputTag("hltFastPVJetTracksAssociator"));
0170 desc.add<double>("minPtRatio", 0.1);
0171 desc.add<double>("minPt", 0.0);
0172 desc.add<bool>("doFilter", false);
0173 desc.add<int>("maxNJetsToCheck", 2);
0174 desc.add<int>("maxNjetsOutput", 2);
0175 desc.add<double>("maxChi2", 20.0);
0176 desc.add<double>("maxTrackPt", 20.0);
0177 desc.add<bool>("newMethod", false);
0178 desc.add<double>("maxETA", 2.4);
0179 desc.add<double>("pvErr_x", 0.0015);
0180 desc.add<double>("pvErr_y", 0.0015);
0181 desc.add<double>("pvErr_z", 1.5);
0182 descriptions.add("jetVertexChecker", desc);
0183 }
0184
0185 DEFINE_FWK_MODULE(JetVertexChecker);