Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #include <boost/bind.hpp>
00012
00013 #include "DataFormats/JetReco/interface/PFJet.h"
00014 #include "DataFormats/Common/interface/Association.h"
00015
00016 #include "RecoTauTag/RecoTau/interface/ConeTools.h"
00017 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00018
00019 #include "FWCore/Framework/interface/EDProducer.h"
00020 #include "FWCore/Framework/interface/EventSetup.h"
00021 #include "FWCore/Framework/interface/ESHandle.h"
00022 #include "FWCore/Framework/interface/Event.h"
00023 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025
00026 class RecoTauJetRegionProducer : public edm::EDProducer {
00027 public:
00028 typedef edm::Association<reco::PFJetCollection> PFJetMatchMap;
00029 explicit RecoTauJetRegionProducer(const edm::ParameterSet& pset);
00030 ~RecoTauJetRegionProducer() {}
00031 void produce(edm::Event& evt, const edm::EventSetup& es);
00032 private:
00033 float deltaR2_;
00034 edm::InputTag inputJets_;
00035 edm::InputTag pfSrc_;
00036 };
00037
00038 RecoTauJetRegionProducer::RecoTauJetRegionProducer(
00039 const edm::ParameterSet& pset) {
00040 deltaR2_ = pset.getParameter<double>("deltaR"); deltaR2_*=deltaR2_;
00041 inputJets_ = pset.getParameter<edm::InputTag>("src");
00042 pfSrc_ = pset.getParameter<edm::InputTag>("pfSrc");
00043 produces<reco::PFJetCollection>("jets");
00044 produces<PFJetMatchMap>();
00045 }
00046
00047 void RecoTauJetRegionProducer::produce(edm::Event& evt,
00048 const edm::EventSetup& es) {
00049
00050 edm::Handle<reco::PFCandidateCollection> pfCandsHandle;
00051 evt.getByLabel(pfSrc_, pfCandsHandle);
00052
00053
00054 typedef edm::Ptr<reco::PFCandidate> PFCandPtr;
00055 std::vector<PFCandPtr> pfCands;
00056 pfCands.reserve(pfCandsHandle->size());
00057 for (size_t icand = 0; icand < pfCandsHandle->size(); ++icand) {
00058 pfCands.push_back(PFCandPtr(pfCandsHandle, icand));
00059 }
00060
00061
00062 edm::Handle<reco::CandidateView> jetView;
00063 evt.getByLabel(inputJets_, jetView);
00064
00065 reco::PFJetRefVector jets =
00066 reco::tau::castView<reco::PFJetRefVector>(jetView);
00067 size_t nJets = jets.size();
00068
00069
00070
00071 edm::ProductID originalId = jets.id();
00072 edm::Handle<reco::PFJetCollection> originalJets;
00073 size_t nOriginalJets = 0;
00074
00075
00076 if (nJets) {
00077 try {
00078 evt.get(originalId, originalJets);
00079 } catch(const cms::Exception &e) {
00080 edm::LogError("MissingOriginalCollection")
00081 << "Can't get the original jets that made: " << inputJets_
00082 << " that have product ID: " << originalId
00083 << " from the event!!";
00084 throw e;
00085 }
00086 nOriginalJets = originalJets->size();
00087 }
00088
00089 std::auto_ptr<reco::PFJetCollection> newJets(new reco::PFJetCollection);
00090
00091
00092
00093 std::vector<int> matchInfo(nOriginalJets, -1);
00094 newJets->reserve(nJets);
00095 for (size_t ijet = 0; ijet < nJets; ++ijet) {
00096
00097 reco::PFJetRef jetRef = jets[ijet];
00098
00099 newJets->emplace_back(*jetRef);
00100 reco::PFJet & newJet = newJets->back();
00101
00102 newJet.clearDaughters();
00103
00104 for ( auto cand : pfCands )
00105 if ( reco::deltaR2(*jetRef,*cand)<deltaR2_ ) newJet.addDaughter(cand);
00106
00107
00108 matchInfo[jetRef.key()] = ijet;
00109 }
00110
00111
00112 edm::OrphanHandle<reco::PFJetCollection> newJetsInEvent =
00113 evt.put(newJets, "jets");
00114
00115
00116 std::auto_ptr<PFJetMatchMap> matching(new PFJetMatchMap(newJetsInEvent));
00117 if (nJets) {
00118 PFJetMatchMap::Filler filler(*matching);
00119 filler.insert(originalJets, matchInfo.begin(), matchInfo.end());
00120 filler.fill();
00121 }
00122 evt.put(matching);
00123 }
00124
00125 #include "FWCore/Framework/interface/MakerMacros.h"
00126 DEFINE_FWK_MODULE(RecoTauJetRegionProducer);