File indexing completed on 2024-04-06 12:24:11
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024 #include "FWCore/Framework/interface/stream/EDProducer.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0027 #include "FWCore/Utilities/interface/InputTag.h"
0028
0029 #include "DataFormats/Math/interface/deltaR.h"
0030 #include "DataFormats/Common/interface/ValueMap.h"
0031 #include "DataFormats/Common/interface/View.h"
0032
0033 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0034 #include "DataFormats/Candidate/interface/Candidate.h"
0035 #include "DataFormats/MuonReco/interface/Muon.h"
0036 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0037 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0038 #include "DataFormats/JetReco/interface/Jet.h"
0039 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0040
0041 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0042
0043 template <typename T>
0044 class DeltaRNearestObjectComputer : public edm::stream::EDProducer<> {
0045 public:
0046 explicit DeltaRNearestObjectComputer(const edm::ParameterSet& iConfig);
0047 ~DeltaRNearestObjectComputer() override;
0048
0049 void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0050
0051 private:
0052 edm::EDGetTokenT<edm::View<reco::Candidate>> probesToken_;
0053 edm::EDGetTokenT<edm::View<T>> objectsToken_;
0054 StringCutObjectSelector<T, true> objCut_;
0055 };
0056
0057 template <typename T>
0058 DeltaRNearestObjectComputer<T>::DeltaRNearestObjectComputer(const edm::ParameterSet& iConfig)
0059 : probesToken_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("probes"))),
0060 objectsToken_(consumes<edm::View<T>>(iConfig.getParameter<edm::InputTag>("objects"))),
0061 objCut_(
0062 iConfig.existsAs<std::string>("objectSelection") ? iConfig.getParameter<std::string>("objectSelection") : "",
0063 true) {
0064 produces<edm::ValueMap<float>>();
0065 }
0066
0067 template <typename T>
0068 DeltaRNearestObjectComputer<T>::~DeltaRNearestObjectComputer() {}
0069
0070 template <typename T>
0071 void DeltaRNearestObjectComputer<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0072 using namespace edm;
0073
0074
0075 Handle<View<reco::Candidate>> probes;
0076 iEvent.getByToken(probesToken_, probes);
0077
0078 Handle<View<T>> objects;
0079 iEvent.getByToken(objectsToken_, objects);
0080
0081
0082 std::vector<float> values;
0083
0084
0085 View<reco::Candidate>::const_iterator probe, endprobes = probes->end();
0086 for (probe = probes->begin(); probe != endprobes; ++probe) {
0087 double dr2min = 10000;
0088 for (unsigned int iObj = 0; iObj < objects->size(); iObj++) {
0089 const T& obj = objects->at(iObj);
0090 if (!objCut_(obj))
0091 continue;
0092 double dr2 = deltaR2(*probe, obj);
0093 if (dr2 < dr2min) {
0094 dr2min = dr2;
0095 }
0096 }
0097 values.push_back(sqrt(dr2min));
0098 }
0099
0100
0101 auto valMap = std::make_unique<ValueMap<float>>();
0102 ValueMap<float>::Filler filler(*valMap);
0103 filler.insert(probes, values.begin(), values.end());
0104 filler.fill();
0105 iEvent.put(std::move(valMap));
0106 }
0107
0108
0109
0110
0111
0112 typedef DeltaRNearestObjectComputer<reco::Candidate> DeltaRNearestCandidateComputer;
0113 typedef DeltaRNearestObjectComputer<reco::Muon> DeltaRNearestMuonComputer;
0114 typedef DeltaRNearestObjectComputer<reco::Electron> DeltaRNearestElectronComputer;
0115 typedef DeltaRNearestObjectComputer<reco::GsfElectron> DeltaRNearestGsfElectronComputer;
0116 typedef DeltaRNearestObjectComputer<reco::Photon> DeltaRNearestPhotonComputer;
0117 typedef DeltaRNearestObjectComputer<reco::Jet> DeltaRNearestJetComputer;
0118
0119 #include "FWCore/Framework/interface/MakerMacros.h"
0120 DEFINE_FWK_MODULE(DeltaRNearestCandidateComputer);
0121 DEFINE_FWK_MODULE(DeltaRNearestMuonComputer);
0122 DEFINE_FWK_MODULE(DeltaRNearestElectronComputer);
0123 DEFINE_FWK_MODULE(DeltaRNearestGsfElectronComputer);
0124 DEFINE_FWK_MODULE(DeltaRNearestPhotonComputer);
0125 DEFINE_FWK_MODULE(DeltaRNearestJetComputer);