Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023
00024
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDFilter.h"
00027
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "DataFormats/TrackReco/interface/Track.h"
00033 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
00034 #include "FWCore/Utilities/interface/InputTag.h"
00035
00036 #include "DataFormats/JetReco/interface/CaloJet.h"
00037 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00038
00039 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00040 #include "DataFormats/VertexReco/interface/Vertex.h"
00041
00042 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00043
00044
00045
00046
00047
00048 class JetVertexChecker : public edm::EDFilter {
00049 public:
00050 explicit JetVertexChecker(const edm::ParameterSet&);
00051 ~JetVertexChecker();
00052
00053 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
00054
00055 private:
00056 virtual bool filter(edm::Event&, const edm::EventSetup&) override;
00057
00058
00059 edm::InputTag m_associator;
00060 edm::InputTag m_primaryVertexProducer;
00061 edm::InputTag m_beamSpot;
00062 bool m_doFilter;
00063 double m_cutMinPt;
00064 double m_cutMinPtRatio;
00065 int32_t m_maxNjets;
00066
00067 };
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080 JetVertexChecker::JetVertexChecker(const edm::ParameterSet& iConfig)
00081 {
00082
00083 m_beamSpot = iConfig.getParameter<edm::InputTag>("beamSpot");
00084 m_associator = iConfig.getParameter<edm::InputTag>("jetTracks");
00085 m_doFilter = iConfig.getParameter<bool>("doFilter");
00086 m_cutMinPt = iConfig.getParameter<double>("minPt");
00087 m_cutMinPtRatio = iConfig.getParameter<double>("minPtRatio");
00088 m_maxNjets = iConfig.getParameter<int32_t>("maxNJetsToCheck");
00089 produces<std::vector<reco::CaloJet> >();
00090 produces<reco::VertexCollection >();
00091 }
00092
00093
00094 JetVertexChecker::~JetVertexChecker()
00095 {
00096
00097
00098
00099
00100 }
00101
00102
00103
00104
00105
00106
00107
00108 bool
00109 JetVertexChecker::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00110 {
00111 using namespace edm;
00112 Handle<reco::JetTracksAssociationCollection> jetTracksAssociation;
00113 iEvent.getByLabel(m_associator, jetTracksAssociation);
00114 std::auto_ptr<std::vector<reco::CaloJet> > pOut(new std::vector<reco::CaloJet> );
00115
00116 bool result=true;
00117 int i = 0;
00118
00119 for(reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin();
00120 it != jetTracksAssociation->end() && i < m_maxNjets; it++, i++) {
00121 if(fabs(it->first->eta()) < 2.4)
00122 {
00123 reco::TrackRefVector tracks = it->second;
00124 math::XYZVector jetMomentum = it->first->momentum();
00125 math::XYZVector trMomentum;
00126 for(reco::TrackRefVector::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack)
00127 {
00128 trMomentum += (*itTrack)->momentum();
00129 }
00130 if(trMomentum.rho()/jetMomentum.rho() < m_cutMinPtRatio || trMomentum.rho() < m_cutMinPt)
00131 {
00132
00133 pOut->push_back(* dynamic_cast<const reco::CaloJet *>(&(*it->first)));
00134 result=false;
00135 }
00136 }
00137 }
00138
00139 iEvent.put(pOut);
00140
00141 edm::Handle<reco::BeamSpot> beamSpot;
00142 iEvent.getByLabel(m_beamSpot,beamSpot);
00143
00144 reco::Vertex::Error e;
00145 e(0, 0) = 0.0015 * 0.0015;
00146 e(1, 1) = 0.0015 * 0.0015;
00147 e(2, 2) = 1.5 * 1.5;
00148 reco::Vertex::Point p(beamSpot->x0(), beamSpot->y0(), beamSpot->z0());
00149 reco::Vertex thePV(p, e, 0, 0, 0);
00150 std::auto_ptr<reco::VertexCollection> pOut2(new reco::VertexCollection);
00151 pOut2->push_back(thePV);
00152 iEvent.put(pOut2);
00153
00154 if(m_doFilter) return result;
00155 else
00156 return true;
00157 }
00158
00159
00160 void
00161 JetVertexChecker::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00162
00163
00164 edm::ParameterSetDescription desc;
00165 desc.setUnknown();
00166 descriptions.addDefault(desc);
00167 }
00168
00169 DEFINE_FWK_MODULE(JetVertexChecker);